Cygrid#

class cygrid.Cygrid(*args, **kwargs)#

Bases: object

Fast Cython-powered gridding software, base class.

DO NOT USE DIRECTLY, BUT THE DERIVED CLASSES WcsGrid OR SlGrid.

The underlying algorithm is a based on serialized convolution with finite gridding kernels. Currently, only Gaussian kernels are provided (which has the drawback of slight degradation of the effective resolution). The algorithm has very small memory footprint, allows easy parallelization, and is very fast.

Look into the set_kernel and grid methods help for more information on how to use this.

Internally, we make use of the HEALPix representation for book-keeping. The idea is the following: for each input point we query which HEALPix pixels are located within the required convolution kernel radius (using HEALPix query_disc function). Likewise, for the target map pixels (any WCS projection) we calculate the HEALPix index they live in. By a simple cross- matching (hash-map based) we can thus easily find out which input pixels contribute to which output pixels. In practice it is a little more complicated, because world pixels could share the same HEALPix index. We use lists (or rather C++ vectors) to account for this.

Parameters:
Optional keyword arguments:
dbg_messagesBoolean (default: False)

Do debugging output

Methods Summary

clear_cache(self)

Clear all internal caches/dictionaries.

clear_data_and_weights(self)

Reset data and weights arrays.

get_datacube(self)

Return final data.

get_unweighted_datacube(self)

Return final unweighted data.

get_weights(self)

Return final weights.

grid(self, ndarray lons, ndarray lats, ...)

Grid irregularly positioned data points (spectra) into the data cube.

set_kernel(self, kernel_type, kernel_params, ...)

Set the gridding kernel type and parameters.

set_num_threads(self, int nthreads)

Change maximum number of threads to use.

Methods Documentation

clear_cache(self)#

Clear all internal caches/dictionaries.

Notes

To speed-up processing, cygrid uses internal caches (see Cygrid paper for details). This is usually a minor contribution to the total memory usage, but there could be scenarios, e.g., when one grids into very large maps having very small pixel sizes, that the internal cache grows too much. In this case, we recommend to sort the input data by latitude, grid the data in chunks, and call this function every now and then.

clear_data_and_weights(self)#

Reset data and weights arrays. (Caches stay filled!)

get_datacube(self)#

Return final data.

Returns:
dataarray [nD] of float32 or float64

The gridded spectral data.

Notes

  • The actual shape of the returned array depends on the input raw data signal shape and the target pixel array size.

get_unweighted_datacube(self)#

Return final unweighted data. (For debugging only.)

Returns:
unweighted_dataarray [nD] of float32 or float64

The gridded spectral unweighted data. The gridded data is the quotient of the unweighted data and the weights.

Notes

  • The actual shape of the returned array depends on the input raw data signal shape and the target pixel array size.

get_weights(self)#

Return final weights.

Returns:
weightsarray [nD] of float32 or float64

The gridded spectral weights.

Notes

  • The actual shape of the returned array depends on the input raw data signal shape and the target pixel array size.

grid(self, ndarray lons, ndarray lats, ndarray data, ndarray weights=None)#

Grid irregularly positioned data points (spectra) into the data cube.

After successful gridding, you can obtain the resulting datacube with the get_datacube method. The associated weight cube is accessible with get_weights.

Parameters:
lons, latsarray [1D] of float64

Flat lists/arrays of input coordinates \((l, b)\).

data/weightsarray [nD] of float32 or float64

The spectra and their weights (optional) for each of the given coordinate pairs, \((l, b)\). First axis must match lons/lats size. The shape of the input data will determine the output datacube shape: the first axes of the datacube will match the ( raw) data (aka input signal) shape - only the first dimension is being stripped as it is associated with the number of coordinate pairs -, while the last axes match the shape of the target pixel array (e.g., the shape of the target map).

Raises:
ShapeError

Input coordinates/data points length mismatch. Number of spectral channels mismatch.

Notes

  • Internally, cygrid always works with a 3D representation of the data cube, i.e., with (possibly redundant) spectral axis. However, as the original shape is stored, the gridded data returned by get_datacube will always be consistent with the input raw data signal (in terms of dimensions).

  • All input parameters need to be C-contiguous (cygrid will re-cast if necessary).

set_kernel(self, kernel_type, kernel_params, double sphere_radius, double hpx_max_resolution)#

Set the gridding kernel type and parameters.

Parameters:
kernel_typestr

Set the kernel type.

The following names/types are available: ‘gauss1d’, ‘gauss2d’, ‘tapered_sinc’ (see Notes for details)

kernel_paramsarray

Set the kernel parameters for the chosen type (see Notes for details)

sphere_radiusdouble

Kernel sphere radius.

This is controls out to which distance the kernel is computed for. For Gaussian kernels, values much larger than \(3\ldots4 \sigma_\mathrm{k}\) do not make much sense.

hpx_max_resolutiondouble

Maximum acceptable HPX resolution (\(\sigma_\mathrm{k}^{[\mathrm{maj}]} / 2\) is a reasonable value).

Notes

Below you find a list of kernel-names and required parameters:

'gauss1d', (kernel_sigma,)
'gauss2d', (kernel_sigma_maj, kernel_sigma_min, PA)
'tapered_sinc', (kernel_sigma, param_a, param_b)

All numbers are in units of degrees, except for PA, param_a and param_b. PA (the position angle) is in units of radians (for efficiency). param_a and param_b should be 2.52 and 1.55, respectively, for optimal results!

The kernel size (sigma) defines the amount of “smoothness” applied to the data. If in doubt a good value is about 50% of the true/input angular resolution of the data (this will result in about 10% degradation of the final angular resolution.)

set_num_threads(self, int nthreads)#

Change maximum number of threads to use.

This is a convenience function, to call omp_set_num_threads(), which is otherwise not possible during runtime from Python.