woden_struct_defs¶
API documentation for woden_struct_defs.h.
Typedefs
-
typedef struct _components_t components_t¶
A struct to contain COMPONENT information for either multiple POINT, GAUSSIAN, or SHAPELET COMPONENTs. Three of these are included in
source_tto specify all COMPONENTs for a single source. This allows common values like RA,Dec,Flux etc information to be stored the same way for different COMPONENT types
-
typedef struct _source_catalogue_t source_catalogue_t¶
A struct to contain multiple
source_ttype sky models andbeam_settings_tprimary beam settings, to be iterated over bycalculate_visibilities
-
typedef struct _beam_settings_t beam_settings_t¶
A struct to contain settings pertaining to the primary beam
-
typedef struct _visibility_set_t visibility_set_t¶
Struct to contain simulation parameters and visibility outputs
-
typedef struct _woden_settings_t woden_settings_t¶
Struct to contain user defined settings for simulation
-
typedef struct _array_layout_t array_layout_t¶
Struct to contain array layout values. Here, a single receiving element is sometimes called an antenna, sometimes called a ‘tile’ (MWA lingo). This is equivalent to a ‘station’ in SKA_LOW talk.
Enums
-
enum e_component_type¶
Values:
-
enumerator POINT¶
Point source type component
-
enumerator GAUSSIAN¶
Gaussian type component
-
enumerator SHAPELET¶
Shapelet type component
-
enumerator POINT¶
-
enum e_beamtype¶
Values:
-
enumerator NO_BEAM¶
Do not use a primary beam in the simulation
-
enumerator GAUSS_BEAM¶
Use a analytic Gaussian primary beam
-
enumerator FEE_BEAM¶
Use the RTS MWA FEE primary beam code
-
enumerator ANALY_DIPOLE¶
Use an analytic MWA dipole primary beam
-
enumerator FEE_BEAM_INTERP¶
Use the RTS MWA FEE primary beam at all frequencies. Should be using an hdf5 file that has been frequency interpolated
-
enumerator MWA_ANALY¶
Use an analytic MWA tile primary beam
-
enumerator NO_BEAM¶
-
struct _components_t¶
- #include <woden_struct_defs.h>
A struct to contain COMPONENT information for either multiple POINT, GAUSSIAN, or SHAPELET COMPONENTs. Three of these are included in
source_tto specify all COMPONENTs for a single source. This allows common values like RA,Dec,Flux etc information to be stored the same way for different COMPONENT typesPublic Members
-
double *ras¶
COMPONENT right ascensions (radians)
-
double *decs¶
COMPONENT declinations (radians)
-
double *ref_freqs¶
COMPONENT Flux density reference frequencies (Hz)
-
user_precision_t *ref_stokesI¶
COMPONENT Stokes I reference flux density (Jy)
-
user_precision_t *ref_stokesQ¶
COMPONENT Stokes Q reference flux density (Jy)
-
user_precision_t *ref_stokesU¶
COMPONENT Stokes U reference flux density (Jy)
-
user_precision_t *ref_stokesV¶
COMPONENT Stokes V reference flux density (Jy)
-
user_precision_t *SIs¶
COMPONENT spectral indexes
-
user_precision_t *shape_coeffs¶
Scaling coefficients for SHAPELET basis functions
-
user_precision_t *n1s¶
1st basis function order for SHAPELET basis functions
-
user_precision_t *n2s¶
2nd basis function order for SHAPELET basis functions
-
user_precision_t *majors¶
GAUSSIAN/SHAPELET major axis (beta1, radians)
-
user_precision_t *minors¶
GAUSSIAN/SHAPELET minor axis (beta2, radians)
-
user_precision_t *pas¶
GAUSSIAN/SHAPELET position angles (radians)
-
user_precision_t *param_indexes¶
An index value to match each coeff, n1, and n2 to the correct ra, dec, major, minor, pa for a SHAPELET
-
user_precision_t *azs¶
SHAPELET source azimuth angles for all time steps
-
user_precision_t *zas¶
SHAPELET source zenith angles for all time steps
-
double *beam_has¶
Hour angle of COMPONENTs for all time steps, used for beam calculations
-
double *beam_decs¶
Declinations of COMPONENTs for all time steps, used for beam calculations
-
int num_primarybeam_values¶
Number of beam calculations needed for COMPONENTs
-
user_precision_complex_t **gxs¶
North-South Beam gain values for all directions, frequencies, and times for these COMPONENTS
-
user_precision_complex_t **Dxs¶
North-South Beam leakage values for all directions, frequencies, and times for these COMPONENTS
-
user_precision_complex_t **Dys¶
East-West Beam leakage values for all directions, frequencies, and times for these COMPONENTS
-
user_precision_complex_t **gys¶
East-West Beam gain values for all directions, frequencies, and times for these COMPONENTS
-
double *ls¶
Device memory l cosine direction coords for these COMPONENTs
-
double *ms¶
Device memory m cosine direction coords for these COMPONENTs
-
double *ns¶
Device memory n cosine direction coords for these COMPONENTs
-
int HELP¶
-
double *ras¶
-
struct _source_t¶
- #include <woden_struct_defs.h>
A struct to contain sky model values for a single SOURCE
Public Members
-
char name[32]¶
Source name
-
int n_comps¶
Total number of COMPONENTs in source
-
int n_points¶
Number of POINT source COMPONENTs
-
int n_gauss¶
Number of GAUSSIAN source COMPONENTs
-
int n_shapes¶
Number of SHAPELET source COMPONENTs
-
int n_shape_coeffs¶
Total number of SHAPELET coefficients
-
components_t point_components¶
components_tholding component information for all POINT COMPONENTs in this SOURCE.
-
components_t gauss_components¶
components_tholding component information for all GAUSSIAN COMPONENTs in this SOURCE.
-
components_t shape_components¶
components_tholding component information for all SHAPELET COMPONENTs in this SOURCE.
-
components_t d_point_components¶
components_tholding component information for all POINT COMPONENTs in this SOURCE.
-
components_t d_gauss_components¶
components_tholding component information for all GAUSSIAN COMPONENTs in this SOURCE.
-
components_t d_shape_components¶
components_tholding component information for all SHAPELET COMPONENTs in this SOURCE.
-
char name[32]¶
-
struct _source_catalogue_t¶
- #include <woden_struct_defs.h>
A struct to contain multiple
source_ttype sky models andbeam_settings_tprimary beam settings, to be iterated over bycalculate_visibilities
-
struct _beam_settings_t¶
- #include <woden_struct_defs.h>
A struct to contain settings pertaining to the primary beam
Public Members
-
user_precision_t gauss_sdec¶
Sine of the declination of the pointing for a Gaussian primary beam
-
user_precision_t gauss_cdec¶
Cosine of the declination of the pointing for a Gaussian primary beam
-
double gauss_ha¶
Hour angle of the pointing for a Gaussian primary beam
-
user_precision_t beam_FWHM_rad¶
FWHM of requested Gaussian primary beam, at reference frequnecy
-
double beam_ref_freq¶
Reference frequency for the given FWHM of Gaussian primary beam
-
e_beamtype beamtype¶
What type of primary beam to simulate - see
e_beamtype
-
user_precision_t *MWAFEE_freqs¶
The frequencies of the initialised MWAFEE beams in FEE_beams
-
int num_MWAFEE¶
Number of MWAFEE beam instances to cover all desired frequencies
-
struct FEEBeamCUDA *cuda_fee_beam¶
Single initialised hyperbeam device model for desired pointing
-
struct FEEBeam *fee_beam¶
Single initialised hyperbeam host model for desired pointing
-
char hyper_error_str[100]¶
-
double base_middle_freq¶
-
uint32_t *hyper_delays¶
-
user_precision_t gauss_sdec¶
-
struct _visibility_set_t¶
- #include <woden_struct_defs.h>
Struct to contain simulation parameters and visibility outputs
Public Members
-
user_precision_t *us_metres¶
Output \(u\) for all time steps, frequency steps, and baselines
-
user_precision_t *vs_metres¶
Output \(v\) for all time steps, frequency steps, and baselines
-
user_precision_t *ws_metres¶
Output \(w\) for all time steps, frequency steps, and baselines
-
double *allsteps_sha0s¶
Sine of hour angle of phase centre for all time steps, frequency steps, and baselines
-
double *allsteps_cha0s¶
Cosine of hour angle of phase centre for all time steps, frequency steps, and baselines
-
double *allsteps_lsts¶
Local sidereal time for all time steps, frequency steps, and baselines (radians)
-
user_precision_t *allsteps_wavelengths¶
Wavelengths for all time steps, frequency steps, and baselines (metres)
-
double *channel_frequencies¶
Frequencies for a frequency steps (Hz)
-
user_precision_t *sum_visi_XX_real¶
Real values for XX polarisation for all time steps, frequency steps, and baselines
-
user_precision_t *sum_visi_XX_imag¶
Imaginary values for XX polarisation for all time steps, frequency steps, and baselines
-
user_precision_t *sum_visi_XY_real¶
Real values for XY polarisation for all time steps, frequency steps, and baselines
-
user_precision_t *sum_visi_XY_imag¶
Imaginary values for XY polarisation for all time steps, frequency steps, and baselines
-
user_precision_t *sum_visi_YX_real¶
Real values for YX polarisation for all time steps, frequency steps, and baselines
-
user_precision_t *sum_visi_YX_imag¶
Imaginary values for YX polarisation for all time steps, frequency steps, and baselines
-
user_precision_t *sum_visi_YY_real¶
Real values for YY polarisation for all time steps, frequency steps, and baselines
-
user_precision_t *sum_visi_YY_imag¶
Imaginary values for YY polarisation for all time steps, frequency steps, and baselines
-
user_precision_t *us_metres¶
-
struct _woden_settings_t¶
- #include <woden_struct_defs.h>
Struct to contain user defined settings for simulation
Public Members
-
double lst_base¶
Local sidereal time for first time step (radians)
-
double ra0¶
Right ascension of phase centre (radians)
-
double dec0¶
Declination of phase centre (radians)
-
double sdec0¶
Sine of Declination of phase centre (radians)
-
double cdec0¶
Cosine of Declination of phase centre (radians)
-
int num_baselines¶
Number of baselines this array layout has
-
int num_freqs¶
Number of frequencies per coarse band
-
double frequency_resolution¶
Frequency resolution of a fine channel (Hz)
-
double base_low_freq¶
The lowest fine channel frequency of band 1
-
int num_time_steps¶
Number of time steps to simulate
-
double time_res¶
Time resolution of simulation (seconds)
-
const char *cat_filename¶
Path to WODEN-style sky model
-
int num_bands¶
Number of coarse frequency bands to simulate
-
int *band_nums¶
Which number coarse bands to simulate (e.g 1,4,6)
-
int sky_crop_type¶
Whether to crop sky models by SOURCE or COMPONENT
-
e_beamtype beamtype¶
What type of primary beam to simulate with
-
user_precision_t gauss_beam_FWHM¶
FWHM of Gaussian primary beam (degrees)
-
double gauss_beam_ref_freq¶
Reference frequency for given Gaussian primary beam FWHM
-
long int chunking_size¶
Maximum number of COMPONENTs to include in a single chunk
-
const char *hdf5_beam_path¶
Path to *.hf file containing MWA FEE beam spherical harmonic information
-
double jd_date¶
Julian date at beginning of simulation
-
bool array_layout_file¶
Do we have a path to the array layout or not
-
const char *array_layout_file_path¶
Path to file containing E,N,H coords of array layout
-
double latitude¶
Latitude of the array to simulate (radians)
-
user_precision_t longitude¶
Longitude of the array to simulate (radians)
-
user_precision_t FEE_ideal_delays[16]¶
Delay values specifying the pointing for the MWA FEE beam model
-
user_precision_t coarse_band_width¶
Frequency bandwidth of a single coarse band (Hz)
-
double gauss_ra_point¶
The initial Right Ascension to point the Gaussian beam at (radians)
-
double gauss_dec_point¶
The initial Declination to point the Gaussian beam at (radians)
-
int num_visis¶
Total number of visiblities to simulate, so freqs*times*baselines
-
double base_band_freq¶
The lowest fine channel frequency in the current band being simulated
-
int do_precession¶
Boolean of whether to apply precession to the array layout or not
-
double lst_base¶
-
struct _array_layout_t¶
- #include <woden_struct_defs.h>
Struct to contain array layout values. Here, a single receiving element is sometimes called an antenna, sometimes called a ‘tile’ (MWA lingo). This is equivalent to a ‘station’ in SKA_LOW talk.
Public Members
-
double *ant_X¶
Local \(X\) location of all antenna/tiles
-
double *ant_Y¶
Local \(Y\) location of all antenna/tiles
-
double *ant_Z¶
Local \(Z\) location of all antenna/tiles
-
double *X_diff_metres¶
The length of all baselines in \(X\) (metres)
-
double *Y_diff_metres¶
The length of all baselines in \(Y\) (metres)
-
double *Z_diff_metres¶
The length of all baselines in \(Z\) (metres)
-
double *ant_east¶
Local east location of all antenna/tiles
-
double *ant_north¶
Local north location of all antenna/tiles
-
double *ant_height¶
Local height location of all antenna/tiles
-
double latitude¶
Latitude of the array (radians)
-
int num_baselines¶
Number of baselines in the array
-
int num_tiles¶
Number of antenna/tiles in the array
-
double lst_base¶
Local sidereal time of the first time step (radians)
-
double *ant_X¶