Fisher¶
Base Class¶
Definition of the FisherForecast interface and a number of general utility functions.
-
class
fisher.fisher_forecast.
FisherForecast
[source]¶ Class that collects and contains information for making Fisher-matrix type forecasts for observations of various types. Forms a base class for fully analytical versions that may be application specific.
-
covar
¶ Internal space for the computation of the covariance matrix.
- Type
-
argument_hash
¶ MD5 hash object indicating last state of covariance computation. Used to determine if the covariance needs to be recomputed.
- Type
-
add_gaussian_prior
(pindex, sigma, verbosity=0)[source]¶ Add a Gaussian prior on the parameter with the specified index.
-
add_gaussian_prior_list
(sigma_list, verbosity=0)[source]¶ Add a Gaussian priors on all of the model parameters.
-
check_gradients
(obs, p, h=None, verbosity=0)[source]¶ Numerically evaluates the gradients using the visibilities function and compares them to the gradients returned by the visibility_gradients function. Numerical differentiation is 2nd-order, centered finite differences. Results are output to standard out. For best results, the user should pass some option for the step size h to use for numerical differentiation.
- Parameters
obs (ehtim.obsdata.Obsdata) – An ehtim Obsdata object with a particular set of observation details (u,v positions, times, etc.) at which to check gradients..
p (list) – List of parameter values.
h (float,list) – List of small steps which define the finite differences. If None, uses 1e-4*p. If a float, uses h*p. If a list, sets each step size separately.
verbosity (int) – Verbosity level. Default: 0.
-
display_image
(p, limits=None, shape=None, verbosity=0, **kwargs)[source]¶ Generate and plot an image for the parent visibility model evaluated at a given set of parameter values. The user is responsible for setting the plot limits and shape intelligently. Note that this uses FFTs and assumes that the resolution and field of view are sufficient to adequately resolve model image features.
- Parameters
p (list) – List of parameter values.
limits (float,list) – Limits on the field of view in which to construct the image in uas. Either a single float, in which case the limits are symmetric and set to [-limits,limits,-limits,limits], or a list of floats. Default: 100.
shape (int,list) – Dimensions of the image. If an int, the dimensions are equal. If a two-element list of ints, sets the dimensions to the two values. Default: 256.
verbosity (int) – Verbosity level. Default: 0.
- Returns
Handles to the figure, axes, and colorbar.
- Return type
(matplotlib.pyplot.figure,matplotlib.pyplot.axes,matplotlib.pyplot.colorbar)
-
fisher_covar
(obs, p, verbosity=0, **kwargs)[source]¶ Returns the Fisher matrix as defined in the accompanying documentation. Intelligently avoids recomputation if the observation and parameters are unchanged.
- Parameters
obs (ehtim.obsdata.Obsdata) – An Obsdata object containing the observation particulars.
p (list) – List of model parameters at which to compute the Fisher matrix.
verbosity (int) – Verbosity level. Default: 0.
- Returns
The Fisher matrix.
- Return type
-
generate_image
(p, limits=None, shape=None, verbosity=0)[source]¶ Generate and return an image for the parent visibility model evaluated at a given set of parameter values. The user is responsible for setting the plot limits and shape intelligently. Note that this uses FFTs and assumes that the resolution and field of view are sufficient to adequately resolve model image features.
- Parameters
p (list) – List of parameter values.
limits (float,list) – Limits on the field of view in which to construct the image in uas. Either a single float, in which case the limits are symmetric and set to [-limits,limits,-limits,limits], or a list of floats. Default: 100.
shape (int,list) – Dimensions of the image. If an int, the dimensions are equal. If a two-element list of ints, sets the dimensions to the two values. Default: 256.
verbosity (int) – Verbosity level. Default: 0.
- Returns
x positions of the pixel centers in uas, y positions of the pixel centers in uas, intensities at pixel centers in Jy/uas^2.
- Return type
-
inverse_fisher_covar
(obs, p, verbosity=0, **kwargs)[source]¶ Returns the inverse of the Fisher matrix as defined in the accompanying documentation. Intelligently avoids recomputation if the observation and parameters are unchanged.
- Parameters
obs (ehtim.obsdata.Obsdata) – An Obsdata object containing the observation particulars.
p (list) – List of model parameters at which to compute the Fisher matrix.
verbosity (int) – Verbosity level. Default: 0.
- Returns
The inverse of the Fisher matrix.
- Return type
-
joint_biparameter_chisq
(obs, p, i1, i2, kind='marginalized', verbosity=0, **kwargs)[source]¶ Computes the ensemble-averaged 2nd-order contribution to the chi-square for two parameters after fixing or marginalizing over all others.
- Parameters
obs (ehtim.obsdata.Obsdata) – An Obsdata object containing the observation particulars.
p (list) – List of parameter values at which to compute uncertainties.
i1 (int) – Index of first parameter.
i2 (int) – Index of second parameter.
kind (str) – Choice of what to do with other parameters. Choices are ‘marginalized’, ‘fixed’. Default: ‘marginalized’.
verbosity (int) – Verbosity level. Default: 0.
- Returns
\(p_{i1}\), \(p_{i2}\), and \(\chi^2\) on a grid of points that covers the inner \(4.5\Sigma\) region.
- Return type
-
marginalized_covariance
(obs, p, ilist=None, verbosity=0, **kwargs)[source]¶ Computes the covariance for a subset of model parameters, marginalized over all parameters outside of the subset.
- Parameters
obs (ehtim.obsdata.Obsdata) – An Obsdata object containing the observation particulars.
p (list) – List of parameter values at which to compute uncertainties.
ilist (int,list) – Index or list of indicies for which to compute the marginalized uncertainties. If None will return a list of all marginalized uncertainties. Default: None.
verbosity (int) – Verbosity level. Default: 0.
- Returns
Marginalized covariance for desired parameter subset.
- Return type
-
marginalized_uncertainties
(obs, p, ilist=None, verbosity=0, **kwargs)[source]¶ Computes the uncertainties on a subset of model parameters, marginalized over all others. Note that this is not the marginalized covariance; each uncertainty is for a single parameter, marginalizing out all others.
- Parameters
obs (ehtim.obsdata.Obsdata) – An Obsdata object containing the observation particulars.
p (list) – List of parameter values at which to compute uncertainties.
ilist (int,list) – Index or list of indicies for which to compute the marginalized uncertainties. If None will return a list of all marginalized uncertainties. Default: None.
verbosity (int) – Verbosity level. Default: 0.
- Returns
Marginalized parameter uncertainty or list of marginalized uncertainties of desired parameters.
- Return type
(float/list)
-
parameter_labels
(verbosity=0)[source]¶ User-defined function in child classes that returns a list of of parameter labels to be used in plotting.
-
plot_1d_forecast
(obs, p, i1, kind='marginalized', labels=None, fig=None, axes=None, colors=None, alphas=0.25, verbosity=0, **kwargs)[source]¶ Plot the marginalized posterior for a single parameter for a given observation.
- Parameters
obs (list,Obsdata) – An Obsdata object or list of Obsdata objects containing the observation particulars.
p (list) – List of parameter values at which to compute uncertainties.
i1 (int) – Index of parameter to be plotted.
kind (str) – Choice of what to do with other parameters. Choices are ‘marginalized’, ‘fixed’. Default: ‘marginalized’.
labels (list,str) – A list of labels for each Obsdata object. When fewer labels than observations are provided, the first set of observations are labeled. Default: None.
fig (matplotlib.figure.Figure) – Figure on which to place plot. If None, a new figure will be created. Default: None.
axes (matplotlib.axes.Axes) – Axes on which to place plot. If None, a new axes will be created. Default: None.
colors (list,str) – A color or list of colors for the plots. When fewer colors than observations are provided, they will be cycled through. If None, a default list will be used. Default: None.
alphas (list,float) – An alpha value or list of alpha values for the filled portion of the plots. Default: 0.25.
verbosity (int) – Verbosity level. Default: 0.
- Returns
Handles to the figure and axes objects in the plot.
- Return type
-
plot_2d_forecast
(obs, p, i1, i2, kind='marginalized', labels=None, fig=None, axes=None, colors=None, cmaps=None, alphas=0.75, verbosity=0, **kwargs)[source]¶ Plot the joint marginalized posterior for a pair of parameters for a given observation.
- Parameters
obs (list,Obsdata) – An Obsdata object or list of Obsdata objects containing the observation particulars.
p (list) – List of parameter values at which to compute uncertainties.
i1 (int) – Index of first parameter to be plotted.
i2 (int) – Index of second parameter to be plotted.
kind (str) – Choice of what to do with other parameters. Choices are ‘marginalized’, ‘fixed’. Default: ‘marginalized’.
labels (list,str) – A list of labels for each Obsdata object. When fewer labels than observations are provided, the first set of observations are labeled. Default: None.
fig (matplotlib.figure.Figure) – Figure on which to place plot. If None, a new figure will be created. Default: None.
axes (matplotlib.axes.Axes) – Axes on which to place plot. If None, a new axes will be created. Default: None.
colors (list,str) – A color or list of colors for the plots. When fewer colors than observations are provided, they will be cycled through. If None, a default list will be used. Default: None.
cmaps (list,matplotlib.colors.Colormap) – A colormap or list of colormaps for the plots. When fewer colormaps than observations are provided, they will be cycled through. If None, a default list will be used. Default: None.
alphas (list,float) – An alpha value or list of alpha values for the filled portion of the plots. Default: 0.75.
verbosity (int) – Verbosity level. Default: 0.
- Returns
Handles to the figure and axes objects in the plot.
- Return type
-
plot_triangle_forecast
(obs, p, ilist=None, kind='marginalized', labels=None, colors=None, cmaps=None, alphas=0.75, fig=None, axes=None, figsize=None, axis_location=None, verbosity=0, **kwargs)[source]¶ Generate a triangle plot of joint posteriors for a selected list of parameter indexes.
- Parameters
obs (list,Obsdata) – An Obsdata object or list of Obsdata objects containing the observation particulars.
p (list) – List of parameter values at which to compute uncertainties.
ilist (list) – List of indicies of parameters to include in triangle plot. If None, includes all parameters. Default: None.
kind (str) – Choice of what to do with other parameters. Choices are ‘marginalized’, ‘fixed’. Default: ‘marginalized’.
labels (list,str) – A list of labels for each Obsdata object. When fewer labels than observations are provided, the first set of observations are labeled. Default: None.
colors (list,str) – A color or list of colors for the plots. When fewer colors than observations are provided, they will be cycled through. If None, a default list will be used. Default: None.
cmaps (list,matplotlib.colors.Colormap) – A colormap or list of colormaps for the plots. When fewer colormaps than observations are provided, they will be cycled through. If None, a default list will be used. Default: None.
alphas (list,float) – An alpha value or list of alpha values for the filled portion of the plots. Default: 0.75.
fig (matplotlib.figure.Figure) – Figure on which to place plot. If None, a new figure will be created. Default: None.
axes (matplotlib.axes.Axes) – Axes on which to place plot. If None, a new axes will be created. Default: None.
figsize (list) – Figure size in inches. If None, attempts a guess. Default: None.
axis_location (list) – Location parameters for triangle plot region. If None, attempts a guess. Default: None.
verbosity (int) – Verbosity level. Default: 0.
- Returns
Handles to the figure and dictionary of axes objects in the plot.
- Return type
-
sample_posterior
(obs, p, N, ilist=None, verbosity=0, **kwargs)[source]¶ Generates a chain of samples appropriate for post-processing as done in MCMC analyses. Such chains may be useful for computing the posterior on complicated non-linear quantities.
- Parameters
obs (ehtim.obsdata.Obsdata) – An Obsdata object containing the observation particulars.
p (list) – List of parameter values at which to compute uncertainties.
N (int) – Number of samples to generate.
ilist (list) – List of parameters to sample, marginalizing over all others. If None, all parameters are sampled. Default: None.
- Returns
Two-dimensional array containing chain of parameter samples.
- Return type
-
uncertainties
(obs, p, ilist=None, verbosity=0, **kwargs)[source]¶ Computes the uncertainties on a subset of model parameters, both fixing and marginalizing over all others.
- Parameters
obs (ehtim.obsdata.Obsdata) – An Obsdata object containing the observation particulars.
p (list) – List of parameter values at which to compute uncertainties.
ilist (int,list) – Index or list of indicies for which to compute the marginalized uncertainties. If None will return a list of all marginalized uncertainties. Default: None.
verbosity (int) – Verbosity level. Default: 0.
- Returns
Parameter uncertainty or list of marginalized uncertainties of desired parameters, marginalized parameter uncertainty or list of marginalized uncertainties of desired parameters.
- Return type
(float/list,float/list)
-
uniparameter_uncertainties
(obs, p, ilist=None, verbosity=0, **kwargs)[source]¶ Computes the uncertainties on a subset of model parameters, fixing all others. Note that this is not a covariance; each uncertainty is for a single parameter, fixing all others. This is probably not what is wanted for forecasting, see marginalized_uncertainties.
- Parameters
obs (ehtim.obsdata.Obsdata) – An Obsdata object containing the observation particulars.
p (list) – List of parameter values at which to compute uncertainties.
ilist (int,list) – Index or list of indicies for which to compute the uncertainties. If None will return a list of all uncertainties. Default: None.
verbosity (int) – Verbosity level. Default: 0.
- Returns
Parameter uncertainty or list of marginalized uncertainties of desired parameters.
- Return type
(float/list)
-
visibilities
(obs, p, verbosity=0)[source]¶ User-defined function in child classes that generates visibilities associated with a given model image object.
- Parameters
obs (ehtim.obsdata.Obsdata) – An ehtim Obsdata object with a particular set of observation details (u,v positions, times, etc.)
p (numpy.ndarray) – list of parameters for the model image used to create object.
verbosity (int) – Verbosity level. Default: 0.
- Returns
List of complex visibilities computed at observations.
- Return type
-
visibility_gradients
(obs, p, verbosity=0, **kwargs)[source]¶ User-defined function in child classes that generates visibility gradients associated with a given model image object.
- Parameters
obs (ehtim.obsdata.Obsdata) – An ehtim Obsdata object with a particular set of observation details (u,v positions, times, etc.)
p (numpy.ndarray) – list of parameters for the model image used to create object.
verbosity (int) – Verbosity level. Default: 0.
- Returns
List of complex visibilities gradients computed at observations.
- Return type
-
Station Gain Mitigation¶
A set of FisherForecast childe classes that facilitate including gain uncertainties.
-
class
fisher.ff_complex_gains.
FF_complex_gains
(ff, marg_method='vMv')[source]¶ FisherForecast with complex gain reconstruction with multiple epochs. This is usually the FisherForecast object that should be used to investigate the impact of uncertain station gains.
- Parameters
ff (FisherForecast) – A FisherForecast object to which we wish to add gains.
marg_method (str) – Method used for intermediate marginalization. Options are ‘covar’ and ‘vMv’. Default: ‘covar’.
-
ff
¶ The FisherForecast object before gain reconstruction.
- Type
-
gain_epochs
¶ 2D array containing the start and end times for each gain solution epoch.
- Type
-
gain_amplitude_priors
¶ list of standard deviations of the log-normal priors on the gain amplitudes. Default: 10.
- Type
-
gain_phase_priors
¶ list of standard deviations on the normal priors on the gain phases. Default: 30.
- Type
-
gain_ratio_amplitude_priors
¶ For polarized gains, list of the log-normal priors on the gain amplitude ratios. Default: 1e-10.
- Type
-
gain_ratio_phase_priors
¶ For polarized gains, list of the normal priors on the gain phase differences. Default: 1e-10.
- Type
-
marg_method
¶ Method used for intermediate marginalization. Options are ‘covar’ and ‘vMv’. Based on preliminary tests, ‘vMv’ is more accurate for the gain marginalization, though this may not be the case for models with very many parameters.
- Type
-
fisher_covar
(obs, p, verbosity=0, **kwargs)[source]¶ Returns the Fisher matrix as defined in the accompanying documentation, marginalized over the complex station gains. Intelligently avoids recomputation if the observation and parameters are unchanged.
- Parameters
obs (ehtim.obsdata.Obsdata) – An Obsdata object containing the observation particulars.
p (list) – List of model parameters at which to compute the Fisher matrix.
verbosity (int) – Verbosity level. Default: 0.
- Returns
The Fisher matrix.
- Return type
-
set_gain_amplitude_prior
(sigma, station=None, verbosity=0)[source]¶ Sets the log-normal priors on the gain amplitudes, either for a specified station or for all stations. Identical priors are assummed across the entire observation.
-
set_gain_epochs
(scans=False, gain_epochs=None)[source]¶ Sets the gain solution intervals (gain epochs) to be used. If neither scans nor gain_epochs selected, will solve for gains on each unique timestamp in the data.
- Parameters
scans (bool) – If True, solves for independent gains by scan. Overrides explicit specification of gain solution intervals. Default: False.
gain_epochs (nd.array) – 2D array containing the start and end times for each gain solution epoch. Default: None.
-
set_gain_phase_prior
(sigma, station=None, verbosity=0)[source]¶ Sets the normal priors on the gain phases, either for a specified station or for all stations. Usually, this should be the default (uniformative). Identical priors are assummed across the entire observation.
-
set_gain_ratio_amplitude_prior
(sigma, station=None, verbosity=0)[source]¶ Sets the log-normal priors on the R/L gain amplitude ratios, either for a specified station or for all stations. Only relevant for polarized models. Default: 1e-10 for all stations. Identical priors are assummed across the entire observation.
-
set_gain_ratio_phase_prior
(sigma, station=None, verbosity=0)[source]¶ Sets the normal priors on the R/L gain phase differences, either for a specified station or for all stations. Only relevant for polarized models. Default: 1e-10 for all stations. Identical priors are assummed across the entire observation.
-
visibilities
(obs, p, verbosity=0, **kwargs)[source]¶ Complex visibilities associated with the underlying FisherForecast object evaluated at the data points in the given Obsdata object for the model with the given parameter values.
- Parameters
obs (ehtim.obsdata.Obsdata) – An ehtim Obsdata object with a particular set of observation details (u,v positions, times, etc.)
p (numpy.ndarray) – list of parameters for the model at which visiblities are desired.
verbosity (int) – Verbosity level. Default: 0.
- Returns
List of complex visibilities computed at observations.
- Return type
-
visibility_gradients
(obs, p, verbosity=0, **kwargs)[source]¶ Gradients of the complex visibilities associated with the underlying FisherForecast object with respect to the model pararmeters evaluated at the data points in the given Obsdata object for the model with the given parameter values.
- Parameters
obs (ehtim.obsdata.Obsdata) – An ehtim Obsdata object with a particular set of observation details (u,v positions, times, etc.)
p (numpy.ndarray) – list of parameters for the model image used to create object.
verbosity (int) – Verbosity level. Default: 0.
- Returns
List of complex visibility gradients computed at observations.
- Return type
-
class
fisher.ff_complex_gains.
FF_complex_gains_single_epoch
(ff)[source]¶ FisherForecast with complex gain reconstruction assuming a single epoch. This is a helper class, and probably not the FisherForecast object that should be used to study the impact of gains. See FF_complex_gains.
- Parameters
ff (FisherForecast) – A FisherForecast object to which we wish to add gains.
-
ff
¶ The FisherForecast object before gain reconstruction.
- Type
-
gain_amplitude_priors
¶ list of standard deviations of the log-normal priors on the gain amplitudes. Default: 10.
- Type
-
gain_phase_priors
¶ list of standard deviations on the normal priors on the gain phases. Default: 30.
- Type
-
gain_ratio_amplitude_priors
¶ For polarized gains, list of the log-normal priors on the gain amplitude ratios. Default: 1e-10.
- Type
-
gain_ratio_phase_priors
¶ For polarized gains, list of the normal priors on the gain phase differences. Default: 1e-10.
- Type
-
parameter_labels
(verbosity=0)[source]¶ Parameter labels to be used in plotting, including those associated with the gains.
-
set_gain_amplitude_prior
(sigma, station=None)[source]¶ Sets the log-normal priors on the gain amplitudes, either for a specified station or for all stations.
-
set_gain_phase_prior
(sigma, station=None)[source]¶ Sets the normal priors on the gain phases, either for a specified station or for all stations. Usually, this should be the default (uniformative).
-
set_gain_ratio_amplitude_prior
(sigma, station=None)[source]¶ Sets the log-normal priors on the R/L gain amplitude ratios, either for a specified station or for all stations. Only relevant for polarized models. Default: 1e-10 for all stations.
-
set_gain_ratio_phase_prior
(sigma, station=None)[source]¶ Sets the normal priors on the R/L gain phase differences, either for a specified station or for all stations. Only relevant for polarized models. Default: 1e-10 for all stations.
-
visibilities
(obs, p, verbosity=0, **kwargs)[source]¶ Complex visibilities associated with the underlying FisherForecast object evaluated at the data points in the given Obsdata object for the model with the given parameter values.
- Parameters
obs (ehtim.obsdata.Obsdata) – An ehtim Obsdata object with a particular set of observation details (u,v positions, times, etc.)
p (numpy.ndarray) – list of parameters for the model at which visiblities are desired.
verbosity (int) – Verbosity level. Default: 0.
- Returns
List of complex visibilities computed at observations.
- Return type
-
visibility_gradients
(obs, p, verbosity=0, **kwargs)[source]¶ Gradients of the complex visibilities associated with the underlying FisherForecast object with respect to the model pararmeters evaluated at the data points in the given Obsdata object for the model with the given parameter values, including gradients with respect to the gain amplitudes, phases, and if this is a polarized model, R/L gain amplitude ratios and phase differences.
- Parameters
obs (ehtim.obsdata.Obsdata) – An ehtim Obsdata object with a particular set of observation details (u,v positions, times, etc.)
p (numpy.ndarray) – list of parameters for the model image used to create object.
verbosity (int) – Verbosity level. Default: 0.
- Returns
List of complex visibility gradients computed at observations.
- Return type
Models¶
Specific models implemented in ngEHTforecast as FisherForecast child classes. Each inherits all of the functionality of FisherForecast, and may be used anywhere a FisherForecast object is required.
-
class
fisher.ff_models.
FF_asymmetric_gaussian
(stokes='I')[source]¶ FisherForecast object for a noncircular Gaussian. Parameter vector is:
p[0] … Total flux in Jy.
p[1] … Symmetrized mean of the FHWM in the major and minor axes in uas: \({\rm FWHM}^{-2} = {\rm FMHM}_{\rm min}^{-2}+{\rm FWHM}_{\rm max}^{-2}\)
p[2] … Asymmetry parameter, \(A\), expected to be in the range [0,1).
p[3] … Position angle in radians of the major axis E of N.
In terms of these, \({\rm FWHM}_{\rm maj}={\rm FWHM}/\sqrt{1-A}\) and \({\rm FWHM}_{\rm min}={\rm FWHM}/\sqrt{1-A}\).
- Parameters
stokes (str) – Indicates if this is limited to Stokes I (‘I’) or include full polarization (‘full’). Default: ‘I’.
-
visibilities
(obs, p, verbosity=0)[source]¶ Generates visibilities associated with a noncircular Gaussian.
- Parameters
obs (ehtim.Obsdata) – ehtim data object
p (numpy.ndarray) – list of parameters for the model image used to create object.
verbosity (int) – Verbosity level. Default: 0.
- Returns
List of complex visibilities computed at observations.
- Return type
-
visibility_gradients
(obs, p, verbosity=0)[source]¶ Generates visibility gradients associated with a noncircular Gaussian.
- Parameters
obs (ehtim.Obsdata) – ehtim data object
p (numpy.ndarray) – list of parameters for the model image used to create object.
verbosity (int) – Verbosity level. Default: 0.
- Returns
List of complex visibilities computed at observations.
- Return type
-
class
fisher.ff_models.
FF_double_smoothed_delta_ring
(stokes='I')[source]¶ FisherForecast object for a delta-ring convolved with a circular Gaussian. Parameter vector is:
p[0] … Flux of n=1 ring in Jy.
p[1] … Diameter of n=1 ring in uas.
p[2] … Twice the standard deviation of the Gaussian smoothing kernel in uas.
p[3] … Ratio of flux in n=2 and n=1 rings.
p[4] … Ratio of difference in the diameters between the n=1 and n=2 rings to the diameter of the n=1 ring.
p[5] … Ratio of the standard deviations of the n=2 and n=1 rings.
- Parameters
stokes (str) – Indicates if this is limited to Stokes I (‘I’) or include full polarization (‘full’). Default: ‘I’.
-
visibilities
(obs, p, verbosity=0)[source]¶ Generates visibilities associated with Gaussian-convolved delta-ring.
- Parameters
obs (ehtim.Obsdata) – ehtim data object
p (numpy.ndarray) – list of parameters for the model image used to create object.
verbosity (int) – Verbosity level. Default: 0.
- Returns
List of complex visibilities computed at observations.
- Return type
-
visibility_gradients
(obs, p, verbosity=0)[source]¶ Generates visibility gradients associated with Gaussian-convolved delta-ring.
- Parameters
obs (ehtim.Obsdata) – ehtim data object
p (numpy.ndarray) – list of parameters for the model image used to create object.
verbosity (int) – Verbosity level. Default: 0.
- Returns
List of complex visibilities computed at observations.
- Return type
-
class
fisher.ff_models.
FF_model_image
(img, stokes='I')[source]¶ FisherForecast from ThemisyPy model_image objects. Uses FFTs and centered finite-difference to compute the visibilities and gradients. May not always produce sensible behaviors. Has not been extensively tested for some time.
- Parameters
img (model_image) – A ThemisPy model_image object.
stokes (str) – Indicates if this is limited to Stokes I (‘I’) or include full polarization (‘full’). Default: ‘I’.
-
img
¶ The ThemisPy model_image object for which to make forecasts.
- Type
model_image
-
visibilities
(obs, p, limits=None, shape=None, padding=4, verbosity=0)[source]¶ Generates visibilities associated with a given model image object. Note that this uses FFTs to create the visibilities and then interpolates.
- Parameters
obs (ehtim.Obsdata) – ehtim data object
p (numpy.ndarray) – list of parameters for the model image used to create object.
limits (list) – May be single number, list with four elements specifying in \(\mu as\) the symmetric limits or explicit limits in [xmin,xmax,ymin,ymax] order. Default: 100.
shape (list) – May be single number or list with two elements indicating the number of pixels in the two directions. Default: 256.
padding (int) – Factor by which to pad image with zeros prior to FFT. Default: 4.
verbosity (int) – Verbosity level. 0 prints nothing. 1 prints various elements of the plotting process. Passed to
model_image.generate_intensity_map()
. Default: 0.
- Returns
List of complex visibilities computed at observations.
- Return type
-
visibility_gradients
(obs, p, h=0.01, limits=None, shape=None, padding=4, verbosity=0)[source]¶ Generates visibilities associated with a given model image object. Note that this uses FFTs to create the visibilities and then interpolates.
- Parameters
obs (ehtim.Obsdata) – ehtim data object
p (numpy.ndarray) – list of parameters for the model image used to create object.
h (float) – fractional step for finite-difference gradient computation. Default: 0.01.
limits (list) – May be single number, list with four elements specifying in \(\mu as\) the symmetric limits or explicit limits in [xmin,xmax,ymin,ymax] order. Default: 100.
shape (list) – May be single number or list with two elements indicating the number of pixels in the two directions. Default: 256.
padding (int) – Factor by which to pad image with zeros prior to FFT. Default: 4.
verbosity (int) – Verbosity level. 0 prints nothing. 1 prints various elements of the plotting process. Passed to
model_image.generate_intensity_map()
. Default: 0.
- Returns
List of complex visibilities gradients computed at observation.
- Return type
-
class
fisher.ff_models.
FF_smoothed_delta_ring
(stokes='I')[source]¶ FisherForecast object for a delta-ring convolved with a circular Gaussian. Parameter vector is:
p[0] … Total flux in Jy.
p[1] … Diameter in uas.
p[2] … Twice the standard deviation of the Gaussian smoothing kernel in uas.
- Parameters
stokes (str) – Indicates if this is limited to Stokes I (‘I’) or include full polarization (‘full’). Default: ‘I’.
-
visibilities
(obs, p, verbosity=0)[source]¶ Generates visibilities associated with Gaussian-convolved delta-ring.
- Parameters
obs (ehtim.Obsdata) – ehtim data object
p (numpy.ndarray) – list of parameters for the model image used to create object.
verbosity (int) – Verbosity level. Default: 0.
- Returns
List of complex visibilities computed at observations.
- Return type
-
visibility_gradients
(obs, p, verbosity=0)[source]¶ Generates visibility gradients associated with Gaussian-convolved delta-ring.
- Parameters
obs (ehtim.Obsdata) – ehtim data object
p (numpy.ndarray) – list of parameters for the model image used to create object.
verbosity (int) – Verbosity level. Default: 0.
- Returns
List of complex visibilities computed at observations.
- Return type
-
class
fisher.ff_models.
FF_smoothed_delta_ring_themage
(N, fov, stokes='I')[source]¶ FisherForecast object for a splined raster (i.e., themage) plus a Gaussian-convolved delta-ring.
Parameter vector is the log of the intensity at the control points:
p[0] ……. \(\ln(I[0,0])\)
p[1] ……. \(\ln(I[1,0])\)
…
p[N-1] ….. \(\ln(I[N-1,0])\)
p[N] ……. \(\ln(I[0,1])\)
…
p[N*N-1] … \(\ln(I[N-1,N-1])\)
p[N*N+0] … Total flux in Jy.
p[N*N+1] … Diameter in uas.
p[N*N+2] … Twice the standard deviation of the Gaussian smoothing kernel in uas.
where each \(I[i,j]\) is measured in Jy/sr.
- Parameters
-
xcp
¶ x-positions of raster control points.
- Type
-
ycp
¶ y-positions of raster control points.
- Type
-
visibilities
(obs, p, verbosity=0)[source]¶ Generates visibilities associated with splined raster + Gaussian-convolved delta-ring.
- Parameters
obs (ehtim.Obsdata) – ehtim data object
p (numpy.ndarray) – list of parameters for the model image used to create object.
verbosity (int) – Verbosity level. Default: 0.
- Returns
List of complex visibilities computed at observations.
- Return type
-
visibility_gradients
(obs, p, verbosity=0)[source]¶ Generates visibility gradients associated with splined raster + Gaussian-convolved delta-ring.
- Parameters
obs (ehtim.Obsdata) – ehtim data object
p (numpy.ndarray) – list of parameters for the model image used to create object.
verbosity (int) – Verbosity level. Default: 0.
- Returns
List of complex visibilities computed at observations.
- Return type
-
class
fisher.ff_models.
FF_splined_raster
(N, fov, stokes='I')[source]¶ FisherForecast object for a splined raster (i.e., themage). Parameter vector is the log of the intensity at the control points:
p[0] ……. \(\ln(I[0,0])\)
p[1] ……. \(\ln(I[1,0])\)
…
p[N-1] ….. \(\ln(I[N-1,0])\)
p[N] ……. \(\ln(I[0,1])\)
…
p[N*N-1] … \(\ln(I[N-1,N-1])\)
where each \(I[i,j]\) is measured in Jy/sr.
- Parameters
-
xcp
¶ x-positions of raster control points.
- Type
-
ycp
¶ y-positions of raster control points.
- Type
-
generate_parameter_list
(glob, verbosity=0, **kwargs)[source]¶ Utility function to quickly and easily generate a set of raster values associated with various potential initialization schemes. These include: an existing FisherForecast object, a FisherForecast class name, an
ehtim.image.Image
object, or the name of a FITS file. If initializing from a FITS file, the image will be blurred to the raster resolution prior to setting the parameter list.- Parameters
glob (str) – Can be an existing FisherForecast child object, the name of a FisherForecast child class,
ehtim.image.Image
object, or a string with the name of a FITS file.verbosity (int) – Verbosity level. Default: 0.
**kwargs (dict) – Additional key-word arguments necesary to define the object. May include parameter values for the relevant FisherForecast object or arguments to the
fisher.fisher_forecast.FisherForecast.generate_image()
function.
- Returns
Approximate parameter list for associated splined raster object.
- Return type
p (list)
-
visibilities
(obs, p, verbosity=0)[source]¶ Generates visibilities associated with a splined raster model.
- Parameters
obs (ehtim.Obsdata) – ehtim data object
p (numpy.ndarray) – list of parameters for the model image used to create object.
verbosity (int) – Verbosity level. Default: 0.
- Returns
List of complex visibilities computed at observations.
- Return type
-
visibility_gradients
(obs, p, verbosity=0)[source]¶ Generates visibility gradients associated with a splined raster model.
- Parameters
obs (ehtim.Obsdata) – ehtim data object
p (numpy.ndarray) – list of parameters for the model image used to create object.
verbosity (int) – Verbosity level. Default: 0.
- Returns
List of complex visibilities computed at observations.
- Return type
-
class
fisher.ff_models.
FF_symmetric_gaussian
(stokes='I')[source]¶ FisherForecast object for a circular Gaussian. Parameter vector is:
p[0] … Total flux in Jy.
p[1] … FWHM in uas.
- Parameters
stokes (str) – Indicates if this is limited to Stokes I (‘I’) or include full polarization (‘full’). Default: ‘I’.
-
visibilities
(obs, p, verbosity=0)[source]¶ Generates visibilities associated with a circular Gaussian.
- Parameters
obs (ehtim.Obsdata) – ehtim data object
p (numpy.ndarray) – list of parameters for the model image used to create object.
verbosity (int) – Verbosity level. Default: 0.
- Returns
List of complex visibilities computed at observations.
- Return type
-
visibility_gradients
(obs, p, verbosity=0)[source]¶ Generates visibility gradients associated with a circular Gaussian.
- Parameters
obs (ehtim.Obsdata) – ehtim data object
p (numpy.ndarray) – list of parameters for the model image used to create object.
verbosity (int) – Verbosity level. Default: 0.
- Returns
List of complex visibilities computed at observations.
- Return type
-
class
fisher.ff_models.
FF_thick_mring
(m, mp=None, mc=None, stokes='I')[source]¶ FisherForecast object for an m-ring model (based on ehtim). Parameter vector is:
p[0] … DOM, PLEASE FILL THIS IN.
- Parameters
m (int) – Stokes I azimuthal Fourier series order.
mp (int) – Linear polarization azimuthal Fourier series order. Only used if stokes==’full’. Default: None.
mc (int) – Circular polarization azimuthal Fourier series order. Only used if stokes==’full’. Default: None.
stokes (str) – Indicates if this is limited to Stokes I (‘I’) or include full polarization (‘full’). Default: ‘I’.
-
visibilities
(obs, p, verbosity=0)[source]¶ Generates visibilities associated with m-ring model.
- Parameters
obs (ehtim.Obsdata) – ehtim data object
p (numpy.ndarray) – list of parameters for the model image used to create object.
verbosity (int) – Verbosity level. Default: 0.
- Returns
List of complex visibilities computed at observations.
- Return type
-
visibility_gradients
(obs, p, verbosity=0)[source]¶ Generates visibility gradients associated with m-ring model.
- Parameters
obs (ehtim.Obsdata) – ehtim data object
p (numpy.ndarray) – list of parameters for the model image used to create object.
verbosity (int) – Verbosity level. Default: 0.
- Returns
List of complex visibilities computed at observations.
- Return type
-
fisher.ff_models.
W_cubic_spline_1d
(k, a=- 0.5)[source]¶ Fourier domain convolution function for the approximate 1D cubic spline. Useful for the splined raster image classes.
- Parameters
k (numpy.ndarray) – Wavenumber.
a (float) – Cubic spline control parameter. Default: -0.5.
- Returns
1D array of values of the cubic spline weight function evaluated at each value of k.
- Return type
Model Combinations¶
Classes that can be used to construct new FisherForecast objects from others. The prototypical example is fisher.ff_metamodels.FF_Sum
, which constructs a new model by directly summing the intensity maps of many previous models.
-
class
fisher.ff_metamodels.
FF_sum
(ff_list=None, stokes='I')[source]¶ FisherForecast object constructed from the sum of multiple FisherForecast objects. For example, a binary might be generated from the sum of two Gaussians.
The parameter vector is constructed from the concatenation of the parameter vectors from each individual object, with all objects after the first gaining a pair of offset parameters. That is:
p[0] ………… Obj1 p[0]
p[1] ………… Obj1 p[1]
…
p[n1] ……….. Obj2 p[0]
p[n1+1] ……… Obj2 p[0]
…
p[n1+n2] …….. Obj2 dx
p[n1+n2+1] …… Obj2 dy
p[n1+n2+2] …… Obj3 p[0]
…
p[n1+n2+n3+2] … Obj3 dx
p[n1+n2+n3+3] … Obj3 dy
…
Prior lists are constructed similarly.
- Parameters
-
add
(ff)[source]¶ Adds a FisherForecast object to the sum.
- Parameters
ff (FisherForecast) – An existing FisherForecast object to add.
-
visibilities
(obs, p, verbosity=0)[source]¶ Generates visibilities associated with the summed model.
- Parameters
obs (ehtim.Obsdata) – ehtim data object
p (numpy.ndarray) – list of parameters for the model image used to create object.
verbosity (int) – Verbosity level. Default: 0.
- Returns
List of complex visibilities computed at observations.
- Return type
-
visibility_gradients
(obs, p, verbosity=0)[source]¶ Generates visibility gradients associated with the summed model.
- Parameters
obs (ehtim.Obsdata) – ehtim data object
p (numpy.ndarray) – list of parameters for the model image used to create object.
verbosity (int) – Verbosity level. Default: 0.
- Returns
List of complex visibilities computed at observations.
- Return type