source_components_gpu¶
API documentation for source_components_gpu.cpp
GPU methods to extrapolate flux densities, assign and apply primary beam gains, and visibility responses for different sky model COMPONENT types.
Functions
-
__device__ gpuUserComplex calc_measurement_equation_gpu(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
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.u = d_us[iBaseline] v = d_vs[iBaseline] w = d_ws[iBaseline] l = d_ns[iComponent] m = d_ms[iComponent] n = d_ns[iComponent]
- 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, agpuUserComplexof the visibility
-
__device__ void apply_beam_gains_stokesIQUV_on_cardinal_gpu(gpuUserComplex g1x, gpuUserComplex D1x, gpuUserComplex D1y, gpuUserComplex g1y, gpuUserComplex g2x, gpuUserComplex D2x, gpuUserComplex D2y, gpuUserComplex g2y, user_precision_t flux_I, user_precision_t flux_Q, user_precision_t flux_U, user_precision_t flux_V, gpuUserComplex visi_component, gpuUserComplex *visi_XX, gpuUserComplex *visi_XY, gpuUserComplex *visi_YX, gpuUserComplex *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 apply_beam_gains_stokesI_on_cardinal_gpu(gpuUserComplex g1x, gpuUserComplex D1x, gpuUserComplex D1y, gpuUserComplex g1y, gpuUserComplex g2x, gpuUserComplex D2x, gpuUserComplex D2y, gpuUserComplex g2y, user_precision_t flux_I, gpuUserComplex visi_component, gpuUserComplex *visi_XX, gpuUserComplex *visi_XY, gpuUserComplex *visi_YX, gpuUserComplex *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 I parameter 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{eqnarray*} \mathrm{V}^{XX}_{12} = (g_{1x}g_{2x}^{\ast} + D_{1x}D_{2x}^{\ast})\mathrm{V}^{I}_{12} \end{eqnarray*}\]\[\begin{eqnarray*} \mathrm{V}^{XY}_{12} = (g_{1x}D_{2y}^{\ast} + D_{1x}g_{2y}^{\ast})\mathrm{V}^{I}_{12} \end{eqnarray*}\]\[\begin{eqnarray*} \mathrm{V}^{YX}_{12} = (D_{1y}g_{2x}^{\ast} + g_{1y}D_{2x}^{\ast})\mathrm{V}^{I}_{12} \end{eqnarray*}\]\[\begin{eqnarray*} \mathrm{V}^{YY}_{12} = (D_{1y}D_{2y}^{\ast} + g_{1y}g_{2y}^{\ast})\mathrm{V}^{I}_{12} \end{eqnarray*}\]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} \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)
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 apply_beam_gains_stokesIQUV_off_cardinal_gpu(gpuUserComplex g1x, gpuUserComplex D1x, gpuUserComplex D1y, gpuUserComplex g1y, gpuUserComplex g2x, gpuUserComplex D2x, gpuUserComplex D2y, gpuUserComplex g2y, user_precision_t flux_I, user_precision_t flux_Q, user_precision_t flux_U, user_precision_t flux_V, gpuUserComplex visi_component, gpuUserComplex *visi_XX, gpuUserComplex *visi_XY, gpuUserComplex *visi_YX, gpuUserComplex *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-east south-west (45 deg), ‘y’ means south-east north-west (135 deg).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} D_{2x}^{\ast} + D_{1x} g_{2x}^{\ast})\mathrm{V}^{Q}_{12} \\ + (g_{1x} g_{2x}^{\ast} - D_{1x} D_{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} g_{2y}^{\ast} + D_{1x} D_{2y}^{\ast})\mathrm{V}^{Q}_{12} \\ + (g_{1x} D_{2y}^{\ast} - D_{1x} g_{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} D_{2x}^{\ast} + g_{1y} g_{2x}^{\ast})\mathrm{V}^{Q}_{12} \\ + (D_{1y} g_{2x}^{\ast} - G_{1y} D_{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} g_{2y}^{\ast} + g_{1y} D_{2y}^{\ast})\mathrm{V}^{Q}_{12} \\ + (D_{1y} D_{2y}^{\ast} - G_{1y} g_{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-east south-west (45 deg)
D1x – [in] Beam leakage antenna 1 from north-east south-west (45 deg)
D1y – [in] Beam gain antenna 1 in south-east north-west (135 deg)
g1y – [in] Beam leakage antenna 1 from south-east north-west (135 deg)
g2x – [in] Beam gain antenna 2 in north-east south-west (45 deg)
D2x – [in] Beam leakage antenna 2 from north-east south-west (45 deg)
D2y – [in] Beam gain antenna 2 in south-east north-west (135 deg)
g2y – [in] Beam leakage antenna 2 from south-east north-west (135 deg)
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 apply_beam_gains_stokesI_off_cardinal_gpu(gpuUserComplex g1x, gpuUserComplex D1x, gpuUserComplex D1y, gpuUserComplex g1y, gpuUserComplex g2x, gpuUserComplex D2x, gpuUserComplex D2y, gpuUserComplex g2y, user_precision_t flux_I, gpuUserComplex visi_component, gpuUserComplex *visi_XX, gpuUserComplex *visi_XY, gpuUserComplex *visi_YX, gpuUserComplex *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 I parameter of a source, simulate the observed XX,XY,YX,YY instrumental cross-correlated visibilities, where ‘x’ means aligned north-east south-west (45 deg), ‘y’ means south-east north-west (135 deg).Performs the following calculations:
\[\begin{eqnarray*} \mathrm{V}^{XX}_{12} = (g_{1x} g_{2x}^{\ast} + D_{1x} D_{2x}^{\ast})\mathrm{V}^{I}_{12} \end{eqnarray*}\]\[\begin{eqnarray*} \mathrm{V}^{XY}_{12} = (g_{1x} D_{2y}^{\ast} + D_{1x} g_{2y}^{\ast})\mathrm{V}^{I}_{12} \end{eqnarray*}\]\[\begin{eqnarray*} \mathrm{V}^{YX}_{12} = (D_{1y} g_{2x}^{\ast} + G_{1y} D_{2x}^{\ast})\mathrm{V}^{I}_{12} \end{eqnarray*}\]\[\begin{eqnarray*} \mathrm{V}^{YY}_{12} = (D_{1y} D_{2y}^{\ast} + G_{1y} g_{2y}^{\ast})\mathrm{V}^{I}_{12} \end{eqnarray*}\]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-east south-west (45 deg)
D1x – [in] Beam leakage antenna 1 from north-east south-west (45 deg)
D1y – [in] Beam gain antenna 1 in south-east north-west (135 deg)
g1y – [in] Beam leakage antenna 1 from south-east north-west (135 deg)
g2x – [in] Beam gain antenna 2 in north-east south-west (45 deg)
D2x – [in] Beam leakage antenna 2 from north-east south-west (45 deg)
D2y – [in] Beam gain antenna 2 in south-east north-west (135 deg)
g2y – [in] Beam leakage antenna 2 from south-east north-west (135 deg)
flux_I – [in] Stokes I 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_gpu(int iBaseline, int iComponent, int num_freqs, int num_baselines, int num_components, int num_times, int beamtype, gpuUserComplex *d_gxs, gpuUserComplex *d_Dxs, gpuUserComplex *d_Dys, gpuUserComplex *d_gys, gpuUserComplex *g1x, gpuUserComplex *D1x, gpuUserComplex *D1y, gpuUserComplex *g1y, gpuUserComplex *g2x, gpuUserComplex *D2x, gpuUserComplex *D2y, gpuUserComplex *g2y)¶
Given the type of primary beam simulated
beamtype, select the beam gains and leakages from the arraysd_g*, d_D*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_gxs,d_Dxs,d_Dys,d_gysshould contain the primary beam values 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- 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_gxs – [in] Pointer towards array of primary beam J[0,0] (north-south gain)
*d_Dxs – [in] Pointer towards array of primary beam J[0,1] (north-south leakage)
*d_Dys – [in] Pointer towards array of primary beam J[1,0] (east-west leakage)
*d_gys – [in] Pointer towards array of primary beam J[1,1] (east-west gain)
*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 get_beam_gains_multibeams_gpu(int iBaseline, int iComponent, int num_freqs, int num_baselines, int num_components, int num_times, int beamtype, gpuUserComplex *d_gxs, gpuUserComplex *d_Dxs, gpuUserComplex *d_Dys, gpuUserComplex *d_gys, int *d_ant1_to_baseline_map, int *d_ant2_to_baseline_map, gpuUserComplex *g1x, gpuUserComplex *D1x, gpuUserComplex *D1y, gpuUserComplex *g1y, gpuUserComplex *g2x, gpuUserComplex *D2x, gpuUserComplex *D2y, gpuUserComplex *g2y)¶
Given the type of primary beam simulated
beamtype, select the beam gains and leakages from the arraysd_g*, d_D*that match the indexesiBaseline,iComponent. This function assumes the primary beams are different for every antenna, so returns different values for each antenna.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_gxs,d_Dxs,d_Dys,d_gysshould contain the primary beam settings for all antennas, all times, all frequencies, and COPMONENTs on the sky, and so should benum_ants*num_freqs*num_components*num_timeslong. The order elements should increment through COMPONENT (fastest changing), freqeuncy, time, and antenna (slowest changing).iBaselineis a combined index of baseline, frequency, and time.If
beamtype == NO_BEAM, setg1x = g1y = g2x = g2y = 1andD1x = D1y = D2x = D2y = 0- 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_gxs – [in] Pointer towards array of primary beam J[0,0] (north-south gain)
*d_Dxs – [in] Pointer towards array of primary beam J[0,1] (north-south leakage)
*d_Dys – [in] Pointer towards array of primary beam J[1,0] (east-west leakage)
*d_gys – [in] Pointer towards array of primary beam J[1,1] (east-west gain)
*d_ant1_to_baseline_map – [in] The index of antenna 1 in all unique pairs of antennas. Used to map iBaseline to the correct antenna 1
*d_ant2_to_baseline_map – [in] The index of antenna 1 in all unique pairs of antennas. Used to map iBaseline to the correct antenna 2
*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_stokesIQUV_gpu(int iBaseline, int iComponent, int num_freqs, int num_baselines, int num_components, int num_times, int beamtype, int off_cardinal, gpuUserComplex *d_gxs, gpuUserComplex *d_Dxs, gpuUserComplex *d_Dys, gpuUserComplex *d_gys, int *d_ant1_to_baseline_map, int *d_ant2_to_baseline_map, int use_twobeams, gpuUserComplex 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_gains_gpuandapply_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_beamtypeoff_cardinal – [in] Boolean to indicate if the dipoles are off-cardinal i.e. if off_cardinal == 1, then the dipoles are aligned at 45 and 135 degrees, if off_cardinal == 0, then the dipoles are aligned at 0 and 90 degrees
*d_gxs – [in] Pointer towards array of primary beam J[0,0] (north-south gain)
*d_Dxs – [in] Pointer towards array of primary beam J[0,1] (north-south leakage)
*d_Dys – [in] Pointer towards array of primary beam J[1,0] (east-west leakage)
*d_gys – [in] Pointer towards array of primary beam J[1,1] (east-west gain)
*d_ant1_to_baseline_map – [in] The index of antenna 1 in all unique pairs of antennas. Used to map iBaseline to the correct antenna 1 (only used when use_twobeams == 1)
*d_ant2_to_baseline_map – [in] The index of antenna 1 in all unique pairs of antennas. Used to map iBaseline to the correct antenna 2 (only used when use_twobeams == 1)
use_twobeams – [in] If 1 (True), assume all primary beams are different for each antenna. If 0 (False), assume all primary beams are the same for all antennas.
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
-
__device__ void update_sum_visis_stokesI_gpu(int iBaseline, int iComponent, int num_freqs, int num_baselines, int num_components, int num_times, int beamtype, int off_cardinal, gpuUserComplex *d_gxs, gpuUserComplex *d_Dxs, gpuUserComplex *d_Dys, gpuUserComplex *d_gys, int *d_ant1_to_baseline_map, int *d_ant2_to_baseline_map, int use_twobeams, gpuUserComplex visi_component, user_precision_t flux_I, 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 I parameter 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_gains_gpuorget_beam_gains_multibeams_gpu, andapply_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_beamtypeoff_cardinal – [in] Boolean to indicate if the dipoles are off-cardinal i.e. if off_cardinal == 1, then the dipoles are aligned at 45 and 135 degrees, if off_cardinal == 0, then the dipoles are aligned at 0 and 90 degrees
*d_gxs – [in] Pointer towards array of primary beam J[0,0] (north-south gain)
*d_Dxs – [in] Pointer towards array of primary beam J[0,1] (north-south leakage)
*d_Dys – [in] Pointer towards array of primary beam J[1,0] (east-west leakage)
*d_gys – [in] Pointer towards array of primary beam J[1,1] (east-west gain)
*d_ant1_to_baseline_map – [in] The index of antenna 1 in all unique pairs of antennas. Used to map iBaseline to the correct antenna 1 (only used when use_twobeams == 1)
*d_ant2_to_baseline_map – [in] The index of antenna 1 in all unique pairs of antennas. Used to map iBaseline to the correct antenna 2 (only used when use_twobeams == 1)
use_twobeams – [in] If 1 (True), assume all primary beams are different for each antenna. If 0 (False), assume all primary beams are the same for all antennas.
visi_component – [in] Complex visibility across antennas 1 and 2
flux_I – [in] Stokes I 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
-
__global__ void kern_make_zeros_user_precision(user_precision_t *array, int num_arr)¶
Kernel to set all elements of
arrayto zero.Call using something like
dim3 grid, threads;
threads.x = 128;
grid.x = (int)ceil( (float)num_arr / (float)threads.x );
- Parameters:
*array – [inout] Array to set to zeros
num_arr – [in] Number of elements in the array
-
void malloc_extrapolated_flux_arrays_gpu(components_t *d_components, int num_comps, int num_freqs)¶
Allocate device memory of extrapolated Stokes arrays int
d_componentsIt does this if d_components.do_QUV == 1:
If do_QUV == 0, only allocate the StokesI array.d_components->extrap_stokesI = NULL; gpuMalloc( (void**)&d_components->extrap_stokesI, num_comps*num_freqs*sizeof(double) ); d_components->extrap_stokesQ = NULL; gpuMalloc( (void**)&d_components->extrap_stokesQ, num_comps*num_freqs*sizeof(double) ); d_components->extrap_stokesU = NULL; gpuMalloc( (void**)&d_components->extrap_stokesU, num_comps*num_freqs*sizeof(double) ); d_components->extrap_stokesV = NULL; gpuMalloc( (void**)&d_components->extrap_stokesV, num_comps*num_freqs*sizeof(double) );
External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
*d_components – [inout] A populated
components_tstructnum_comps – [in] Number of components
num_freqs – [in] Number of frequencies
-
__device__ void extrap_stokes_power_law_gpu(user_precision_t *d_ref_fluxes, user_precision_t *d_SIs, double *d_extrap_freqs, int iFluxComp, int iFreq, user_precision_t *extrap_flux)¶
Assuming a simple power-law SED of \( S \propto \nu^\alpha \), extrapolate flux density parameters to the requested frequencies.
Assumes a referece frequency of 200MHz, and uses the reference flux density
d_ref_fluxes[iFluxComp], andspectral index d_SIs[iFluxComp]to calculate the flux densityextrap_fluxat the requested frequency d_extrap_freqs[iFreq].- Parameters:
d_ref_fluxes – [in] Reference fluxes for all components
d_SIs – [in] Spectral indexes for all components
d_extrap_freqs – [in] Pointer to array of frequencies to extrapolate to
iFluxComp – [in] Index of which component to extrapolate
iFreq – [in] Index of which frequency to extrapolate to
*extrap_flux – [inout] Pointer to extrapolated flux (Jy)
-
__global__ void kern_extrap_power_laws_stokesI(int num_extrap_freqs, double *d_extrap_freqs, int num_comps, components_t d_components)¶
Kernel to run
extrap_stokes_power_law_gpufor all Stokes I components ind_components.Fills the array
d_components.extrap_stokesIwith the extrapolated Stokes I flux densities.- Parameters:
num_extrap_freqs – [in] The number of frequencies to extrapolate to.
d_extrap_freqs – [in] Pointer to an array of frequencies to extrapolate to on the device.
num_comps – [in] The number of components to extrapolate.
d_components – [inout] The components to extrapolate on the device.
-
__global__ void kern_extrap_power_laws_stokesV(int num_extrap_freqs, double *d_extrap_freqs, int num_comps, components_t d_components)¶
Kernel to run
extrap_stokes_power_law_gpufor all Stokes V components ind_components.Fills the array
d_components.extrap_stokesVwith the extrapolated Stokes V flux densities.- Parameters:
num_extrap_freqs – [in] The number of frequencies to extrapolate to.
d_extrap_freqs – [in] Pointer to an array of frequencies to extrapolate to on the device.
num_comps – [in] The number of components to extrapolate.
d_components – [inout] The components to extrapolate on the device.
-
__global__ void kern_extrap_power_laws_linpol(int num_extrap_freqs, double *d_extrap_freqs, int num_comps, components_t d_components)¶
Kernel to run
extrap_stokes_power_law_gpufor all linear polarisation components ind_components.Temporarily fills the array
d_components.extrap_stokesQwith the extrapolated linear polarisation flux densities. The intention is to the usekern_apply_rotation_measureafter the fact, which will used_components.extrap_stokesQas input flux, do rotations, and then filld_components.extrap_stokesQandd_components.extrap_stokesU.- Parameters:
num_extrap_freqs – [in] The number of frequencies to extrapolate to.
d_extrap_freqs – [in] Pointer to an array of frequencies to extrapolate to on the device.
num_comps – [in] The number of components to extrapolate.
d_components – [inout] The components to extrapolate on the device.
-
__device__ void extrap_stokes_curved_power_law_gpu(user_precision_t *d_ref_fluxes, user_precision_t *d_SIs, user_precision_t *d_qs, double *d_extrap_freqs, int iFluxComp, int iFreq, user_precision_t *extrap_flux)¶
Assuming a curved power-law SED of \( S \propto \nu^\alpha \exp(q \ln (\nu)^2 ) \), extrapolate Stokes I flux density parameters to the requested frequencies.
Assumes a referece frequency of 200MHz, and uses the reference flux density
d_ref_fluxes[iFluxComp],spectral index d_SIs[iFluxComp], curvatureqparamspectral index d_qs[iFluxComp]to calculate the flux densityextrap_fluxat the requested frequency d_extrap_freqs[iFreq].- Parameters:
d_ref_fluxes – [in] Reference fluxes for all components
d_SIs – [in] Spectral indexes for all components
d_qs – [in] Curvature
qparam for all componentsd_extrap_freqs – [in] Pointer to array of frequencies to extrapolate to
iFluxComp – [in] Index of which component to extrapolate
iFreq – [in] Index of which frequency to extrapolate to
*extrap_flux – [inout] Pointer to extrapolated flux (Jy)
-
__global__ void kern_extrap_curved_power_laws_stokesI(int num_extrap_freqs, double *d_extrap_freqs, int num_comps, components_t d_components)¶
Kernel to run
extrap_stokes_curved_power_law_gpufor all Stokes I components ind_components.Fills the array
d_components.extrap_stokesIwith the extrapolated Stokes I flux densities.- Parameters:
num_extrap_freqs – [in] The number of frequencies to extrapolate to.
d_extrap_freqs – [in] Pointer to an array of frequencies to extrapolate to on the device.
num_comps – [in] The number of components to extrapolate.
d_components – [inout] The components to extrapolate on the device.
-
__global__ void kern_extrap_curved_power_laws_stokesV(int num_extrap_freqs, double *d_extrap_freqs, int num_comps, components_t d_components)¶
Kernel to run
extrap_stokes_curved_power_law_gpufor all Stokes I components ind_components.Fills the array
d_components.extrap_stokesVwith the extrapolated Stokes V flux densities.- Parameters:
num_extrap_freqs – [in] The number of frequencies to extrapolate to.
d_extrap_freqs – [in] Pointer to an array of frequencies to extrapolate to on the device.
num_comps – [in] The number of components to extrapolate.
d_components – [inout] The components to extrapolate on the device.
-
__global__ void kern_extrap_curved_power_laws_linpol(int num_extrap_freqs, double *d_extrap_freqs, int num_comps, components_t d_components)¶
Kernel to run
extrap_stokes_curved_power_law_gpufor all linear polarisation components ind_components.Temporarily fills the array
d_components.extrap_stokesQwith the extrapolated linear polarisation flux densities. The intention is to the usekern_apply_rotation_measureafter the fact, which will used_components.extrap_stokesQas input flux, do rotations, and then filld_components.extrap_stokesUandd_components.extrap_stokesV.- Parameters:
num_extrap_freqs – [in] The number of frequencies to extrapolate to.
d_extrap_freqs – [in] Pointer to an array of frequencies to extrapolate to on the device.
num_comps – [in] The number of components to extrapolate.
d_components – [inout] The components to extrapolate on the device.
-
__global__ void kern_polarisation_fraction_stokesV(int num_extrap_freqs, double *d_extrap_freqs, int num_comps, components_t d_components)¶
Kernel to calculate polarisation fractions for all Stokes V components in
d_components.MUST have the
d_components.extrap_stokesIfilled already. Fills the arrayd_components.extrap_stokesVusingd_components.extrap_stokesIandd_components.stokesV_pol_fracs.- Parameters:
num_extrap_freqs – [in] The number of frequencies to extrapolate to.
d_extrap_freqs – [in] Pointer to an array of frequencies to extrapolate to on the device.
num_comps – [in] The number of components to extrapolate.
d_components – [inout] The components to extrapolate on the device.
-
__global__ void kern_polarisation_fraction_linpol(int num_extrap_freqs, double *d_extrap_freqs, int num_comps, components_t d_components)¶
Kernel to calculate polarisation fractions for all linear polarisation components in
d_components.MUST have the
d_components.extrap_stokesIfilled already. Temporarily fills the arrayd_components.extrap_stokesQusingd_components.extrap_stokesIandd_components.linpol_pol_fracs. The intention is to the usekern_apply_rotation_measureafter the fact, which will used_components.extrap_stokesQas input flux, do rotations, and then filld_components.extrap_stokesUandd_components.extrap_stokesV.- Parameters:
num_extrap_freqs – [in] The number of frequencies to extrapolate to.
d_extrap_freqs – [in] Pointer to an array of frequencies to extrapolate to on the device.
num_comps – [in] The number of components to extrapolate.
d_components – [inout] The components to extrapolate on the device.
-
__device__ void extrap_stokes_list_fluxes_gpu(user_precision_t *list_stokes, double *list_freqs, int *arr_num_list_values, int *list_start_indexes, double *d_extrap_freqs, int iFluxComp, int iFreq, user_precision_t *extrap_flux)¶
Assuming a list-type spectral model, extrapolates the flux for a given component and frequency.
- Parameters:
list_stokes – [in] Array containing all the list-fluxes used for the extrapolation.
list_freqs – [in] Array containing all the list-frequencies used for the extrapolation; must match order of
list_stokes.arr_num_list_values – [in] Array containing the number of list values for each component.
list_start_indexes – [in] Array containing the starting index of each component within
list_stokesandlist_freqs.d_extrap_freqs – [in] The frequencies to use for the extrapolation.
iFluxComp – [in] The index of the flux component to store the result in.
iFreq – [in] The index of the frequency component to extrapolate.
extrap_flux – [inout] The extrapolated flux.
-
__global__ void kern_extrap_list_fluxes(user_precision_t *list_stokes, double *list_freqs, int *num_list_values, int *list_start_indexes, int *list_comp_inds, int num_extrap_freqs, double *d_extrap_freqs, int num_comps, user_precision_t *extrap_stokes)¶
Extrapolates list-type spectral model fluxes to given frequencies. To be clear,
list_stokesandlist_freqsare arrays of flux densities and frequencies, respectively, fornum_compsnumber of components. Each component can have any number of flux and freq entries, given bynum_list_values.Fills
extrap_stokeswith the extrapolated stokes flux densities. Runs the functionextrap_stokes_list_fluxes_gpu.- Parameters:
list_stokes – [in] Array containing all the list-fluxes used for the extrapolation.
list_freqs – [in] Array containing all the list-frequencies used for the extrapolation; must match order of
list_stokes.num_list_values – [in] Array containing the number of list values for each component.
list_start_indexes – [in] Array containing the starting index of each component within
list_stokesandlist_freqs.list_comp_inds – [in] Array containing the component index for list component; used to index the extrapolated fluxes to
extrap_stokes.num_extrap_freqs – [in] The number of extrapolation frequencies.
d_extrap_freqs – [in] Pointer to the array of extrapolation frequencies.
num_comps – [in] The number of components.
extrap_stokes – [inout] Output extrapolated fluxes
-
__global__ void kern_apply_rotation_measure(int num_extrap_freqs, double *d_extrap_freqs, int num_comps, components_t d_components)¶
Apply rotation measure to existing fluxes in
components.extrap_stokesQ.Fills
d_components.extrap_stokesQandd_components.extrap_stokesUwith the rotated linear polarisation flux densities. Assumes the model\[\begin{split} Q = \mathrm{P}(\lambda) \cos(2\chi_0 + 2\phi_{\textrm{RM}} \lambda^2), \\ U = \mathrm{P}(\lambda) \sin(2\chi_0 + 2\phi_{\textrm{RM}} \lambda^2). \end{split}\]where \(\mathrm{P}(\lambda)\) is the polarised flux, stored in
components.extrap_stokesQ, \(\chi_0\) is the intrinsic polarisation angle (components.intr_pol_angle), and \(\phi_{\textrm{RM}}\) is the rotation measure (components.rm_values).- Parameters:
num_extrap_freqs – The number of extrapolation frequencies.
d_extrap_freqs – The frequencies to extrapolate to.
num_comps – The number of components.
d_components – The components to apply rotation measure to
-
void extrap_power_laws_stokesI_gpu(components_t d_components, int n_powers, double *d_extrap_freqs, int num_extrap_freqs)¶
Extrapolates power law SEDSs for Stokes I flux densities to given frequencies on the GPU.
Calls
kern_extrap_power_laws_stokesIto filld_components.extrap_stokesIwith the extrapolated flux densities with inputs fromd_components.External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
d_components – The device components structure containing necessary data.
n_powers – The number of power laws to extrapolate.
d_extrap_freqs – Pointer to the device array of extrapolation frequencies.
num_extrap_freqs – The number of extrapolated frequencies.
-
void extrap_curved_power_laws_stokesI_gpu(components_t d_components, int n_curves, double *d_extrap_freqs, int num_extrap_freqs)¶
Extrapolates curved power law SEds for Stokes I flux densities to given frequencies on the GPU.
Calls
kern_extrap_curved_power_laws_stokesIto filld_components.extrap_stokesIwith the extrapolated flux densities with inputs fromd_components.External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
d_components – The device components structure containing necessary data for extrapolation.
n_curves – The number of curves to be extrapolated.
d_extrap_freqs – Pointer to the device array of extrapolation frequencies.
num_extrap_freqs – The number of frequencies in the d_extrap_freqs array.
-
void extrap_list_fluxes_stokesI_gpu(components_t d_components, int n_lists, double *d_extrap_freqs, int num_extrap_freqs)¶
Extrapolate list-type SEDs for Stokes I flux densities to given frequencies on GPU.
Calls
kern_extrap_list_fluxesto filld_components.extrap_stokesIwith the extrapolated flux densities with inputs fromd_components.External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
d_components – The device pointer to the components data structure.
n_lists – The number of lists to process.
d_extrap_freqs – The device pointer to the array of frequencies at which to extrapolate.
num_extrap_freqs – The number of frequencies in the d_extrap_freqs array.
-
void extrap_power_laws_stokesV_gpu(components_t d_components, double *d_extrap_freqs, int num_extrap_freqs)¶
Extrapolates power law SEDSs for Stokes V flux densities to given frequencies on the GPU.
Calls
kern_extrap_power_laws_stokesVto filld_components.extrap_stokesVwith the extrapolated flux densities with inputs fromd_components.External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
d_components – The device components structure containing necessary data.
d_extrap_freqs – Pointer to the device array of extrapolation frequencies.
num_extrap_freqs – The number of extrapolated frequencies.
-
void extrap_curved_power_laws_stokesV_gpu(components_t d_components, double *d_extrap_freqs, int num_extrap_freqs)¶
Extrapolates curved power law SEds for Stokes V flux densities to given frequencies on the GPU.
Calls
kern_extrap_curved_power_laws_stokesVto filld_components.extrap_stokesVwith the extrapolated flux densities with inputs fromd_components.External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
d_components – The device components structure containing necessary data for extrapolation.
d_extrap_freqs – Pointer to the device array of extrapolation frequencies.
num_extrap_freqs – The number of frequencies in the d_extrap_freqs array.
-
void polarisation_fraction_stokesV_gpu(components_t d_components, double *d_extrap_freqs, int num_extrap_freqs)¶
Calculates polarisations fraction for Stokes V flux densities to given frequencies on the GPU.
Calls
kern_polarisation_fraction_stokesVto filld_components.extrap_stokesVwith the extrapolated flux densities with inputs fromd_components. Requiresd_components.extrap_stokesIto be filled first.External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
d_components – The device components structure containing necessary data for extrapolation.
d_extrap_freqs – Pointer to the device array of extrapolation frequencies.
num_extrap_freqs – The number of frequencies in the d_extrap_freqs array.
-
void extrap_list_fluxes_stokesV_gpu(components_t d_components, double *d_extrap_freqs, int num_extrap_freqs)¶
Extrapolate list-type SEDs for Stokes V flux densities to given frequencies on GPU.
Calls
kern_extrap_list_fluxesto filld_components.extrap_stokesVwith the extrapolated flux densities with inputs fromd_components.External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
d_components – The device pointer to the components data structure.
d_extrap_freqs – The device pointer to the array of frequencies at which to extrapolate.
num_extrap_freqs – The number of frequencies in the d_extrap_freqs array.
-
void extrap_power_laws_linpol_gpu(components_t d_components, double *d_extrap_freqs, int num_extrap_freqs)¶
Extrapolates power law SEDSs for linear polarisation flux given frequencies on the GPU.
Calls
kern_extrap_power_laws_linpolto filld_components.extrap_stokesQwith the extrapolated flux densities with inputs fromd_components. Once called, the functionkern_apply_rotation_measureshould be called to rotate the linear polarisation fluxes into bothd_components.extrap_stokesQandd_components.extrap_stokesU.External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
d_components – The device components structure containing necessary data.
d_extrap_freqs – Pointer to the device array of extrapolation frequencies.
num_extrap_freqs – The number of extrapolated frequencies.
-
void extrap_curved_power_laws_linpol_gpu(components_t d_components, double *d_extrap_freqs, int num_extrap_freqs)¶
Extrapolates curved power law SEDs for linear polarisation flux for given frequencies on the GPU.
Calls
kern_extrap_curved_power_laws_linpolto fill to filld_components.extrap_stokesQwith the extrapolated flux densities with inputs fromd_components. Once called, the functionkern_apply_rotation_measureshould be called to rotate the linear polarisation fluxes into bothd_components.extrap_stokesQandd_components.extrap_stokesU.External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
d_components – The device components structure containing necessary data for extrapolation.
d_extrap_freqs – Pointer to the device array of extrapolation frequencies.
num_extrap_freqs – The number of frequencies in the d_extrap_freqs array.
-
void polarisation_fraction_linpol_gpu(components_t d_components, double *d_extrap_freqs, int num_extrap_freqs)¶
Calculates polarisations fraction for linear polarisation flux on the GPU.
Calls
kern_polarisation_fraction_linpolto filld_components.extrap_stokesQwith the extrapolated flux densities with inputs fromd_components. Requiresd_components.extrap_stokesIto be filled first. Once called, the functionkern_apply_rotation_measureshould be called to rotate the linear polarisation fluxes into bothd_components.extrap_stokesQandd_components.extrap_stokesU.External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
d_components – The device components structure containing necessary data for extrapolation.
d_extrap_freqs – Pointer to the device array of extrapolation frequencies.
num_extrap_freqs – The number of frequencies in the d_extrap_freqs array.
-
void extrap_list_fluxes_linpol_gpu(components_t d_components, double *d_extrap_freqs, int num_extrap_freqs)¶
Extrapolate list-type SEDs for linear polarisation flux for given frequencies on the GPU.
Calls
kern_extrap_list_fluxestwice to separately filld_components.extrap_stokesQandd_components.extrap_stokesU. Assumes information for boths Q and U list type fluxes are stored ind_components. Does not required rotation measure to be applied as this model type just requires lists of Q and U fluxes.External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
d_components – The device pointer to the components data structure.
d_extrap_freqs – The device pointer to the array of frequencies at which to extrapolate.
num_extrap_freqs – The number of frequencies in the d_extrap_freqs array.
-
void extrap_p_list_fluxes_linpol_gpu(components_t d_components, double *d_extrap_freqs, int num_extrap_freqs)¶
Extrapolate list-type SEDs for linear polarisation flux for given frequencies on the GPU.
Calls
kern_extrap_list_fluxesto filld_components.extrap_stokesQwith the extrapolated flux densities with inputs fromd_components. Assumeslinpol_ptype list flux information exists ind_components. Once called, the functionkern_apply_rotation_measureshould be called to rotate the linear polarisation fluxes into bothd_components.extrap_stokesQandd_components.extrap_stokesU.External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
d_components – The device components structure containing necessary data.
d_extrap_freqs – Pointer to the device array of extrapolation frequencies.
num_extrap_freqs – The number of extrapolated frequencies.
-
void apply_rotation_measure_gpu(components_t d_components, double *d_extrap_freqs, int num_extrap_freqs)¶
Applies rotation measure to the given components on the GPU.
Calls
kern_apply_rotation_measure. Assumesd_components.extrap_stokesQhas been filled with the total linear polarisation flux. The function will filld_components.extrap_stokesQandd_components.extrap_stokesUwith the rotated linear polarisation fluxes.External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
d_components – The device pointer to the components structure.
d_extrap_freqs – The device pointer to the array of extrapolated frequencies.
num_extrap_freqs – The number of extrapolated frequencies.
-
__global__ void kern_calc_visi_point_or_gauss(components_t d_components, 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, int num_components, int num_baselines, int num_freqs, int num_cross, int num_times, e_beamtype beamtype, e_component_type comptype, int off_cardinal_dipoles)¶
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. Assumes all fluxes have been extrapolated to the requested frequencies.Uses the functions
calc_measurement_equation_gpuandupdate_sum_visis_stokes*gpuas 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
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
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_cross – [in] Overall number of cross-correlations (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
off_cardinal_dipoles – [in] Boolean to indicate if the dipoles are off-cardinal i.e. if off_cardinal == 1, then the dipoles are aligned at 45 and 135 degrees, if off_cardinal == 0, then the dipoles are aligned at 0 and 90 degrees
-
__global__ void kern_calc_visi_shapelets(components_t d_components, 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_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_cross, const int num_coeffs, int num_times, e_beamtype beamtype, int off_cardinal_dipoles)¶
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_gpu,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_gpu::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).A further difference is that the shapelet \( u_{\mathrm{comp}}, v_{\mathrm{comp}} \) used in the visibility envelope equation should be stored in metres only - for every baseline (fastest changing), every time step, and every SHAPELET component (slowest changing). They will be scaled by wavelength inside the kernel, as memory costs become prohibitive to store for every wavelength as well. The visibility \(u,v,w\) are stored for every baseline, wavelength, and time step.
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
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 SHAPELET component in the simulation (metres)
*d_v_shapes – [in] Array of \( v_{\mathrm{comp}} \) for every baseline, frequency, and SHAPELET component in the simulation (metres)
*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_cross – [in] Overall number of cross-correlations (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_beamtypeoff_cardinal_dipoles – [in] Boolean to specify if the dipoles in the beam are off the cardinal axes (45 and 135 degrees) or aligned with north-south and east-west (0 and 90 degrees). Effects what visibilities are calculated.
-
void copy_components_to_GPU(source_t *chunked_source, source_t *d_chunked_source, e_component_type comptype)¶
Copies the specified type of source components from host memory to device memory.
This function copies the specified type of source components from the host memory pointed to by
chunked_sourceto the device memory pointed to byd_chunked_source, allocating the device memory in the process Thecomptypeparameter specifies the type of components to copy.- Parameters:
chunked_source – [in] Pointer to the host memory containing the source components to copy.
d_chunked_source – [inout] Pointer to the device memory where the source components will be copied.
comptype – [in] The type of source components to copy.
-
void free_components_gpu(source_t *d_chunked_source, e_component_type comptype)¶
Frees device memory associated with
d_chunked_source, depending on what type of COMPONENT you are freeing.Either frees all device memory that will have been allocated from either
d_chunked_source->point_components,d_chunked_source->gauss_components, ord_chunked_source->shape_components, depending oncomptype.External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
*d_chunked_source – [inout] A
source_twith device memory to be freedcomptype – [in] Which type of COMPONENT, either POINT, GAUSSIAN, or SHAPELET
-
void free_beam_gains_gpu(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.External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
d_beam_gains – [inout] A
beam_gains_twith device memory to be freedbeamtype – [in] An
e_beamtypedescribing the primary beam type
-
source_t *copy_chunked_source_to_GPU(source_t *chunked_source)¶
Copies a populated
source_tfrom host to device.Depending on what is in chunked_source (e.g. POINT, GAUSSIAN, SHAPELET, POWER_LAW, CURVED_POWER_LAW, LIST) type information, only what needs to be copied across to the GPU is (no empty pointers arrays are copied across)
External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
*chunked_source – [inout] A populated
source_tstruct
-
void copy_CPU_component_gains_to_GPU_beam_gains(components_t *components, beam_gains_t *d_beam_gains, int num_gains)¶
Copies component gains from CPU to GPU beam gains. This function transfers the component gains stored in the CPU memory to the GPU memory, specifically to the beam gains structure. It is used to ensure that the GPU has the necessary data to perform computations involving beam gains.
Specifically, these member arrays
gxs, Dxs, Dys, gyshave memory allocated on GPU ind_beam_gainsand are copied from thecomponentsCPU memory.External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
components – Pointer to the components structure containing the gains on the CPU.
d_beam_gains – Pointer to the beam gains structure on the GPU where the gains will be copied.
num_gains – The number of gains to be copied from the CPU to the GPU.
-
void copy_CPU_beam_gains_to_GPU_beam_gains(beam_gains_t *beam_gains, beam_gains_t *d_beam_gains, int num_gains)¶
Copies beam gains data from CPU memory to GPU memory.
Specifically, these member arrays
gxs, Dxs, Dys, gyshave memory allocated on GPU ind_beam_gainsand are copied from thebeam_gainsCPU memory.External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
beam_gains – Pointer to the beam gains data in CPU memory.
d_beam_gains – Pointer to the beam gains data in GPU memory.
num_gains – The number of beam gains to copy.
-
void free_extrapolated_flux_arrays_gpu(components_t *d_components)¶
Free device memory of extrapolated Stokes arrays from
d_componentsIt does this if d_components.do_QUV == 1:
gpuFree d_components->extrap_stokesI ); gpuFree d_components->extrap_stokesQ ); gpuFree d_components->extrap_stokesU ); gpuFree d_components->extrap_stokesV );
If do_QUV == 0, only free the StokesI array.
External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
*d_components – [inout] A populated
components_tstruct
-
__global__ void kern_calc_autos(components_t d_components, beam_gains_t d_component_beam_gains, int beamtype, int num_components, int num_baselines, int num_freqs, int num_times, int num_ants, 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, int use_twobeams, int *d_ant1_to_auto_map, int *d_ant2_to_auto_map, int off_cardinal_dipoles)¶
Calculate the auto-correlations for all antennas given the fluxes already calculated in
d_componentsand beam gains ind_component_beam_gains. Stores the outputs at the end ofd_sum_visi*, after the cross-correlations.d_component_beam_gainsshould have gains for as many antennas (a.k.a stations/tiles) as detailed bynum_ants. This function then just does the dot product of the component fluxes and beam gains specific to each baseline.When called with
dim3 grid, threads, kernel should be called withgrid.xand grid.y defined, where:grid.x * threads.x >=
num_freqs*num_timesgrid.y * threads.y >=
num_ants
- 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
beam_gains_tstruct as filled bysource_components::source_component_commonbeamtype – [in] Beam type see
woden_struct_defs.e_beamtypenum_components – [in] Number of components in
d_componentsnum_baselines – [in] Number of baselines for one time, one frequency step
num_freqs – [in] Number of frequncies in the simulation
num_times – [in] Number of times in the simulation
num_ants – [in] Number of antennas in the array
*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
use_twobeams – [in] If True, use a two primary beams per visibility. Otherwise, assume all primary beams are identical
*d_ant1_to_auto_map – [in] An index of all primary beams to auto-correlations Currently this is just an index of all antennas. Gets passed to
get_beam_gains_multibeams_gpu*d_ant2_to_auto_map – [in] An index of all primary beams to auto-correlations Currently this is just an index of all antennas. Gets passed to
get_beam_gains_multibeams_gpuoff_cardinal_dipoles – [in] Boolean to specify if the dipoles in the beam are off the cardinal axes (45 and 135 degrees) or aligned with north-south and east-west (0 and 90 degrees). Effects what visibilities are calculated.
-
void fill_ant_to_baseline_mapping_gpu(int num_ants, int *d_ant1_to_baseline_map, int *d_ant2_to_baseline_map)¶
Fill the
d_ant1_to_baseline_mapandd_ant2_to_baseline_mapdevice arrays with indexes corresponding to ant1 and ant2 for all unique baselines in an array ofnum_antsantennas.The
d_ant1_to_baseline_mapandd_ant2_to_baseline_mapshould already have their memory allocated- Parameters:
num_ants – [in] Number of antennas in the array
*d_ant1_to_baseline_map – [inout] Device memory-allocated array of size
((num_ants - 1)*num_ants) / 2*d_ant2_to_baseline_map – [inout] Device memory-allocated array of size
((num_ants - 1)*num_ants) / 2
-
void calc_visi_point_or_gauss_gpu(components_t d_components, beam_gains_t d_component_beam_gains, calc_visi_inouts_t *d_calc_visi_inouts, visibility_set_t *d_visibility_set, int num_components, e_beamtype beamtype, e_component_type comptype, woden_settings_t *woden_settings)¶
Calculate the visibility response to a number
num_componentsof either POINT or GAUSSIAN COMPONENTs, and sum the outputs tosum_visi_*_real,sum_visi_*_imagind_visibility_set. Assumes all fluxes have been extrapolated to the requested frequencies.Runs
kern_calc_visi_point_or_gauss.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
External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
d_components – [in] components Pointer to a populated
components_tstruct as filled bysource_components::source_component_commond_component_beam_gains – [in] Pointer to a populated
beam_gains_tstruct as filled bysource_components::source_component_commond_calc_visi_inouts – The input/output parameters for visibility calculation.
d_visibility_set – The visibility set to update with the results.
num_components – The number of components.
beamtype – Beam type, see
woden_struct_defs.e_beamtypecomptype – Component type, either POINT or GAUSSIAN
woden_settings – The filled WODEN settings struct
-
void calc_visi_shapelets_gpu(components_t d_components, beam_gains_t d_component_beam_gains, calc_visi_inouts_t *d_calc_visi_inouts, visibility_set_t *d_visibility_set, int num_shapes, int num_shape_coeffs, e_beamtype beamtype, woden_settings_t *woden_settings)¶
the visibility response to a number
num_shapesof SHAPELET COMPONENTs, and sum the outputs tosum_visi_*_real,sum_visi_*_imagind_visibility_set. Assumes all fluxes have been extrapolated to the requested frequencies.Uses the functions
calc_measurement_equation_cpuandupdate_sum_visis_stokesI*_cpuas detailed incalc_visi_point_or_gauss_cputo calculate the visibilities in serial on the CPU. 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_gpu::kern_calc_uvw_shapelet.This differs from
calc_visi_point_or_gauss_gpuin 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 arraycomponents.param_indexesis used to match the basis function informationcomponents.n1s,components.n2s,components.shape_coeffsto the COPMONENT information (e.g.components.ls).A further difference is that the shapelet \( u_{\mathrm{comp}}, v_{\mathrm{comp}} \) used in the visibility envelope equation should be stored in metres only - for every baseline (fastest changing), every time step, and every SHAPELET component (slowest changing). They will be scaled by wavelength inside the kernel, as memory costs become prohibitive to store for every wavelength as well. The visibility \(u,v,w\) are stored for every baseline, wavelength, and time step.
sbfis the shapelet basis function lookup table, and should have been generated withcreate_sbfExternal GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
d_components – [in] components Pointer to a populated
components_tstruct as filled bysource_components::source_component_commond_component_beam_gains – [in] Pointer to a populated
beam_gains_tstruct as filled bysource_components::source_component_commond_calc_visi_inouts – The input/output parameters for visibility calculation.
d_visibility_set – The visibility set to update.
num_shapes – The number of shapelets.
num_shape_coeffs – The number of shapelet coefficients.
beamtype – The type of beam.
woden_settings – The WODEN settings.
-
void malloc_beam_gains_gpu(beam_gains_t *d_component_beam_gains, int beamtype, int num_gains)¶
Allocate memory for beam gains on GPU.
Allocates the following memory:
d_component_beam_gains->Dxsd_component_beam_gains->Dysd_component_beam_gains->gxsd_component_beam_gains->gys
depending on what beamtype is given, the leakgages may or may not be allocated.
External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
d_component_beam_gains – The beam gains struct to allocate memory in
beamtype – The type of beam.
num_gains – The number of gains.
-
void calc_autos_gpu(components_t *d_components, beam_settings_t *beam_settings, beam_gains_t *d_component_beam_gains, visibility_set_t *d_visibility_set, woden_settings_t *woden_settings, int num_components, int use_twobeams)¶
Calculate the auto-correlations for all antennas given the fluxes already calculated in
d_componentsand beam gains ind_component_beam_gains. Stores the outputs ind_visibility_set.d_component_beam_gainsshould have gains for as many antennas (a.k.a stations/tiles). This function then just does the dot product of the component fluxes and beam gains specific to each baseline.External GPU code linked in from source_components_gpu.cpp, see source_components_gpu.h for full details.
- Parameters:
d_components – Pointer to the device memory holding the components.
beam_settings – Pointer to the beam settings structure.
d_component_beam_gains – Pointer to the device memory holding the component beam gains.
d_visibility_set – Pointer to the device memory holding the visibility set.
woden_settings – Pointer to the WODEN settings structure.
num_components – The number of components to process.
use_twobeams – If True (1), use two primary beams per visibility. Otherwise, assume all primary beams are identical
-
__global__ void kern_calc_measurement_equation(int num_components, int num_baselines, user_precision_t *d_us, user_precision_t *d_vs, user_precision_t *d_ws, double *d_ls, double *d_ms, double *d_ns, gpuUserComplex *d_visis)¶
Wrapper kernel to test
calc_measurement_equation_gpuin unit tests.
-
__global__ void kern_apply_beam_gains(int num_gains, gpuUserComplex *d_g1xs, gpuUserComplex *d_D1xs, gpuUserComplex *d_D1ys, gpuUserComplex *d_g1ys, gpuUserComplex *d_g2xs, gpuUserComplex *d_D2xs, gpuUserComplex *d_D2ys, gpuUserComplex *d_g2ys, user_precision_t *d_flux_Is, user_precision_t *d_flux_Qs, user_precision_t *d_flux_Us, user_precision_t *d_flux_Vs, gpuUserComplex *d_visi_components, gpuUserComplex *d_visi_XXs, gpuUserComplex *d_visi_XYs, gpuUserComplex *d_visi_YXs, gpuUserComplex *d_visi_YYs, int off_cardinal_dipoles, int do_QUV)¶
Wrapper kernel to test
apply_beam_gains_stokesI*_off_cardinal_gpuandapply_beam_gains_stokesI*_on_cardinal_gpuin unit tests.
-
__global__ void kern_get_beam_gains(int num_components, int num_baselines, int num_freqs, int num_cross, int num_times, int beamtype, gpuUserComplex *d_g1xs, gpuUserComplex *d_D1xs, gpuUserComplex *d_D1ys, gpuUserComplex *d_g1ys, gpuUserComplex *d_recov_g1x, gpuUserComplex *d_recov_D1x, gpuUserComplex *d_recov_D1y, gpuUserComplex *d_recov_g1y, gpuUserComplex *d_recov_g2x, gpuUserComplex *d_recov_D2x, gpuUserComplex *d_recov_D2y, gpuUserComplex *d_recov_g2y, int use_twobeams, int num_ants, int *d_ant1_to_baseline_map, int *d_ant2_to_baseline_map)¶
Wrapper kernel to test
get_beam_gains_gpuandget_beam_gains_multibeams_gpuin unit tests.
-
__global__ void kern_update_sum_visis_stokesIQUV(int num_freqs, int num_baselines, int num_components, int num_times, int beamtype, int off_cardinal_dipoles, gpuUserComplex *d_g1xs, gpuUserComplex *d_D1xs, gpuUserComplex *d_D1ys, gpuUserComplex *d_g1ys, int *d_ant1_to_baseline_map, int *d_ant2_to_baseline_map, int use_twobeams, gpuUserComplex *d_visi_components, user_precision_t *d_flux_I, user_precision_t *d_flux_Q, user_precision_t *d_flux_U, user_precision_t *d_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)¶
Wrapper kernel to test
update_sum_visis_stokesIQUV_gpuin unit tests.