W-tower gridder kernel functions

C/C++

struct sdp_GridderWtowerUVW
#include “ska-sdp-func/grid_data/sdp_gridder_wtower_uvw.h”

Uses 3D w-towers implementation and PSWF for subgrid (de)gridding.

These functions are ported from the implementation provided in Peter Wortmann’s subgrid_imaging Jupyter notebook.

sdp_GridderWtowerUVW *sdp_gridder_wtower_uvw_create(sdp_MemType type, sdp_GridProjection proj, int subgrid_size, int support, int oversampling, int w_support, int w_oversampling, sdp_Error *status)

Create plan for w-towers (de)gridder.

Parameters:
  • type – Number type to use internally (either complex float or double)

  • proj – Image grid projection to use

  • subgrid_size – Sub-grid size in pixels (quadratic).

  • w_step – Spacing between w-planes.

  • support – Kernel support size in (u, v).

  • oversampling – Oversampling factor for uv-kernel.

  • w_support – Support size in w.

  • w_oversampling – Oversampling factor for w-kernel.

  • status – Error status.

Returns:

Handle to gridder plan.

int64_t sdp_gridder_wtower_uvw_degrid(sdp_GridderWtowerUVW *plan, const sdp_Mem *subgrid_image, int subgrid_offset_u, int subgrid_offset_v, int subgrid_offset_w, double freq0_hz, double dfreq_hz, const sdp_Mem *uvws, const sdp_Mem *start_chs, const sdp_Mem *end_chs, const sdp_Mem *conjugate, sdp_Mem *vis, int64_t start_row, int64_t end_row, sdp_Error *status)

Degrid visibilities using w-stacking/towers.

The caller must ensure the output visibility array is sized correctly.

Parameters:
  • plan – Handle to gridder plan.

  • subgrid_image – Fourier transformed subgrid to degrid from. Note that the subgrid could especially span the entire grid, in which case this could simply be the entire (corrected) image.

  • subgrid_offset_u – Offset of subgrid centre relative to grid centre.

  • subgrid_offset_v – Offset of subgrid centre relative to grid centre.

  • subgrid_offset_w – Offset of subgrid centre relative to grid centre.

  • freq0_hz – Frequency of first channel (Hz).

  • dfreq_hz – Channel separation (Hz).

  • uvwsfloat[uvw_count, 3] UVW coordinates of visibilities (in m).

  • start_chsint[uvw_count] First channel to degrid for every uvw.

  • end_chsint[uvw_count] Channel at which to stop degridding for every uvw.

  • conjugatebool[uvw_count] Whether to conjugate visibilities (optional)

  • viscomplex[uvw_count, ch_count] Output degridded visibilities.

  • start_row – Row (uvw index) at which to start processing data.

  • end_row – Row (uvw index) at which to stop processing data (exclusive).

  • status – Error status.

Returns:

Number of visibilites degridded.

void sdp_gridder_wtower_uvw_degrid_correct(sdp_GridderWtowerUVW *plan, sdp_Mem *facet, int facet_offset_l, int facet_offset_m, int w_offset, sdp_Error *status)

Do degrid correction to enable degridding from the FT of the image.

Parameters:
  • plan – Handle to gridder plan.

  • facetcomplex[facet_size_l, facet_size_m] Facet.

  • facet_offset_l – Offset of facet centre in l, relative to image centre.

  • facet_offset_m – Offset of facet centre in m, relative to image centre.

  • w_offset – Offset in w, to allow for w-stacking.

  • status – Error status.

int64_t sdp_gridder_wtower_uvw_grid(sdp_GridderWtowerUVW *plan, const sdp_Mem *vis, const sdp_Mem *uvws, const sdp_Mem *start_chs, const sdp_Mem *end_chs, const sdp_Mem *conjugate, double freq0_hz, double dfreq_hz, sdp_Mem *subgrid_image, int subgrid_offset_u, int subgrid_offset_v, int subgrid_offset_w, int64_t start_row, int64_t end_row, sdp_Error *status)

Grid visibilities using w-stacking/towers.

The caller must ensure the output subgrid image is sized correctly.

Parameters:
  • plan – Handle to gridder plan.

  • viscomplex[uvw_count, ch_count] Input visibilities.

  • uvwsfloat[uvw_count, 3] UVW coordinates of visibilities (in m).

  • start_chsint[uvw_count] First channel to degrid for every uvw.

  • end_chsint[uvw_count] Channel at which to stop degridding for every uvw.

  • conjugatechar[uvw_count] Whether to interpret visibilities as conjugated.

  • freq0_hz – Frequency of first channel (Hz).

  • dfreq_hz – Channel separation (Hz).

  • subgrid_image – Fourier transformed subgrid to be gridded to. Note that the subgrid could especially span the entire grid, in which case this could simply be the entire (corrected) image.

  • subgrid_offset_u – Offset of subgrid centre relative to grid centre.

  • subgrid_offset_v – Offset of subgrid centre relative to grid centre.

  • subgrid_offset_w – Offset of subgrid centre relative to grid centre.

  • start_row – Row (uvw index) at which to start processing data.

  • end_row – Row (uvw index) at which to stop processing data (exclusive).

  • status – Error status.

Returns:

Number of visibilities gridded.

void sdp_gridder_wtower_uvw_grid_correct(sdp_GridderWtowerUVW *plan, sdp_Mem *facet, int facet_offset_l, int facet_offset_m, int w_offset, sdp_Error *status)

Do grid correction after gridding.

Parameters:
  • plan – Handle to gridder plan.

  • facetcomplex[facet_size_l, facet_size_m] Facet.

  • facet_offset_l – Offset of facet centre in l, relative to image centre.

  • facet_offset_m – Offset of facet centre in m, relative to image centre.

  • w_offset – Offset in w, to allow for w-stacking.

  • status – Error status.

void sdp_gridder_wtower_uvw_free(sdp_GridderWtowerUVW *plan)

Destroy plan for w-towers (de)gridder.

Parameters:

plan – Handle to gridder plan.

int sdp_gridder_wtower_uvw_num_w_planes(const sdp_GridderWtowerUVW *plan, int gridding)

Report total number of w-planes processed in the sub-grid stack.

Parameters:
  • plan – Handle to gridder plan.

  • gridding – 0 for w-planes in degridding, 1 for w-planes in gridding.

int sdp_gridder_wtower_uvw_image_size_l(const sdp_GridderWtowerUVW *plan)

Accessor function to return the image size along the l-axis.

Parameters:

plan – Handle to gridder plan.

Returns:

The configured image size along the l-axis.

int sdp_gridder_wtower_uvw_image_size_m(const sdp_GridderWtowerUVW *plan)

Accessor function to return the image size along the m-axis.

Parameters:

plan – Handle to gridder plan.

Returns:

The configured image size along the m-axis.

int sdp_gridder_wtower_uvw_oversampling(const sdp_GridderWtowerUVW *plan)

Accessor function to return the (u, v)-oversampling factor.

Parameters:

plan – Handle to gridder plan.

Returns:

The configured (u,v)-oversampling factor.

double sdp_gridder_wtower_uvw_shear_u(const sdp_GridderWtowerUVW *plan)

Accessor function to return the shear factor in the u-dimension.

Parameters:

plan – Handle to gridder plan.

Returns:

The configured shear factor in the u-dimension.

double sdp_gridder_wtower_uvw_shear_v(const sdp_GridderWtowerUVW *plan)

Accessor function to return the shear factor in the v-dimension.

Parameters:

plan – Handle to gridder plan.

Returns:

The configured shear factor in the v-dimension.

int sdp_gridder_wtower_uvw_subgrid_size(const sdp_GridderWtowerUVW *plan)

Accessor function to return the sub-grid size.

Parameters:

plan – Handle to gridder plan.

Returns:

The configured sub-grid size.

int sdp_gridder_wtower_uvw_support(const sdp_GridderWtowerUVW *plan)

Accessor function to return the kernel support size in (u, v).

Parameters:

plan – Handle to gridder plan.

Returns:

The configured kernel support size in (u, v).

double sdp_gridder_wtower_uvw_theta_l(const sdp_GridderWtowerUVW *plan)

Accessor function to return the padded field of view along the l-axis.

Parameters:

plan – Handle to gridder plan.

Returns:

The configured padded field of view along the l-axis, in direction cosines.

double sdp_gridder_wtower_uvw_theta_m(const sdp_GridderWtowerUVW *plan)

Accessor function to return the padded field of view along the m-axis.

Parameters:

plan – Handle to gridder plan.

Returns:

The configured padded field of view along the m-axis, in direction cosines.

double sdp_gridder_wtower_uvw_shift_l(const sdp_GridderWtowerUVW *plan)

Accessor function to return the image centre position on the l-axis.

Parameters:

plan – Handle to gridder plan.

Returns:

Image centre on the l-axis, in direction cosines.

double sdp_gridder_wtower_uvw_shift_m(const sdp_GridderWtowerUVW *plan)

Accessor function to return the image centre position on the m-axis.

Parameters:

plan – Handle to gridder plan.

Returns:

Image centre on the m-axis, in direction cosines.

double sdp_gridder_wtower_uvw_shift_n(const sdp_GridderWtowerUVW *plan)

Accessor function to return the phase centre position on the n-axis.

Parameters:

plan – Handle to gridder plan.

Returns:

Phase centre on the n-axis, in direction cosines.

int sdp_gridder_wtower_uvw_w_oversampling(const sdp_GridderWtowerUVW *plan)

Accessor function to return the w-oversampling factor.

Parameters:

plan – Handle to gridder plan.

Returns:

The configured w-oversampling factor.

double sdp_gridder_wtower_uvw_w_step(const sdp_GridderWtowerUVW *plan)

Accessor function to return the distance between w-layers.

Parameters:

plan – Handle to gridder plan.

Returns:

The configured distance between w-layers.

int sdp_gridder_wtower_uvw_w_support(const sdp_GridderWtowerUVW *plan)

Accessor function to return the kernel support size in w.

Parameters:

plan – Handle to gridder plan.

Returns:

The configured kernel support size in w.

double sdp_gridder_wtower_uvw_recentre_l(const sdp_GridderWtowerUVW *plan)

Accessor function to return recentring factor l.

Parameters:

plan – Handle to gridder plan.

Returns:

The configured recentring factor

double sdp_gridder_wtower_uvw_recentre_m(const sdp_GridderWtowerUVW *plan)

Accessor function to return recentring factor m.

Parameters:

plan – Handle to gridder plan.

Returns:

The configured recentring factor

Python

class ska_sdp_func.grid_data.GridderWtowerUVW(proj: GridProjection, subgrid_size: int, support: int, oversampling: int, w_support: int, w_oversampling: int, use_float: bool = False)

Uses w-towers / w-stacking for uvw (de)gridding.

These functions are ported from the implementation provided in Peter Wortmann’s subgrid_imaging Jupyter notebook.

__init__(proj: GridProjection, subgrid_size: int, support: int, oversampling: int, w_support: int, w_oversampling: int, use_float: bool = False)

Create plan for w-towers (de)gridder.

Parameters:
  • proj – Image grid projection to use

  • image_size_m – Total image size in pixels on m axis.

  • subgrid_size – Sub-grid size in pixels.

  • theta_l – Total image size on l axis in direction cosines.

  • theta_m – Total image size on m axis in direction cosines.

  • shear_u – Shear parameter in u (use zero for no shear).

  • shear_v – Shear parameter in v (use zero for no shear).

  • support – Kernel support size in (u, v).

  • oversampling – Oversampling factor for uv-kernel.

  • w_support – Support size in w.

  • w_oversampling – Oversampling factor for w-kernel.

  • use_float – Use single precision internally (double otherwise).

degrid(subgrid_image: ndarray, subgrid_offset_u: int, subgrid_offset_v: int, subgrid_offset_w: int, freq0_hz: float, dfreq_hz: float, uvws: ndarray, start_chs: ndarray, end_chs: ndarray, vis: ndarray, conjugate: ndarray = None, start_row: int = -1, end_row: int = -1)

Deprecated since version 1.2.0: Use degrid_subgrid() instead.

Degrid visibilities using w-stacking/towers.

The caller must ensure the output visibility array is sized correctly.

Parameters:
  • subgrid_image – Fourier transformed subgrid to degrid from. Note that the subgrid could especially span the entire grid, in which case this could simply be the entire (corrected) image.

  • subgrid_offset_u – Offset of subgrid centre relative to grid centre.

  • subgrid_offset_v – Offset of subgrid centre relative to grid centre.

  • subgrid_offset_w – Offset of subgrid centre relative to grid centre.

  • freq0_hz – Frequency of first channel (Hz)

  • dfreq_hz – Channel width (Hz)

  • uvwsfloat[uvw_count, 3] UVW coordinates of visibilities (in m)

  • start_chsint[uvw_count] First channel to degrid for every uvw

  • end_chsint[uvw_count] Channel at which to stop degridding for every uvw

  • viscomplex[uvw_count, ch_count] Output visibilities

  • conjugatebool[uvw_count] Whether to conjugate visibilities (optional)

  • start_row – Row (uvw index) at which to start processing data.

  • end_row – Row (uvw index) at which to stop processing data (exclusive).

degrid_subgrid(subgrid_image: ndarray, subgrid_offset: tuple[int, int, int], ch_count: int, freq0_hz: float, dfreq_hz: float, uvws: ndarray, start_chs: ndarray, end_chs: ndarray, vis: ndarray = None, conjugate: ndarray = None, start_row: int = -1, end_row: int = -1)

Degrid visibilities using w-stacking/towers.

This matches the interface in the notebook version.

If supplied, the caller must ensure the output visibility array is sized correctly; otherwise, it will be created and returned.

Parameters:
  • subgrid_image – Fourier transformed subgrid to degrid from. Note that the subgrid could especially span the entire grid, in which case this could simply be the entire (corrected) image.

  • subgrid_offset – Tuple of integers containing offset of subgrid in (u, v, w) relative to grid centre.

  • ch_count – Channel count (determines size of array returned)

  • freq0_hz – Frequency of first channel (Hz)

  • dfreq_hz – Channel width (Hz)

  • uvwsfloat[uvw_count, 3] UVW coordinates of visibilities (in m)

  • start_chsint[uvw_count] First channel to degrid for every uvw

  • end_chsint[uvw_count] Channel at which to stop degridding for every uvw

  • viscomplex[uvw_count, ch_count] Output visibilities

  • conjugatebool[uvw_count] Whether to conjugate visibilities (optional)

  • start_row – Row (uvw index) at which to start processing data.

  • end_row – Row (uvw index) at which to stop processing data (exclusive).

degrid_correct(facet: ndarray, facet_offset_l: int, facet_offset_m: int, w_offset: int = 0)

Do degrid correction to enable degridding from the FT of the image.

Parameters:
  • facetcomplex[facet_size,facet_size] Facet.

  • facet_offset_l – Offset of facet centre relative to image centre.

  • facet_offset_m – Offset of facet centre relative to image centre.

  • w_offset – Offset in w, to allow for w-stacking.

grid(vis: ndarray, uvw: ndarray, start_chs: ndarray, end_chs: ndarray, freq0_hz: float, dfreq_hz: float, subgrid_image: ndarray, subgrid_offset_u: int, subgrid_offset_v: int, subgrid_offset_w: int, conjugate: ndarray = None, start_row: int = -1, end_row: int = -1)

Deprecated since version 1.2.0: Use grid_subgrid() instead.

Grid visibilities using w-stacking/towers.

The caller must ensure the output subgrid_image is sized correctly.

Parameters:
  • viscomplex[uvw_count, ch_count] Input visibilities

  • uvwfloat[uvw_count, 3] UVW coordinates of visibilities (in m)

  • start_chsint[uvw_count] First channel to grid for every uvw

  • end_chsint[uvw_count] Channel at which to stop gridding for every uvw

  • freq0_hz – Frequency of first channel (Hz)

  • dfreq_hz – Channel width (Hz)

  • subgrid_image – Fourier transformed subgrid to be gridded to. Note that the subgrid could especially span the entire grid, in which case this could simply be the entire (corrected) image.

  • subgrid_offset_u – Offset of subgrid centre relative to grid centre.

  • subgrid_offset_v – Offset of subgrid centre relative to grid centre.

  • subgrid_offset_w – Offset of subgrid centre relative to grid centre.

  • conjugatechar[uvw_count] Whether to interpret visibilities as conjugated (optional)

  • start_row – Row (uvw index) at which to start processing data.

  • end_row – Row (uvw index) at which to stop processing data (exclusive).

grid_subgrid(vis: ndarray, uvw: ndarray, start_chs: ndarray, end_chs: ndarray, ch_count: int, freq0_hz: float, dfreq_hz: float, subgrid_image: ndarray, subgrid_offset: tuple[int, int, int], conjugate: ndarray = None, start_row: int = -1, end_row: int = -1)

Grid visibilities using w-stacking/towers.

This matches the interface in the notebook version.

The caller must ensure the output subgrid_image is sized correctly.

Parameters:
  • viscomplex[uvw_count, ch_count] Input visibilities

  • uvwfloat[uvw_count, 3] UVW coordinates of visibilities (in m)

  • start_chsint[uvw_count] First channel to grid for every uvw

  • end_chsint[uvw_count] Channel at which to stop gridding for every uvw

  • ch_count – Channel count

  • freq0_hz – Frequency of first channel (Hz)

  • dfreq_hz – Channel width (Hz)

  • subgrid_image – Fourier transformed subgrid to be gridded to. Note that the subgrid could especially span the entire grid, in which case this could simply be the entire (corrected) image.

  • subgrid_offset – Tuple of integers containing offset of subgrid in (u, v, w) relative to grid centre.

  • conjugatechar[uvw_count] Whether to interpret visibilities as conjugated (optional)

  • start_row – Row (uvw index) at which to start processing data.

  • end_row – Row (uvw index) at which to stop processing data (exclusive).

grid_correct(facet: ndarray, facet_offset_l: int, facet_offset_m: int, w_offset: int = 0)

Do grid correction after gridding.

Parameters:
  • facetcomplex[facet_size,facet_size] Facet.

  • facet_offset_l – Offset of facet centre relative to image centre.

  • facet_offset_m – Offset of facet centre relative to image centre.

  • w_offset – Offset in w, to allow for w-stacking.

property image_size_l

Returns the image size on l axis.

property image_size_m

Returns the image size on m axis.

property oversampling

Returns the (u,v)-oversampling factor.

property shear_u

Returns the shear factor in the u-dimension.

property shear_v

Returns the shear factor in the v-dimension.

property subgrid_size

Returns the sub-grid size.

property support

Returns the kernel support size in (u, v).

property theta_l

Returns the padded field of view on l axis, in direction cosines.

property theta_m

Returns the padded field of view on m axis, in direction cosines.

property shift_l

Returns the l coordinate of the image centre.

property shift_m

Returns the m coordinate of the image centre.

property shift_n

Returns the n shift applied to phase centre in gridder.

property recentre_l

Returns the l recentering factor

property recentre_m

Returns the m recentering factor

property grid_projection

Returns the image grid projection.

property w_oversampling

Returns the w-oversampling factor.

property w_step

Returns the distance between w-layers.

property w_support

Returns the kernel support size in w.