source_components¶
API documentation for source_components.c
Device methods to extrapolate flux densities, assign and apply primary beam gains, and visibility responses for different sky model COMPONENT types.
Typedefs
-
typedef struct _d_beam_gains_t d_beam_gains_t¶
A struct to contain primary beam values for a give COMPONENT
Functions
-
__device__ void extrap_stokes(user_precision_t *d_allsteps_wavelengths, double *d_ref_freqs, user_precision_t *d_ref_stokesI, user_precision_t *d_ref_stokesQ, user_precision_t *d_ref_stokesU, user_precision_t *d_ref_stokesV, user_precision_t *d_SIs, int iComponent, int iBaseline, user_precision_t *flux_I, user_precision_t *flux_Q, user_precision_t *flux_U, user_precision_t *flux_V)¶
Assuming a simple power-law SED of \( S \propto \nu^\alpha \), extrapolate Stokes flux density parameters at the requested allsteps_wavelengths.
The 4 arrays of reference flux densities (
d_ref_stokesI,d_ref_stokesQ,d_ref_stokesU,d_ref_stokesV), reference frequenciesd_ref_freqs, and spectral indexesd_SIsshould all be the same length.iComponentindexes these 6 arrays.iBaselineindexes the arrayd_allsteps_wavelengths, containging the wavelengths to extrapolate the reference fluxes to. Returns the flux density values extrapolated to the selected wavelength and reference data (flux_I,flux_Q,flux_U,flux_V).- Parameters
d_allsteps_wavelengths – [in] Wavelengths array to index with iBaseline (metres)
d_ref_freqs – [in] Component reference frequencies to index with iBaseline (Hz)
d_ref_stokesI – [in] Stokes I reference flux densities to reference with iComponent (Jy)
d_ref_stokesQ – [in] Stokes Q reference flux densities to reference with iComponent (Jy)
d_ref_stokesU – [in] Stokes U reference flux densities to reference with iComponent (Jy)
d_ref_stokesV – [in] Stokes V reference flux densities to reference with iComponent (Jy)
d_SIs – [in] Spectral indicies to index with iBaseline
iComponent – [in] Component index
iBaseline – [in] Baseline index
*flux_I – [inout] Pointer to extrapolated Stokes I (Jy)
*flux_Q – [inout] Pointer to extrapolated Stokes Q (Jy)
*flux_U – [inout] Pointer to extrapolated Stokes U (Jy)
*flux_V – [inout] Pointer to extrapolated Stokes V (Jy)
-
__device__ cuUserComplex calc_measurement_equation(user_precision_t *d_us, user_precision_t *d_vs, user_precision_t *d_ws, double *d_ls, double *d_ms, double *d_ns, const int iBaseline, const int iComponent)¶
Calculate the complex exponential phase part of the visibility function from the given u,v,w and l,m,n coord arrays.
Calculates the following complex exponential:
\[ \exp \left( 2\pi i\left( ul + vm + w(n-1) \right) \right) \]where
u = d_us[iBaseline] v = d_vs[iBaseline] w = d_ws[iBaseline] l = d_ns[iComponent] m = d_ms[iComponent] n = d_ns[iComponent]
Note there is no negative sign here - testing with the current implementation through imaging software, and from comparisons to the OSKAR and RTS simulators this yields the correct output visibilities.
- Parameters
*d_us – [in] Pointer to \(u\) coodinates (wavelengths)
*d_vs – [in] Pointer to \(v\) coodinates (wavelengths)
*d_ws – [in] Pointer to \(w\) coodinates (wavelengths)
*d_ls – [in] Pointer to \(l\) coodinates
*d_ms – [in] Pointer to \(m\) coodinates
*d_ns – [in] Pointer to \(n\) coodinates
iBaseline – [in] Index for \(u,v,w\)
iComponent – [in] Index for \(l,m,n\)
- Returns
visi, acuUserComplexof the visibility
-
__device__ void apply_beam_gains(cuUserComplex g1x, cuUserComplex D1x, cuUserComplex D1y, cuUserComplex g1y, cuUserComplex g2x, cuUserComplex D2x, cuUserComplex D2y, cuUserComplex g2y, user_precision_t flux_I, user_precision_t flux_Q, user_precision_t flux_U, user_precision_t flux_V, cuUserComplex visi_component, cuUserComplex *visi_XX, cuUserComplex *visi_XY, cuUserComplex *visi_YX, cuUserComplex *visi_YY)¶
Given primary beam gains and leakage terms for antenna 1
g1x, D1x, D1y, gyand antenna 2g1x, D1x, D1y, gy,the complex visibility phase across those two antennasvisi, and the Stokes parameters of a source, simulate the observed XX,XY,YX,YY instrumental cross-correlated visibilities, where ‘x’ means aligned north-south, ‘y’ means east-west.Performs the following calculations:
\[\begin{split}\begin{eqnarray*} \mathrm{V}^{XX}_{12} = (g_{1x}g_{2x}^{\ast} + D_{1x}D_{2x}^{\ast})\mathrm{V}^{I}_{12} + (g_{1x}g_{2x}^{\ast} - D_{1x}D_{2x}^{\ast})\mathrm{V}^{Q}_{12} \\ + (g_{1x}D_{2x}^{\ast} + D_{1x}g_{2x}^{\ast})\mathrm{V}^{U}_{12} + i(g_{1x}D_{2x}^{\ast} - D_{1x}g_{2x}^{\ast})\mathrm{V}^{V}_{12} \end{eqnarray*}\end{split}\]\[\begin{split}\begin{eqnarray*} \mathrm{V}^{XY}_{12} = (g_{1x}D_{2y}^{\ast} + D_{1x}g_{2y}^{\ast})\mathrm{V}^{I}_{12} + (g_{1x}D_{2y}^{\ast} - D_{1x}g_{2y}^{\ast})\mathrm{V}^{Q}_{12} \\ + (g_{1x}g_{2y}^{\ast} + D_{1x}D_{2y}^{\ast})\mathrm{V}^{U}_{12} + i(g_{1x}g_{2y}^{\ast} - D_{1x}D_{2y}^{\ast})\mathrm{V}^{V}_{12} \end{eqnarray*}\end{split}\]\[\begin{split}\begin{eqnarray*} \mathrm{V}^{YX}_{12} = (D_{1y}g_{2x}^{\ast} + g_{1y}D_{2x}^{\ast})\mathrm{V}^{I}_{12} + (D_{1y}g_{2x}^{\ast} - g_{1y}D_{2x}^{\ast})\mathrm{V}^{Q}_{12} \\ + (D_{1y}D_{2x}^{\ast} + g_{1y}g_{2x}^{\ast})\mathrm{V}^{U}_{12} + i(D_{1y}D_{2x}^{\ast} - g_{1y}g_{2x}^{\ast})\mathrm{V}^{V}_{12} \end{eqnarray*}\end{split}\]\[\begin{split}\begin{eqnarray*} \mathrm{V}^{YY}_{12} = (D_{1y}D_{2y}^{\ast} + g_{1y}g_{2y}^{\ast})\mathrm{V}^{I}_{12} + (D_{1y}D_{2y}^{\ast} - g_{1y}g_{2y}^{\ast})\mathrm{V}^{Q}_{12} \\ + (D_{1y}g_{2y}^{\ast} + g_{1y}D_{2y}^{\ast})\mathrm{V}^{U}_{12} + i(D_{1y}g_{2y}^{\ast} - g_{1y}D_{2y}^{\ast})\mathrm{V}^{V}_{12} \end{eqnarray*}\end{split}\]where \({\ast}\) means complex conjugate, and
\[\begin{split}\begin{eqnarray*} \mathrm{V}_{12} &=& \mathrm{V}_{\mathrm{env}}\exp \left( 2\pi i\left( u_{12}l + v_{12}m + w_{12}(n-1) \right) \right) \\ \mathrm{V}^{I}_{12} &=& I(l,m) \mathrm{V}_{12} \\ \mathrm{V}^{Q}_{12} &=& Q(l,m) \mathrm{V}_{12} \\ \mathrm{V}^{U}_{12} &=& U(l,m) \mathrm{V}_{12} \\ \mathrm{V}^{V}_{12} &=& V(l,m) \mathrm{V}_{12} \end{eqnarray*}\end{split}\]where \( \mathrm{V}_{\mathrm{env}} \) is an visibility envelope that turns a component into a GAUSSIAN or SHAPELET component, and has been applied in
visi_component.- Parameters
g1x – [in] Beam gain antenna 1 in north-south
D1x – [in] Beam leakage antenna 1 from north-south
D1y – [in] Beam gain antenna 1 in east-west
g1y – [in] Beam leakage antenna 1 from east-west
g2x – [in] Beam gain antenna 2 in north-south
D2x – [in] Beam leakage antenna 2 from north-south
D2y – [in] Beam gain antenna 2 in east-west
g2y – [in] Beam leakage antenna 2 from east-west
flux_I – [in] Stokes I flux density (Jy)
flux_Q – [in] Stokes Q flux density (Jy)
flux_U – [in] Stokes U flux density (Jy)
flux_V – [in] Stokes V flux density (Jy)
visi_component – [in] Complex visibility across antennas 1 and 2
visi_XX – [inout] Output XX instrumental visibility
visi_XY – [inout] Output XY instrumental visibility
visi_YX – [inout] Output YX instrumental visibility
visi_YY – [inout] Output YY instrumental visibility
-
__device__ void get_beam_gains(int iBaseline, int iComponent, int num_freqs, int num_baselines, int num_components, int num_times, int beamtype, cuUserComplex *d_primay_beam_J00, cuUserComplex *d_primay_beam_J01, cuUserComplex *d_primay_beam_J10, cuUserComplex *d_primay_beam_J11, cuUserComplex *g1x, cuUserComplex *D1x, cuUserComplex *D1y, cuUserComplex *g1y, cuUserComplex *g2x, cuUserComplex *D2x, cuUserComplex *D2y, cuUserComplex *g2y)¶
Given the type of primary beam simulated
beamtype, select the beam gains and leakages from the arraysd_primay_beam_J*that match the indexesiBaseline,iComponent. NOTE that currently, primary beams are assumed to be the same for all recieving elements, and so this function only takes in one set of primary beam arrays.This function is built to return the correct beam gain for a given component on the sky, at a given time, at a given frequency, for a given baseline. The 4 arrays
d_primay_beam_J00,d_primay_beam_J01,d_primay_beam_J10,d_primay_beam_J11should contain the primary beam settings for all times, all frequencies, and COPMONENTs on the sky, and so should benum_freqs*num_components*num_timeslong. The order elements should increment through COMPONENT (fastest changing), freqeuncy, and time (slowest changing).iBaselineis a combined index of baseline, frequency, and time.If
beamtype == NO_BEAM, setg1x = g1y = g2x = g2y = 1andD1x = D1y = D2x = D2y = 0- Todo:
Prepare this function for two different primary beam responses
- Parameters
iBaseline – [in] Index of which baseline, freq, and time we are on
iComponent – [in] COMPONENT index
num_freqs – [in] Number of frequencies in simulation
num_baselines – [in] Number of baselines for one time and one frequency step
num_components – [in] Number of COMPONENTs
num_times – [in] Number of times in simulation
beamtype – [in] Beam type see
woden_struct_defs.e_beamtype*d_primay_beam_J00 – [in] Pointer towards array of primary beam J[0,0] (north-south gain)
*d_primay_beam_J01 – [in] Pointer towards array of primary beam J[0,1] (north-south leakage)
*d_primay_beam_J10 – [in] Pointer towards array of primary beam J[1,0] (east-west gain)
*d_primay_beam_J11 – [in] Pointer towards array of primary beam J[1,1] (east-west leakage)
*g1x – [inout] Beam gain antenna 1 in north-south
*D1x – [inout] Beam leakage antenna 1 from north-south
*D1y – [inout] Beam gain antenna 1 in east-west
*g1y – [inout] Beam leakage antenna 1 from east-west
*g2x – [inout] Beam gain antenna 2 in north-south
*D2x – [inout] Beam leakage antenna 2 from north-south
*D2y – [inout] Beam gain antenna 2 in east-west
*g2y – [inout] Beam leakage antenna 2 from east-west
-
__device__ void update_sum_visis(int iBaseline, int iComponent, int num_freqs, int num_baselines, int num_components, int num_times, int beamtype, cuUserComplex *d_primay_beam_J00, cuUserComplex *d_primay_beam_J01, cuUserComplex *d_primay_beam_J10, cuUserComplex *d_primay_beam_J11, cuUserComplex visi_component, user_precision_t flux_I, user_precision_t flux_Q, user_precision_t flux_U, user_precision_t flux_V, user_precision_t *d_sum_visi_XX_real, user_precision_t *d_sum_visi_XX_imag, user_precision_t *d_sum_visi_XY_real, user_precision_t *d_sum_visi_XY_imag, user_precision_t *d_sum_visi_YX_real, user_precision_t *d_sum_visi_YX_imag, user_precision_t *d_sum_visi_YY_real, user_precision_t *d_sum_visi_YY_imag)¶
Given the visibility between two recieving elements for a COMPONENT
visi_component, at the baseline/freq/time index given byiBaselineand COMPONENT indexiComponent, apply the instrumental primary beam and COMPONENT Stokes parameters to create instrumental XX,XY,YX,YY visibilities, and sum them into real and imaginary XX,XY,YX,YY visibilities arraysd_sim_visi_*_realandd_sim_visi_*_imag.Uses
get_beam_gainsandapply_beam_gainsas described above to apply the gains - see descriptions for what should be the arguments to them.- Parameters
iBaseline – [in] Index of which baseline, freq, and time we are on
iComponent – [in] COMPONENT index
num_freqs – [in] Number of frequencies in simulation
num_baselines – [in] Number of baselines for one time and one frequency step
num_components – [in] Number of COMPONENTs
num_times – [in] Number of times in simulation
beamtype – [in] Beam type see
woden_struct_defs.e_beamtype*d_primay_beam_J00 – [in] Pointer towards array of primary beam J[0,0] (north-south gain)
*d_primay_beam_J01 – [in] Pointer towards array of primary beam J[0,1] (north-south leakage)
*d_primay_beam_J10 – [in] Pointer towards array of primary beam J[1,0] (east-west gain)
*d_primay_beam_J11 – [in] Pointer towards array of primary beam J[1,1] (east-west leakage)
visi_component – [in] Complex visibility across antennas 1 and 2
flux_I – [in] Stokes I flux density (Jy)
flux_Q – [in] Stokes Q flux density (Jy)
flux_U – [in] Stokes U flux density (Jy)
flux_V – [in] Stokes V flux density (Jy)
*d_sum_visi_XX_real – [inout] Pointer to array to sum real XX visibility into
*d_sum_visi_XX_imag – [inout] Pointer to array to sum imaginary XX visibility into
*d_sum_visi_XY_real – [inout] Pointer to array to sum real XY visibility into
*d_sum_visi_XY_imag – [inout] Pointer to array to sum imaginary XY visibility into
*d_sum_visi_YX_real – [inout] Pointer to array to sum real YX visibility into
*d_sum_visi_YX_imag – [inout] Pointer to array to sum imaginary YX visibility into
*d_sum_visi_YY_real – [inout] Pointer to array to sum real YY visibility into
*d_sum_visi_YY_imag – [inout] Pointer to array to sum imaginary YY visibility into
-
void source_component_common(int num_components, int num_shape_coeffs, components_t *components, double *d_freqs, woden_settings_t *woden_settings, beam_settings_t *beam_settings, e_component_type comptype, components_t *d_components, d_beam_gains_t *d_component_beam_gains)¶
Performs necessary calculations that are common to all POINT, GAUSSIAN, and SHAPELET component types, calculating the \(l,m,n\) coordinates and primary beam values.
componentsshould be a populatedcomponents_tstruct, with appropriate fields filled depending on whether the COMPONENTS are POINT, GAUSSIAN, or SHAPELET (set bycomptype).This function uses
fundamental_coords::kern_calc_lmnto calculate \(l,m,n\), and stores the ouputs in thecomponents_tstructd_components. Uses thed_beam_gains_tstructd_component_beam_gainsto store the calculated primary beam responses. Each gain and leakage will havenum_components*woden_settings->num_time_steps*woden_settings->num_freqselements of memory allocated in the process.If the primary beam requested (set via beam_settings->beamtype ) is either GAUSS_BEAM or MWA_ANALY, both
components->beam_hasandcomponents->beam_decsneed to be set.Device arrays
d_primay_beam_J*are used to store primary beam outputs, and should have \(l,m,n\) coordinate arraysd_ls,d_ms,d_nsshould havenum_componentselements.sin_para_angsandcos_para_angsare sine and cosine of the parallactic angle, and only needed ifbeam_settings.beamtype == FEE_BEAM. They should benum_components*woden_settings->num_time_stepslong.beam_has,beam_decsare hour angles and declinations of components at all time steps, and only used whenbeam_settings.beamtype == GAUSS_BEAM. They should be also benum_components*woden_settings->num_time_steps.- Parameters
num_components – [in] Number of copmponents stored in
componentsnum_shape_coeffs – [in] Total number of shapelet coefficients in
componentscomponents – [in] A populated
components_tstruct as filled bycreate_sky_model::crop_sky_modelorchunk_sky_model::create_chunked_sky_modelsd_freqs – [in] Frequencies to calculate beam responses to (stored on the device)
woden_settings – [in] Populated
woden_settings_tstructbeam_settings – [in] Populated
beam_settings_tstructcomptype – [in] Either POINT, GAUSSIAN, SHAPELET
d_components – [inout] Pointer to
components_tstruct in which to malloc and store results on the device ind_component_beam_gains – [inout] Pointer to
d_beam_gains_tstruct in which to malloc and store results on the device in
-
__global__ void kern_calc_visi_point_or_gauss(components_t d_components, d_beam_gains_t d_component_beam_gains, user_precision_t *d_us, user_precision_t *d_vs, user_precision_t *d_ws, user_precision_t *d_sum_visi_XX_real, user_precision_t *d_sum_visi_XX_imag, user_precision_t *d_sum_visi_XY_real, user_precision_t *d_sum_visi_XY_imag, user_precision_t *d_sum_visi_YX_real, user_precision_t *d_sum_visi_YX_imag, user_precision_t *d_sum_visi_YY_real, user_precision_t *d_sum_visi_YY_imag, user_precision_t *d_allsteps_wavelengths, int num_components, int num_baselines, int num_freqs, int num_visis, int num_times, e_beamtype beamtype, e_component_type comptype)¶
Kernel to calculate the visibility response to a number
num_componentsof either POINT or GAUSSIAN COMPONENTs, and sum the outputs tod_sum_visi_*_real,d_sum_visi_*_imag.Uses the functions
extrap_stokes,calc_measurement_equation,update_sum_visisas detailed above to calculate the visibilities. Sets off a thread for each visibility to be calculated, with each thread looping over all COMPONENTs. This seems to keep the memory access low enough to be faster than having a second dimension over COMPONENT. Specify if you want POINT or GAUSSIAN usingcomptype.The code is almost exactly the same for POINT or GAUSSIANs, except when running with a Gaussian, a visiblity envelope is calculated as
\[\begin{split}\begin{eqnarray*} \mathrm{V}_{\mathrm{env}}&=& \exp\left( -\frac{\pi^2}{4\ln(2)} \left( k_x^2\theta_{\mathrm{maj}}^2 + k_y^2\theta_{\mathrm{min}}^2\right) \right) \\ k_x &=& \cos(\phi_{PAj})v + \sin(\phi_{PAj})u \\ k_y &=& -\sin(\phi_{PAj})v + \cos(\phi_{PAj})u \end{eqnarray*}\end{split}\]where \( \theta_{\mathrm{maj}}, \theta_{\mathrm{min}}, \phi_{PA} \) are the major axis, minor axis, and position angle. The visiblity equation is multipled by this envelope to modify from a POINT source into a GAUSSIAN
When called with
dim3 grid, threads, kernel should be called withgrid.xdefined, where:grid.x * threads.x >=
num_visi
- Parameters
d_components – [in] d_components Pointer to a populated
components_tstruct as filled bysource_components::source_component_commond_component_beam_gains – [in] Pointer to a populated
d_beam_gains_tstruct as filled bysource_components::source_component_common*d_us – [in] Visibility coord \(u\) for every baseline, frequency, and time step in the simulation (wavelengths)
*d_vs – [in] Visibility coord \(v\) for every baseline, frequency, and time step in the simulation (wavelengths)
*d_ws – [in] Visibility coord \(w\) for every baseline, frequency, and time step in the simulation (wavelengths)
*d_sum_visi_XX_real – [inout] Pointer to array to sum real XX visibility into
*d_sum_visi_XX_imag – [inout] Pointer to array to sum imaginary XX visibility into
*d_sum_visi_XY_real – [inout] Pointer to array to sum real XY visibility into
*d_sum_visi_XY_imag – [inout] Pointer to array to sum imaginary XY visibility into
*d_sum_visi_YX_real – [inout] Pointer to array to sum real YX visibility into
*d_sum_visi_YX_imag – [inout] Pointer to array to sum imaginary YX visibility into
*d_sum_visi_YY_real – [inout] Pointer to array to sum real YY visibility into
*d_sum_visi_YY_imag – [inout] Pointer to array to sum imaginary YY visibility into
*d_allsteps_wavelengths – [in] Wavelength for every baseline, frequency, and time step in the simulation
num_components – [in] Either number of POINT of GAUSSIANS components
num_baselines – [in] Number of baselines for one time, one frequency step
num_freqs – [in] Number of frequencies in simulation
num_visis – [in] Overall number of visibilities (baselines*freqs*times) in simulation
num_times – [in] Number of time steps in simulation
beamtype – [in] Beam type see
woden_struct_defs.e_beamtypecomptype – [in] Component type, either POINT or GAUSSIAN
-
__global__ void kern_calc_visi_shapelets(components_t d_components, d_beam_gains_t d_component_beam_gains, user_precision_t *d_us, user_precision_t *d_vs, user_precision_t *d_ws, user_precision_t *d_allsteps_wavelengths, user_precision_t *d_u_shapes, user_precision_t *d_v_shapes, user_precision_t *d_w_shapes, user_precision_t *d_sum_visi_XX_real, user_precision_t *d_sum_visi_XX_imag, user_precision_t *d_sum_visi_XY_real, user_precision_t *d_sum_visi_XY_imag, user_precision_t *d_sum_visi_YX_real, user_precision_t *d_sum_visi_YX_imag, user_precision_t *d_sum_visi_YY_real, user_precision_t *d_sum_visi_YY_imag, user_precision_t *d_sbf, int num_shapes, int num_baselines, int num_freqs, int num_visis, const int num_coeffs, int num_times, e_beamtype beamtype)¶
Kernel to calculate the visibility response to a number
num_shapesof SHAPELET COMPONENTs, and sum the outputs tod_sum_visi_*_real,d_sum_visi_*_imag.Uses the functions
extrap_stokes,calc_measurement_equation,update_sum_visisas detailed above to calculate the visibilities. Furthermore calculates the visibility envelope \(\mathrm{V}_{\mathrm{env}}\) to convert the basic visibility into a SHAPELET:\[\begin{split}\begin{eqnarray*} \mathrm{V}_{\mathrm{env}} &=& \sum^{n_k +n_l < n_\mathrm{max}}_{k,l} C_{n_k,n_l} B_{n_k,n_l}(k_x,k_y)\\ k_x &=& \frac{\pi}{\sqrt{2\ln(2)}} \left[\cos(\phi_{PA})v_{\mathrm{comp}} + \sin(\phi_{PA})u_{\mathrm{comp}} \right] \\ k_y &=& \frac{\pi}{\sqrt{2\ln(2)}} \left[-\sin(\phi_{PA})v_{\mathrm{comp}} + \cos(\phi_{PA})u_{\mathrm{comp}} \right] \end{eqnarray*}\end{split}\]where \( B_{n_k,n_l} \) is a stored 2D shapelet basis function from a look up table, of order \( n_k,n_l \) and scaled by the major and minor axis, \(C_{n_k,n_l}\) is a scaling coefficient, and \( u_{\mathrm{comp}}, v_{\mathrm{comp}} \) are visibility coordinates with the SHAPELET ra,dec as the phase centre. These should have been calculated using
fundamental_coords::kern_calc_uvw_shapelet.This kernel differs from
kern_calc_visi_point_or_gaussin that each SHAPELET component is built up from multiple shapelet basis functions, and so the kernel (and sometimes the sky model) is split over the shapelet basis functions, meaning multiple basis function calculations will use the same COMPONENT \(l,m,n\). The arrayd_components.param_indexesis used to match the basis function informationd_components.n1s,d_components.n2s,d_components.shape_coeffsto the COPMONENT information (e.g.d_components.ls).Sets off a thread for each visibility to be calculate, with each thread looping over all cofficients. This seems to keep the memory access low enough to be faster than having a second dimension over coefficient.
d_sbfis the shapelet basis function lookup table, and should have been generated withshapelet_basis.create_sbfand copied into device memory.When called with
dim3 grid, threads, kernel should be called withgrid.xdefined, where:grid.x * threads.x >=
num_visi
- Parameters
d_components – [in] d_components Pointer to a populated
components_tstruct as filled bysource_components::source_component_commond_component_beam_gains – [in] Pointer to a populated
d_beam_gains_tstruct as filled bysource_components::source_component_common*d_us – [in] Visibility coord \(u\) for every baseline, frequency, and time step in the simulation (wavelengths)
*d_vs – [in] Visibility coord \(v\) for every baseline, frequency, and time step in the simulation (wavelengths)
*d_ws – [in] Visibility coord \(w\) for every baseline, frequency, and time step in the simulation (wavelengths)
*d_allsteps_wavelengths – [in] Wavelength for every baseline, frequency, and time step in the simulation
*d_u_shapes – [in] Array of \( u_{\mathrm{comp}} \) for every baseline, frequency, and time step in the simulation (wavelengths)
*d_v_shapes – [in] Array of \( v_{\mathrm{comp}} \) for every baseline, frequency, and time step in the simulation (wavelengths)
*d_w_shapes – [in] Array of \( w_{\mathrm{comp}} \) for every baseline, frequency, and time step in the simulation (wavelengths)
*d_sum_visi_XX_real – [inout] Pointer to array to sum real XX visibility into
*d_sum_visi_XX_imag – [inout] Pointer to array to sum imaginary XX visibility into
*d_sum_visi_XY_real – [inout] Pointer to array to sum real XY visibility into
*d_sum_visi_XY_imag – [inout] Pointer to array to sum imaginary XY visibility into
*d_sum_visi_YX_real – [inout] Pointer to array to sum real YX visibility into
*d_sum_visi_YX_imag – [inout] Pointer to array to sum imaginary YX visibility into
*d_sum_visi_YY_real – [inout] Pointer to array to sum real YY visibility into
*d_sum_visi_YY_imag – [inout] Pointer to array to sum imaginary YY visibility into time step in the simulation
*d_sbf – [in] Shapelet basis function lookup table
num_shapes – [in] Number of SHAPELETs components
num_baselines – [in] Number of baselines for one time, one frequency step
num_freqs – [in] Number of frequencies in simulation
num_visis – [in] Overall number of visibilities (baselines*freqs*times) in simulation
num_coeffs – [in] Number of shapelet basis functions and coefficents
num_times – [in] Number of time steps in simulation
beamtype – [in] Beam type see
woden_struct_defs.e_beamtype
-
void free_d_components(components_t d_components, e_component_type comptype)¶
Frees device memory associated with
d_components, depending on what type of COMPONENT you are freeing.Frees the device memory set by
source_components::source_component_common. Different COMPONENT types use different parameters so need to knowcomptypeto free the correct memory locations- Parameters
d_components – [inout] A
components_twith device memory to be freedcomptype – [in] Which type of COMPONENT, either POINT, GAUSSIAN, or SHAPELET
-
void free_beam_gains(d_beam_gains_t d_beam_gains, e_beamtype beamtype)¶
Frees the device memory on
d_beam_gains, given the type of primary beambeamtype.Only certain primary beam models have leakage terms, so need to know the
beamtypeto free the correct locations.- Parameters
d_beam_gains – [inout] A
d_beam_gains_twith device memory to be freedcomptype – [in] An
e_beamtypedescribing the primary beam type
-
struct _d_beam_gains_t¶
- #include <source_components.h>
A struct to contain primary beam values for a give COMPONENT
Public Members
-
cuUserComplex *d_gxs = NULL¶
Device copy of North-South Beam gain values for all directions, frequencies, and times for these COMPONENTS
-
cuUserComplex *d_Dxs = NULL¶
Device copy of North-South Beam leakage values for all directions, frequencies, and times for these COMPONENTS
-
cuUserComplex *d_Dys = NULL¶
Device copy of East-West Beam leakage values for all directions, frequencies, and times for these COMPONENTS
-
cuUserComplex *d_gys = NULL¶
Device copy of East-West Beam gain values for all directions, frequencies, and times for these COMPONENTS
-
cuUserComplex *d_gxs = NULL¶