model_image_refractive_scattering Class Reference

Defines the interface for models that generates refractive scattered images. More...

#include "model/model_movie_refractive_scattering.h"

Inheritance diagram for model_image_refractive_scattering:
Collaboration diagram for model_image_refractive_scattering:

Classes

class  _P_phi
 Wandering magnetic field angular model. More...
 

Public Member Functions

 model_image_refractive_scattering (model_image &model, size_t nModes, double tobs, double frequency=230e9, std::string scattering_model="dipole", double observer_screen_distance=2.82 *3.086e21, double source_screen_distance=5.53 *3.086e21, double theta_maj_mas_cm=1.38, double theta_min_ma_cm=0.703, double POS_ANG=81.9, double scatt_alpha=1.38, double r_in=800e5, double r_out=1e20, double vs_ss_kms=50.0, double vy_ss_kms=0.0)
 
virtual size_t size () const
 
virtual std::string model_tag () const
 Return a string that contains a unique identifying tag for use with the ThemisPy plotting features.
 
virtual void generate_model (std::vector< double > parameters)
 Takes model parameters and generates scattered observables e.g. VA.
 
void set_image_resolution (size_t nray)
 
void set_screen_size (double fov)
 
void get_unscattered_image (std::vector< std::vector< double > > &alpha, std::vector< std::vector< double > > &beta, std::vector< std::vector< double > > &I) const
 
void get_ensemble_average_image (std::vector< std::vector< double > > &alpha, std::vector< std::vector< double > > &beta, std::vector< std::vector< double > > &I) const
 
void get_image (std::vector< std::vector< double > > &alpha, std::vector< std::vector< double > > &beta, std::vector< std::vector< double > > &I) const
 
virtual void set_mpi_communicator (MPI_Comm comm)
 
- Public Member Functions inherited from model_image
virtual void generate_complex_visibilities ()
 A one-time generate function that will generate the complex visibilities and store them. This must be called after generate_model has been called.
 
virtual std::complex< double > visibility (datum_visibility &d, double accuracy)
 Returns complex visibility in Jy computed from the image given a datum_visibility_amplitude object, containing all of the accoutrements. While this provides access to the actual data value, the two could be separated if necessary. Also takes an accuracy parameter with the same units as the data, indicating the accuracy with which the model must generate a comparison value. Note that this can be redefined in child classes.
 
virtual double visibility_amplitude (datum_visibility_amplitude &d, double accuracy)
 Returns visibility ampitudes in Jy computed from the image given a datum_visibility_amplitude object, containing all of the accoutrements. While this provides access to the actual data value, the two could be separated if necessary. Also takes an accuracy parameter with the same units as the data, indicating the accuracy with which the model must generate a comparison value. Note that this can be redefined in child classes.
 
virtual double closure_phase (datum_closure_phase &d, double accuracy)
 Returns closure phase in degrees computed from the image given a datum_closure_phase object, containing all of the accoutrements. While this provides access to the actual data value, the two could be separated if necessary. Also takes an accuracy parameter with the same units as the data, indicating the accuracy with which the model must generate a comparison value. Note that this can be redefined in child classes.
 
virtual double closure_amplitude (datum_closure_amplitude &d, double accuracy)
 Returns closure amplitude computed from the image given a datum_closure_phase object, containing all of the accoutrements. While this provides access to the actual data value, the two could be separated if necessary. Also takes an accuracy parameter with the same units as the data, indicating the accuracy with which the model must generate a comparison value. Note that this can be redefined in child classes.
 
void output_image (std::string fname, bool rotate=false)
 
void get_image (std::vector< std::vector< double > > &alpha, std::vector< std::vector< double > > &beta, std::vector< std::vector< double > > &I) const
 Provides direct access to the constructed image. Sets a 2D grid of angles (alpha, beta) in radians and intensities in Jy per steradian.
 
void get_visibilities (std::vector< std::vector< double > > &u, std::vector< std::vector< double > > &v, std::vector< std::vector< std::complex< double > > > &V) const
 Provides direct access to the complex visibilities. Sets a 2D grid of baselines (u,v) in lambda, and visibilites in Jy.
 
void get_visibility_amplitudes (std::vector< std::vector< double > > &u, std::vector< std::vector< double > > &v, std::vector< std::vector< double > > &V) const
 Provides direct access to the visibility amplitudes. Sets a 2D grid of baselines (u,v) in lambda, and visibilites in Jy.
 
void use_spline_interp (bool use_spline)
 Provides ability to use bicubic spline interpolator (true) instead of regular bicubic. Code defaults to false.
 
void write_model_tag_file (std::string tagfilename="model_image.tag") const
 Write a unique identifying tag for use with the ThemisPy plotting features. This calls the overloaded version with the outstream, which is the only function that need be rewritten in child classes.
 
virtual void write_model_tag_file (std::ofstream &tagout) const
 Write a unique identifying tag for use with the ThemisPy plotting features. For most child classes, the default implementation is suffcient. However, should that not be the case, this is the only function that need be rewritten in child classes.
 

Private Member Functions

virtual void generate_image (std::vector< double > parameters, std::vector< std::vector< double > > &I, std::vector< std::vector< double > > &alpha, std::vector< std::vector< double > > &beta)
 A user-supplied function that generates and returns rectalinear grid of intensities in Jy/str located at pixels centered on angular positions alpha and beta, both specified in radians and aligned with a fiducial direction. Note that the parameter vector has had the position removed. Note that it will be assumed that alpha and beta are defined as the image appears on the sky, e.g., beta running from S to N and alpha running from E to W.
 
void ensemble_blur_image ()
 < Blurs the input model image using the ensemble averaged blurring kernel. More...
 
void compute_kphase_screen (std::vector< double > screen_params)
 
void generate_model_visibilities ()
 
double Dphi (double r, double phi) const
 Power spectrum of phase screen.
 
double Q (double qx, double qy) const
 Generate random Gaussian epsilon screen from parameter list.
 
void make_epsilon_screen (std::vector< double > screen_params)
 
std::vector< std::vector< double > > ifft_2d (const std::vector< std::vector< std::complex< double > > > &V)
 
std::vector< std::vector< std::complex< double > > > ifft_shift (const std::vector< std::vector< std::complex< double > > > &V)
 

Private Attributes

model_image_model
 
std::vector< double > _current_model_params
 Current parameters for the intrinsic model.
 
const size_t _nModes
 Number of modes to include in scattering screen.
 
const double _tobs
 Observing time since start in seconds.
 
const double _frequency
 frequency of observation
 
const double _wavelength
 
const std::string _scattering_model
 Which scattering model to use.
 
const double _observer_screen_distance
 
const double _source_screen_distance
 
const double _theta_maj_mas_cm
 
const double _theta_min_mas_cm
 Semi-major and minor axis length of Gaussian component of scattering screen in mas at reference wavelength of 1 cm.
 
const double _phi0
 (90-POS_ANG)
 
const double _M
 Magnification factor.
 
const double _rF
 Fresnel scale.
 
const double _scatt_alpha
 Scattering turbulence and dissipation scale power law parameter. Kolmogorov scattering has alpha=5/3.
 
const double _rin
 Inner scale for the turbulence in the scattering screen.
 
const double _vx_ss_kms
 
const double _vy_ss_kms
 Velocity of scattering screen with respect to earth in km/s.
 
std::vector< std::vector< std::complex< double > > > _epsilon
 
size_t _nray
 
double _fov
 
double _zeta0
 
double _Qbar
 Power spectrum normalization.
 
double _C_scatt_0
 Scattering coefficent for phase structure function.
 
double _kZeta
 Scattering kernel phase structure coefficient for different wandering magnetic field models.
 
double _Bmaj
 
double _Bmin
 
std::valarray< double > _i2drxK
 
std::valarray< double > _i2dryK
 
std::valarray< double > _i2dDxPhi
 
std::valarray< double > _i2dDyPhi
 
Interpolator2D _i2D_Dxphi
 
Interpolator2D _i2D_Dyphi
 2D interpolator object for the gradient of the phase screen.
 
std::vector< std::vector< double > > _Iea
 2D grid of intensities of ensemble averaged image at pixel locations in Jy/str.
 
std::vector< std::vector< double > > _u
 2D grid of horizontal baseline locations in lambda, relative to the fiducial direction of the image (i.e., unrotated by the position angle).
 
std::vector< std::vector< double > > _v
 2D grid of vertical baseline locations in lambda, relative to the fiducial direction of the image (i.e., unrotated by the position angle).
 
std::vector< std::vector< std::complex< double > > > _Vsrc
 complex visibilities of intrinsic source More...
 
std::valarray< double > _i2drx
 
std::valarray< double > _i2dry
 
std::valarray< double > _i2dIea
 
Interpolator2D _i2D_Iea
 2D interpolator object to estimate intensity of ensemble average position at arbitrary alpha,beta coordinates. model_image uses bicubic interpolation. More...
 
_P_phi _P_phi_func
 Normalization for wandering magnetic field model.
 
double _P_phi_norm
 

Additional Inherited Members

- Protected Attributes inherited from model_image
MPI_Comm _comm
 
bool _generated_model
 True when a model is generated with generate_model.
 
bool _generated_visibilities
 True when model visibilities have been computed.
 
bool _use_spline
 True when want to use bicubic spline interpolator.
 
double _position_angle
 Position angle of the image, assumed to be the last element of the parameter list passed to generate_model. Assumed to be in radians and defined E of N.
 
std::vector< double > _current_parameters
 Last set of parameters passed to generate_image, useful to determine if it is necessary to recompute the model. Useful, e.g., if we are varying only position angle, or if recomputing for a number of different data sets at the same set of parameters.
 
std::vector< std::vector< double > > _alpha
 2D grid of horizonal pixel locations in radians, relative to the fiducial direction of the image (i.e., unrotated by the position angle).
 
std::vector< std::vector< double > > _beta
 2D grid of vertical pixel locations in radians, relative to the fiducial direction of the image (i.e., unrotated by the position angle).
 
std::vector< std::vector< double > > _I
 2D grid of intensities at pixel locations in Jy/str.
 
- Protected Attributes inherited from model_visibility
MPI_Comm _comm
 
- Protected Attributes inherited from model_visibility_amplitude
MPI_Comm _comm
 
- Protected Attributes inherited from model_closure_phase
MPI_Comm _comm
 
- Protected Attributes inherited from model_closure_amplitude
MPI_Comm _comm
 

Detailed Description

Defines the interface for models that generates refractive scattered movies (it is templated off the model image class that you will scatter)

Author
Paul Tiede
Date
Oct. 2018

Scattering implementation assumes we are in average strong-scattering regime where diffractive scintillation has been averaged over. Parameter list:

  • parameters[0...n-1] ... Model parameters for the source sans position angle, i.e. if a RIAF then spin inclination and so on.
  • parameters[n] ... position angle for the image.
  • parameters[n+1...m] ... nModes^2-1=m-n-1 normalized Fourier modes of the scattering screen where nModes is passed in the constructor. ... note the normalized means the distribution of these modes is given by a Gaussian with zero mean and unit variance.
Author
Paul Tiede
Date
Oct. 2018

Scattering implementation assumes we are in average strong-scattering regime where diffractive scintillation has been averaged over. Parameter list:

  • parameters[0...n-1] ... Model parameters for the source sans position angle, i.e. if a RIAF then spin inclination and so on.
  • parameters[n] ... position angle for the images.
  • parameters[n+1...m] ... nModes^2-1=m-n-1 normalized Fourier modes of the scattering screen where nModes is passed in the constructor. ... note the normalized means the distribution of these modes is given by a Gaussian with zero mean and unit variance.

Constructor & Destructor Documentation

model_image_refractive_scattering ( model_image model,
size_t  nModes,
double  tobs,
double  frequency = 230e9,
std::string  scattering_model = "dipole",
double  observer_screen_distance = 2.82*3.086e21,
double  source_screen_distance = 5.53*3.086e21,
double  theta_maj_mas_cm = 1.38,
double  theta_min_ma_cm = 0.703,
double  POS_ANG = 81.9,
double  scatt_alpha = 1.38,
double  r_in = 800e5,
double  r_out = 1e20,
double  vs_ss_kms = 50.0,
double  vy_ss_kms = 0.0 
)

Constructor for refractive scattering class. Depends on whether the source you want to scatter is a model image class or one of the eht observables model classes e.g. model_visibility_amplitude.

nModes defines the number of modes to use when constructing the scattering screen. The rest of the parameters define the ensemble average scattering screen, where the default values follow Johnson et al. 2018.

Here is the call graph for this function:

Member Function Documentation

void ensemble_blur_image ( )
private

Utility function. Generates the gradient of the phase screen since this is what we need.

Here is the call graph for this function:

Here is the caller graph for this function:

void set_image_resolution ( size_t  nray)

Sets the image scattered image resolution to be used resolution of the image is nrayxnray, the default is 128,128 which is probably too much

void set_mpi_communicator ( MPI_Comm  comm)
virtual

Defines a set of processors provided to the model for parallel computation via an MPI communicator. Only facilates code parallelization if the model computation is parallelized via MPI.

Reimplemented from model_image.

Here is the call graph for this function:

void set_screen_size ( double  fov)

Sets the fov size of the image in units of radians The current default is 100uas.

virtual size_t size ( ) const
inlinevirtual

Size of the supplied model image plus the number of modes to include in the screen. Gives the number of parameters to fit.

Reimplemented from model_image.

Here is the call graph for this function:

Member Data Documentation

Interpolator2D _i2D_Iea
private

Dphi function for the ensemble averaged image.

std::vector<std::vector<std::complex<double> > > _Vsrc
private

Arrays for the interpolator objects


The documentation for this class was generated from the following files: