source_components

API documentation for source_components.cpp

Device methods to extrapolate flux densities, assign and apply primary beam gains, and visibility responses for different sky model COMPONENT types.

Typedefs

typedef struct _d_beam_gains_t d_beam_gains_t

A struct to contain primary beam values for a give COMPONENT. d_gxs,d_Dxs,d_Dys,d_gys should be used when all antennas have the same primary beam, and d_gxs,d_Dxs,d_Dys,d_gys used when all primary beams are different.

Functions

__device__ gpuUserComplex 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, a gpuUserComplex of the visibility

__device__ void apply_beam_gains_stokesIQUV_on_cardinal(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, gy and antenna 2 g1x, D1x, D1y, gy,the complex visibility phase across those two antennas visi, 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(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, gy and antenna 2 g1x, D1x, D1y, gy,the complex visibility phase across those two antennas visi, 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(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, gy and antenna 2 g1x, D1x, D1y, gy,the complex visibility phase across those two antennas visi, 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(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, gy and antenna 2 g1x, D1x, D1y, gy,the complex visibility phase across those two antennas visi, 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(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 arrays d_primay_beam_J* that match the indexes iBaseline, 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_J11 should contain the primary beam settings for all times, all frequencies, and COPMONENTs on the sky, and so should be num_freqs*num_components*num_times long. The order elements should increment through COMPONENT (fastest changing), freqeuncy, and time (slowest changing). iBaseline is a combined index of baseline, frequency, and time.

If beamtype == NO_BEAM, set g1x = g1y = g2x = g2y = 1 and D1x = 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(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 arrays d_primay_beam_J* that match the indexes iBaseline, iComponent. This function assumes the primary beams are different for every antenna, so returns different values for each antenna.

Todo:

Currently this is only set up to work with the MWA_FEE and MWA_FEE_INTERP primary beams

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_J11 should contain the primary beam settings for all antennas, all times, all frequencies, and COPMONENTs on the sky, and so should be num_ants*num_freqs*num_components*num_times long. The order elements should increment through COMPONENT (fastest changing), freqeuncy, time, and antenna (slowest changing). iBaseline is a combined index of baseline, frequency, and time.

If beamtype == NO_BEAM, set g1x = g1y = g2x = g2y = 1 and D1x = 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(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 by iBaseline and COMPONENT index iComponent, 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 arrays d_sim_visi_*_real and d_sim_visi_*_imag.

Uses get_beam_gains and apply_beam_gains as 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

  • off_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(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 by iBaseline and COMPONENT index iComponent, 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 arrays d_sim_visi_*_real and d_sim_visi_*_imag.

Uses get_beam_gains or get_beam_gains_multibeams, and apply_beam_gains as 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

  • off_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

void malloc_extrapolated_flux_arrays(components_t *d_components, int num_comps, int num_freqs)

Allocate device memory of extrapolated Stokes arrays int d_components

It does this if d_components.do_QUV == 1:

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) );
If do_QUV == 0, only allocate the StokesI array.

Parameters:
  • *d_components[inout] A populated components_t struct

  • num_comps[in] Number of components

  • num_freqs[in] Number of frequencies

__device__ void extrap_stokes_power_law(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], and spectral index d_SIs[iFluxComp] to calculate the flux density extrap_flux at 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 for all Stokes I components in d_components.

Fills the array d_components.extrap_stokesI with 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 for all Stokes V components in d_components.

Fills the array d_components.extrap_stokesV with 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 for all linear polarisation components in d_components.

Temporarily fills the array d_components.extrap_stokesQ with the extrapolated linear polarisation flux densities. The intention is to the use kern_apply_rotation_measure after the fact, which will use d_components.extrap_stokesQ as input flux, do rotations, and then fill d_components.extrap_stokesU and d_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_curved_power_law(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], curvature q param spectral index d_qs[iFluxComp] to calculate the flux density extrap_flux at 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 q param 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_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 for all Stokes I components in d_components.

Fills the array d_components.extrap_stokesI with 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 for all Stokes I components in d_components.

Fills the array d_components.extrap_stokesV with 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 for all linear polarisation components in d_components.

Temporarily fills the array d_components.extrap_stokesQ with the extrapolated linear polarisation flux densities. The intention is to the use kern_apply_rotation_measure after the fact, which will use d_components.extrap_stokesQ as input flux, do rotations, and then fill d_components.extrap_stokesU and d_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_stokesI filled already. Fills the array d_components.extrap_stokesV using d_components.extrap_stokesI and d_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_stokesI filled already. Temporarily fills the array d_components.extrap_stokesQ using d_components.extrap_stokesI and d_components.linpol_pol_fracs. The intention is to the use kern_apply_rotation_measure after the fact, which will use d_components.extrap_stokesQ as input flux, do rotations, and then fill d_components.extrap_stokesU and d_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(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_stokes and list_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.

Fills extrap_stokes with the extrapolated stokes flux densities. Runs the function extrap_stokes_list_fluxes.

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_stokes and list_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

void extrapolate_Stokes(source_t *d_chunked_source, double *d_extrap_freqs, int num_extrap_freqs, e_component_type comptype)

Extrapolate the flux densities of certain spectral model type in a source to a set of frequencies.

This function is a wrapper for the various extrapolation functions. It will extrapolate all spectral model types for a give component type as specified by comptype. It will extrapolate to all the frequencies specified in d_extrap_freqs, and store the results within d_chunked_source

Parameters:
  • d_chunked_source[inout] Pointer to the source data.

  • d_extrap_freqs[in] Pointer to an array of frequencies to extrapolate to.

  • num_extrap_freqs[in] Number of frequencies in the d_extrap_freqs array.

  • comptype[in] The type of component to extrapolate (e.g. POINT, GAUSSIAN, SHAPELET).

void source_component_common(woden_settings_t *woden_settings, beam_settings_t *beam_settings, double *d_freqs, source_t *chunked_source, source_t *d_chunked_source, d_beam_gains_t *d_component_beam_gains, e_component_type comptype, visibility_set_t *d_visibility_set)

Performs necessary calculations that are common to all POINT, GAUSSIAN, and SHAPELET component types: calculating the \(l,m,n\) coordinates; extrapolating flux densities to given frequencies; calculating primary beam values.

*chunked_source should be a populated source_t struct as filled by chunk_sky_model::create_chunked_sky_models, and then remapped by chunk_sky_model::remap_source_for_gpu. The remapping has a memory cost, hence should be done iteratively over chunks (this is done inside calculate_visibilities.cu).

This function uses fundamental_coords::kern_calc_lmn to calculate \(l,m,n\), and stores the ouputs in the components_t struct d_components. Uses the d_beam_gains_t struct d_component_beam_gains to store the calculated primary beam responses. Each gain and leakage will have num_components*woden_settings->num_time_steps*woden_settings->num_freqs elements of memory allocated in the process.

If the primary beam requested (set via beam_settings->beamtype ) is either GAUSS_BEAM or MWA_ANALY, both chunked_source->components->beam_has and chunked_source->components->beam_decs need to be set.

If woden_settings->do_autos is True, will use source_components::kern_calc_autos to calculate all the auto-correlations, and stores them in d_visibility_set.

Parameters:
  • woden_settings[in] Populated woden_settings_t struct

  • beam_settings[in] Populated beam_settings_t struct

  • d_freqs[in] Frequencies to calculate beam responses to (stored on the device)

  • *chunked_source[in] A populated source_t struct as filled by chunk_sky_model::create_chunked_sky_models and then remapped by chunk_sky_model::remap_source_for_gpu

  • *d_chunked_source[inout] A populated source_t struct with device memory as filled by source_components::copy_chunked_source_to_GPU

  • d_component_beam_gains[inout] Pointer to d_beam_gains_t struct in which to malloc and store results on the device in

  • comptype[in] Either POINT, GAUSSIAN, SHAPELET - selects which set of COMPONENTS stored in chunked_source and chunked_source to work with

  • d_visibility_set[in] If woden_settings->do_autos is True, use this to store auto correlations in. Needs to have device memory allocated. Otherwise, d_visibility_set is unused, so can just be an empty pointer

__global__ void kern_calc_visi_point_or_gauss(components_t d_components, d_beam_gains_t d_component_beam_gains, user_precision_t *d_us, user_precision_t *d_vs, user_precision_t *d_ws, user_precision_t *d_sum_visi_XX_real, user_precision_t *d_sum_visi_XX_imag, user_precision_t *d_sum_visi_XY_real, user_precision_t *d_sum_visi_XY_imag, user_precision_t *d_sum_visi_YX_real, user_precision_t *d_sum_visi_YX_imag, user_precision_t *d_sum_visi_YY_real, user_precision_t *d_sum_visi_YY_imag, 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_components of either POINT or GAUSSIAN COMPONENTs, and sum the outputs to d_sum_visi_*_real, d_sum_visi_*_imag.

Uses the functions extrap_stokes, calc_measurement_equation, update_sum_visis as 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 using comptype.

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 with grid.x defined, where:

  • grid.x * threads.x >= num_visi

Parameters:
  • d_components[in] d_components Pointer to a populated components_t struct as filled by source_components::source_component_common

  • d_component_beam_gains[in] Pointer to a populated d_beam_gains_t struct as filled by source_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_beamtype

  • comptype[in] Component type, either POINT or GAUSSIAN

  • off_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.

__global__ void kern_calc_visi_shapelets(components_t d_components, d_beam_gains_t d_component_beam_gains, user_precision_t *d_us, user_precision_t *d_vs, user_precision_t *d_ws, user_precision_t *d_allsteps_wavelengths, user_precision_t *d_u_shapes, user_precision_t *d_v_shapes, user_precision_t *d_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_shapes of SHAPELET COMPONENTs, and sum the outputs to d_sum_visi_*_real, d_sum_visi_*_imag.

Uses the functions extrap_stokes, calc_measurement_equation, update_sum_visis as detailed above to calculate the visibilities. Furthermore calculates the visibility envelope \(\mathrm{V}_{\mathrm{env}}\) to convert the basic visibility into a SHAPELET:

\[\begin{split}\begin{eqnarray*} \mathrm{V}_{\mathrm{env}} &=& \sum^{n_k +n_l < n_\mathrm{max}}_{k,l} C_{n_k,n_l} B_{n_k,n_l}(k_x,k_y)\\ k_x &=& \frac{\pi}{\sqrt{2\ln(2)}} \left[\cos(\phi_{PA})v_{\mathrm{comp}} + \sin(\phi_{PA})u_{\mathrm{comp}} \right] \\ k_y &=& \frac{\pi}{\sqrt{2\ln(2)}} \left[-\sin(\phi_{PA})v_{\mathrm{comp}} + \cos(\phi_{PA})u_{\mathrm{comp}} \right] \end{eqnarray*}\end{split}\]

where \( B_{n_k,n_l} \) is a stored 2D shapelet basis function from a look up table, of order \( n_k,n_l \) and scaled by the major and minor axis, \(C_{n_k,n_l}\) is a scaling coefficient, and \( u_{\mathrm{comp}}, v_{\mathrm{comp}} \) are visibility coordinates with the SHAPELET ra,dec as the phase centre. These should have been calculated using fundamental_coords::kern_calc_uvw_shapelet.

This kernel differs from kern_calc_visi_point_or_gauss in 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 array d_components.param_indexes is used to match the basis function information d_components.n1s, d_components.n2s, d_components.shape_coeffs to 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_sbf is the shapelet basis function lookup table, and should have been generated with shapelet_basis.create_sbf and copied into device memory.

When called with dim3 grid, threads, kernel should be called with grid.x defined, where:

  • grid.x * threads.x >= num_visi

Parameters:
  • d_components[in] d_components Pointer to a populated components_t struct as filled by source_components::source_component_common

  • d_component_beam_gains[in] Pointer to a populated d_beam_gains_t struct as filled by source_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_beamtype

  • off_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_source to the device memory pointed to by d_chunked_source, allocating the device memory in the process The comptype parameter 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_d_components(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, or d_chunked_source->shape_components, depending on comptype.

Parameters:
  • *d_chunked_source[inout] A source_t with device memory to be freed

  • comptype[in] Which type of COMPONENT, either POINT, GAUSSIAN, or SHAPELET

void free_beam_gains(d_beam_gains_t d_beam_gains, e_beamtype beamtype)

Frees the device memory on d_beam_gains, given the type of primary beam beamtype.

Only certain primary beam models have leakage terms, so need to know the beamtype to free the correct locations.

Parameters:
  • d_beam_gains[inout] A d_beam_gains_t with device memory to be freed

  • beamtype[in] An e_beamtype describing the primary beam type

source_t *copy_chunked_source_to_GPU(source_t *chunked_source)

Copies a populated source_t from 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)

Parameters:

*chunked_source[inout] A populated source_t struct

void free_extrapolated_flux_arrays(components_t *d_components)

Free device memory of extrapolated Stokes arrays from d_components

It 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.

Parameters:

*d_components[inout] A populated components_t struct

__global__ void kern_calc_autos(components_t d_components, d_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)

Calculate the auto-correlations for all antennas given the fluxes already calculated in d_components and beam gains in d_component_beam_gains. Stores the outputs at the end of d_sum_visi*, after the cross-correlations.

Currently, the primary beam for all antennas (or what the MWA calls tiles) are identical. So d_component_beam_gains should have gains for one primary beam. This function then just does the dot product of the component fluxes and beam gains to produce linear stokes polarisation auto-correlations.

When called with dim3 grid, threads, kernel should be called with grid.x and grid.y defined, where:

  • grid.x * threads.x >= num_freqs*num_times

  • grid.y * threads.y >= num_ants

Parameters:
  • d_components[in] d_components Pointer to a populated components_t struct as filled by source_components::source_component_common

  • d_component_beam_gains[in] Pointer to a populated d_beam_gains_t struct as filled by source_components::source_component_common

  • beamtype[in] Beam type see woden_struct_defs.e_beamtype

  • num_components[in] Number of components in d_components

  • num_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

  • *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

void fill_ant_to_baseline_mapping(int num_ants, int *d_ant1_to_baseline_map, int *d_ant2_to_baseline_map)

Fill the d_ant1_to_baseline_map and d_ant2_to_baseline_map device arrays with indexes corresponding to ant1 and ant2 for all unique baselines in an array of num_ants antennas.

The d_ant1_to_baseline_map and d_ant2_to_baseline_map should 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

struct _d_beam_gains_t
#include <source_components.h>

A struct to contain primary beam values for a give COMPONENT. d_gxs,d_Dxs,d_Dys,d_gys should be used when all antennas have the same primary beam, and d_gxs,d_Dxs,d_Dys,d_gys used when all primary beams are different.

Public Members

gpuUserComplex *d_gxs = NULL

Device copy of North-South Beam gain values for all beams, directions, frequencies, and times for these COMPONENTS

gpuUserComplex *d_Dxs = NULL

Device copy of North-South Beam leakage values for all beams, directions, frequencies, and times for these COMPONENTS

gpuUserComplex *d_Dys = NULL

Device copy of East-West Beam leakage values for all beams, directions, frequencies, and times for these COMPONENTS

gpuUserComplex *d_gys = NULL

Device copy of East-West Beam gain values for all beams, directions, frequencies, and times for these COMPONENTS

int *d_ant1_to_baseline_map = NULL

The index of antenna 1 in all unique pairs of antennas. Used to map iBaseline to the correct antenna 1

int *d_ant2_to_baseline_map = NULL

The index of antenna 2 in all unique pairs of antennas. Used to map iBaseline to the correct antenna 2

int use_twobeams

The beam gains were made with unique primary beams so should use two antenna patterns per visibility