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.
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, cuUserComplex *d_primay_beam_J00, cuUserComplex *d_primay_beam_J01, cuUserComplex *d_primay_beam_J10, cuUserComplex *d_primay_beam_J11, double *d_freqs, double *d_ls, double *d_ms, double *d_ns, double *d_ras, double *d_decs, user_precision_t *azs, user_precision_t *zas, user_precision_t *sin_para_angs, user_precision_t *cos_para_angs, double *beam_has, double *beam_decs, woden_settings_t *woden_settings, beam_settings_t *beam_settings)¶
Performs necessary calculations that are common to all POINT, GAUSSIAN, and SHAPELET component types, calculating the \(l,m,n\) coordinates and primary beam values.
Uses
kern_calc_lmnto calculate \(l,m,n\). Uses eitherprimary_beam_cuda.calculate_gaussian_beam,FEE_primary_beam_cuda.calc_CUDA_FEE_beam, orprimary_beam_cuda.calculate_analytic_dipole_beamto calculate the primary beam.Device arrays
d_primay_beam_J*are used to store primary beam outputs, and should havenum_components*woden_settings->num_time_steps*woden_settings->num_freqselements memory allocated. \(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 COMPONENTs
*d_primay_beam_J00 – [in] Pointer towards array of primary beam J[0,0] (north-south gain) to be populated
*d_primay_beam_J01 – [in] Pointer towards array of primary beam J[0,1] (north-south leakage) to be populated
*d_primay_beam_J10 – [in] Pointer towards array of primary beam J[1,0] (east-west gain) to be populated
*d_primay_beam_J11 – [in] Pointer towards array of primary beam J[1,1] (east-west leakage) to be populated
d_freqs – [in] Frequencies used in the simulation
d_ls – [inout] Pointer to array to store l coords in
d_ms – [inout] Pointer to array to store m coords in
d_ns – [inout] Pointer to array to store d coords in
d_ras – [in] Array of Right Ascensions of COPMONENTs (radians)
d_decs – [in] Array of Declinations of COPMONENTs (radians)
azs – [in] Azimuth angles for all time steps for all COMPONENTs (radians)
zas – [in] Zenith Angles for all time steps for all COMPONENTs
sin_para_angs – [in] Sine of parallactic angle for all time steps for all COMPONENTs
cos_para_angs – [in] Sine of parallactic angle for all time steps for all COMPONENTs
beam_has – [in] Hour angles for all time steps for all COMPONENTs (radians)
beam_decs – [in] Declinations for all time steps for all COMPONENTs (radians)
woden_settings – [in] A populated
woden_settings_tstructbeam_settings – [in] A populated
beam_settings_tstruct
-
__global__ void kern_calc_visi_point(double *d_point_freqs, user_precision_t *d_point_stokesI, user_precision_t *d_point_stokesQ, user_precision_t *d_point_stokesU, user_precision_t *d_point_stokesV, user_precision_t *d_point_SIs, 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, double *d_ls, double *d_ms, double *d_ns, int num_points, int num_baselines, int num_freqs, int num_visis, 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)¶
Kernel to calculate the visibility response to a number
num_pointof POINT 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 calculate, 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.When called with
dim3 grid, threads, kernel should be called withgrid.xdefined, where:grid.x * threads.x >=
num_visi
- Parameters
*d_point_freqs – [in] Frequencies in the simulation (Hz)
*d_point_stokesI – [in] Array of reference Stokes I flux densities for all POINTs (Jy)
*d_point_stokesQ – [in] Array of reference Stokes Q flux densities for all POINTs (Jy)
*d_point_stokesU – [in] Array of reference Stokes U flux densities for all POINTs (Jy)
*d_point_stokesV – [in] Array of reference Stokes V flux densities for all POINTs (Jy)
d_point_SIs – [in] Array of spectral indicies for all POINTs
*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
*d_ls – [in] \(l\) sky coord for all time steps for all POINTs
*d_ms – [in] \(m\) sky coord for all time steps for all POINTs
*d_ns – [in] \(n\) sky coord for all time steps for all POINTs
num_points – [in] Number of POINTs
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_beamtype*d_primay_beam_J00 – [in] Array of primary beam J[0,0] (north-south gain)
*d_primay_beam_J01 – [in] Array of primary beam J[0,1] (north-south leakage)
*d_primay_beam_J10 – [in] Array of primary beam J[1,0] (east-west leakage)
*d_primay_beam_J11 – [in] Array of primary beam J[1,1] (east-west gain)
-
__global__ void kern_calc_visi_gaussian(double *d_gauss_freqs, user_precision_t *d_gauss_stokesI, user_precision_t *d_gauss_stokesQ, user_precision_t *d_gauss_stokesU, user_precision_t *d_gauss_stokesV, user_precision_t *d_gauss_SIs, 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, double *d_ls, double *d_ms, double *d_ns, user_precision_t *d_gauss_pas, user_precision_t *d_gauss_majors, user_precision_t *d_gauss_minors, int num_gauss, int num_baselines, int num_freqs, int num_visis, 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)¶
Kernel to calculate the visibility response to a number
num_gaussof GAUSS 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 GAUSSIAN:\[\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.
Sets off a thread for each visibility to be calculate, 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.
When called with
dim3 grid, threads, kernel should be called withgrid.xdefined, where:grid.x * threads.x >=
num_visi
- Parameters
*d_gauss_freqs – [in] Frequencies in the simulation (Hz)
*d_gauss_stokesI – [in] Array of reference Stokes I flux densities for all GAUSSIANs (Jy)
*d_gauss_stokesQ – [in] Array of reference Stokes Q flux densities for all GAUSSIANs (Jy)
*d_gauss_stokesU – [in] Array of reference Stokes U flux densities for all GAUSSIANs (Jy)
*d_gauss_stokesV – [in] Array of reference Stokes V flux densities for all GAUSSIANs (Jy)
d_gauss_SIs – [in] Array of spectral indicies for all GAUSSIANs
*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
*d_ls – [in] \(l\) sky coord for all time steps for all GAUSSIANs
*d_ms – [in] \(m\) sky coord for all time steps for all GAUSSIANs
*d_ns – [in] \(n\) sky coord for all time steps for all GAUSSIANs
*d_gauss_pas – [in] Position angles for all GAUSSIANs (radians)
*d_gauss_majors – [in] Major axis for all GAUSSIANs (FWHM, radians)
*d_gauss_minors – [in] Minor axis for all GAUSSIANs (FWHM, radians)
num_gauss – [in] Number 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_beamtype*d_primay_beam_J00 – [in] Array of primary beam J[0,0] (north-south gain)
*d_primay_beam_J01 – [in] Array of primary beam J[0,1] (north-south leakage)
*d_primay_beam_J10 – [in] Array of primary beam J[1,0] (east-west leakage)
*d_primay_beam_J11 – [in] Array of primary beam J[1,1] (east-west gain)
-
__global__ void kern_calc_visi_shapelets(double *d_shape_freqs, user_precision_t *d_shape_stokesI, user_precision_t *d_shape_stokesQ, user_precision_t *d_shape_stokesU, user_precision_t *d_shape_stokesV, user_precision_t *d_shape_SIs, 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_shape_pas, user_precision_t *d_shape_majors, user_precision_t *d_shape_minors, user_precision_t *d_shape_n1s, user_precision_t *d_shape_n2s, user_precision_t *d_shape_coeffs, user_precision_t *d_shape_param_indexes, double *d_ls, double *d_ms, double *d_ns, user_precision_t *d_sbf, int num_shapes, int num_baselines, int num_freqs, int num_visis, const int num_coeffs, 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)¶
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_pointandkern_calc_visi_gaussianin 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_shape_param_indexesis used to match the basis function informationd_shape_n1s,d_shape_n2s,d_shape_coeffsto the COPMONENT information (e.g.d_shape_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_shape_freqs – [in] Frequencies in the simulation (Hz)
*d_shape_stokesI – [in] Array of reference Stokes I flux densities for all SHAPELETs (Jy)
*d_shape_stokesQ – [in] Array of reference Stokes Q flux densities for all SHAPELETs (Jy)
*d_shape_stokesU – [in] Array of reference Stokes U flux densities for all SHAPELETs (Jy)
*d_shape_stokesV – [in] Array of reference Stokes V flux densities for all SHAPELETs (Jy)
d_shape_SIs – [in] Array of spectral indicies for all SHAPELETs
*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_shape_pas – [in] Position angles for all SHAPELETs (radians)
*d_shape_majors – [in] \(\beta_1 \) scaling (major axis) for all SHAPELETs (radians)
*d_shape_minors – [in] \(\beta_2 \) scaling (minor axis for all SHAPELETs (radians)
*d_shape_n1s – [in] Array of first order of 2D shapelet basis function
*d_shape_n2s – [in] Array of second order of 2D shapelet basis function
*d_shape_coeffs – [in] Array of scaling coefficients for basis functions
*d_shape_param_indexes – [in] Array of index mapping basis functions to COMPONENT parameters
*d_ls – [in] \(l\) sky coord for all time steps for all SHAPELETs
*d_ms – [in] \(m\) sky coord for all time steps for all SHAPELETs
*d_ns – [in] \(n\) sky coord for all time steps for all SHAPELETs
*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*d_primay_beam_J00 – [in] Array of primary beam J[0,0] (north-south gain)
*d_primay_beam_J01 – [in] Array of primary beam J[0,1] (north-south leakage)
*d_primay_beam_J10 – [in] Array of primary beam J[1,0] (east-west leakage)
*d_primay_beam_J11 – [in] Array of primary beam J[1,1] (east-west gain)