pymaster.utils module
- class pymaster.utils.NmtParams[source]
Bases:
objectClass that holds the values of several parameters controlling the default behaviour of different NaMaster methods. Currently these are:
sht_calculator: the software package to use to calculate spherical harmonic transforms. Defaults to ducc, with the other choice being healpy.n_iter_default: number of iterations to use when computing map2alm spherical harmonic transforms of most maps. Defaults to 3.n_iter_mask_default: number of iterations to use when computing map2alm spherical harmonic transforms of masks. Defaults to 3.tol_pinv_default: relative eigenvalue threshold to use when computing matrix pseudo-inverses.
All variables can be changed using the
set_methods described below, and their current values can be checked withget_default_params(). Note that, except forsht_calculator, all of these variables can be tweaked when calling various NaMaster functions. The values stored in this object only hold the values they default to if they are not set in those function calls.
- pymaster.utils.set_sht_calculator(calc_name)[source]
Select default spherical harmonic transform calculator.
- Parameters:
calc_name (
str) – Calculator name. Allowed options are'ducc'or'healpy'.
- pymaster.utils.set_n_iter_default(n_iter, mask=False)[source]
Select the number of Jacobi iterations to use when computing map2alm spherical harmonic transforms.
- pymaster.utils.set_tol_pinv_default(tol_pinv)[source]
Select the relative eigenvalue threshold to use when computing matrix pseudo-inverses. Check the docstring of
moore_penrose_pinvh().- Parameters:
tol_pinv (
float) – Relative eigenvalue threshold.
- pymaster.utils.get_default_params()[source]
Returns a dictionary with the current values of all default parameters.
- class pymaster.utils.NmtMapInfo(wcs, axes)[source]
Bases:
objectThis object contains information about the pixelization of a curved-sky map.
NmtMapInfoobjects can be compared with one another using==and!=to check for fields with compatible pixelizations.- Parameters:
wcs (WCS) – a WCS object (see the astropy documentation). If
None, HEALPix pixelization is assumed, andaxesshould be a 1-element sequence with the number of pixels of the map.axes (array) – shape of the maps (length-2 for CAR maps, length-1 for HEALPix).
- reform_map(maps)[source]
Modifies a map to make it compatible with the standards used by NaMaster for map manipulation (e.g. spherical harmonic transforms). This includes flattening 2D maps, and flipping their coordinate axes if required by their associated WCS information. HEALPix maps are unmodified.
- Parameters:
maps (array) – 2D (for HEALPix) or 3D (for CAR) array containing a set of maps.
- Returns:
Reformed map.
- Return type:
(array)
- get_lmax()[source]
Returns the default maximum multipole \(\ell_{\rm max}\) associated with this pixelization scheme. This will be \(3 N_{\rm side}-1\) for HEALPix or \(\pi/{\rm min}(\Delta\theta,\Delta\varphi)\) for CAR maps (with \(\Delta\theta\) and \(\Delta\varphi\) the constant intervals of colatitude and azimuth in radians).
- class pymaster.utils.NmtAlmInfo(lmax)[source]
Bases:
objectObject containing information useful to manipulate sets of spherical harmonic coefficients \(a_{\ell m}\).
- Parameters:
lmax (
int) – Maximum multipole \(\ell_{\rm max}\) to which the \(a_{\ell m}\) s are calculated.
- pymaster.utils.mask_apodization(mask_in, aposize, apotype='C1')[source]
Apodizes a mask with an given apodization scale using different methods. A given pixel is determined to be “masked” if its value is 0. This method only works for HEALPix maps. Three apodization methods are currently implemented:
“C1” apodization: all pixels are multiplied by a factor \(f\) given by
\[\begin{split}f=\left\{ \begin{array}{cc} x-\sin(2\pi x)/(2\pi) & x<1\\ 1 & {\rm otherwise} \end{array} \right..\end{split}\]where \(x=\sqrt{(1-\cos\theta)/(1-\cos\theta_*)}\), with \(\theta_*\) the apodization scale, and \(\theta\) the separation between a pixel and its nearest masked pixel (i.e. where the mask takes a zero value).
“C2” apodization: similar to “C1”, but using the apodization function:
\[\begin{split}f=\left\{ \begin{array}{cc} \frac{1}{2}\left[1-\cos(\pi x)\right] & x<1\\ 1 & {\rm otherwise} \end{array} \right..\end{split}\]“Smooth” apodization: this is carried out in three steps:
All pixels within a disk of radius \(2.5\theta_*\) of a masked pixel are masked.
The resulting map is smoothed with a Gaussian window function with standard deviation \(\theta_*\).
The resulting map is multiplied by the original mask to ensure that all pixels that were previously masked are still masked.
Note
Note that, confusingly, the definition of the “C1” and “C2” apodization windows above is the opposite of the similar \(C^1\) and \(C^2\) functions in Grain et al. 2009. This is due to a typo in the early development of NaMaster, which may be too disruptive to fix at this stage. This may be modified in a future release. Until then, users should rely on the documentation above, and not the defitions in Grain et al. 2009.
- pymaster.utils.mask_apodization_flat(mask_in, lx, ly, aposize, apotype='C1')[source]
Apodizes a flat-sky mask. See the docstrings of
mask_apodization()for a description of the different methods implemented.- Parameters:
- Returns:
Apodized mask.
- Return type:
(array)
- pymaster.utils.synfast_spherical(nside, cls, spin_arr, beam=None, seed=-1, wcs=None, lmax=None)[source]
Generates full-sky Gaussian random fields according to a given power spectrum. This function should produce outputs similar to healpy’s synfast.
- Parameters:
nside (
int) – HEALpix resolution parameter. If you want rectangular pixel maps, ignore this parameter and pass a WCS object aswcs(see below).cls (array) – Array contiaining power spectra. Shape should be
(n_cls, n_ell), wheren_clsis the number of power spectra needed to define all the fields. This should ben_cls = n_maps * (n_maps + 1) / 2, wheren_mapsis the total number of maps required (1 for each spin-0 field, 2 for each spin>0 field). Power spectra must be provided only for the upper-triangular part of the full power spectrum matrix in row-major order (e.g. ifn_maps=3, there will be 6 power spectra ordered as [1x1, 1x2, 1x3, 2x2, 2x3, 3x3]).spin_arr (array) – Array containing the spins of all the fields to generate.
beam (array) – 2D array containing the instrumental beam of each field to simulate (the output map(s) will be convolved with it).
seed (
int) – Seed for the random number generator. If negative, a random seed will be used.wcs (WCS) – A WCS object if using rectangular pixels (see the astropy documentation).
lmax (
int) – Maximum multipole up to which the spherical harmonic coefficients of the maps will be generated.
- Returns:
A set of full-sky maps (1 for each spin-0 field, 2 for each spin-s field).
- Return type:
(array)
- pymaster.utils.synfast_flat(nx, ny, lx, ly, cls, spin_arr, beam=None, seed=-1)[source]
Generates flat-sky Gaussian random fields according to a given power spectrum. This is the flat-sky equivalent of
synfast_spherical().- Parameters:
nx (
int) – Number of pixels in the x axis.ny (
int) – Number of pixels in the y axis.lx (
float) – Patch size in the x axis (radians).ly (
float) – Patch size in the y axis (radians).cls (array) – Array contiaining power spectra. Shape should be
(n_cls, n_ell), wheren_clsis the number of power spectra needed to define all the fields. This should ben_cls = n_maps * (n_maps + 1) / 2, wheren_mapsis the total number of maps required (1 for each spin-0 field, 2 for each spin>0 field). Power spectra must be provided only for the upper-triangular part of the full power spectrum matrix in row-major order (e.g. ifn_maps=3, there will be 6 power spectra ordered as [1x1, 1x2, 1x3, 2x2, 2x3, 3x3]). The power spectra are assumed to be sampled at all integer multipoles from \(\ell=0\) ton_ell-1.spin_arr (array) – Array containing the spins of all the fields to generate.
beam (array) – 2D array containing the instrumental beam of each field to simulate (the output map(s) will be convolved with it).
seed (
int) – Seed for the random number generator. If negative, a random seed will be used.
- Returns:
An array of flat-sky maps (1 for each spin-0 field, 2 for each spin-s field) with shape
(nmaps, ny, nx).- Return type:
(array)
- pymaster.utils.moore_penrose_pinvh(mat, tol_pinv)[source]
Compute the Moore-Penrose pseudo-inverse of a Hermitian matrix. This is done by diagonalising the matrix, setting the inverse of all eigenvales below a given threshold to zero, and reconstructing the inverse matrix from the modified eigenvalues. The inverse of all eigenvalues smaller than a factor
tol_pinvtimes the largest eigenvalue will be set to zero. Iftol_pinv <= 0, the standard inverse is computed.- Parameters:
mat (array) – 2D Hermitian matrix.
tol_pinv (
float) – Relative eigenvalue threshold.
- Returns:
Pseudo-inverse of input matrix.
- Return type:
(array)
- pymaster.utils.map2alm(map, spin, map_info, alm_info, *, n_iter)[source]
Computes the spherical harmonic transform (SHT) for a set of input maps. The SHT implementation to be used can be selected via
set_sht_calculator()(see the documentation ofNmtParams).- Parameters:
map (array) – 2D array with shape
(nmaps, npix)wherenmapsis either 1 (for spin-0 fields) or 2 (for spin-s fields), andnpixis the number of pixels. If using CAR rectangular pixels, the map should be provided flattened.spin (
int) – Field spin.map_info (
NmtMapInfo) – Object describing the pixelization of the input map.alm_info (
NmtAlmInfo) – Object describing the structure of the output spherical harmonic coefficients.n_iter (
int) – Number of Jacobi iterations used to improve the accuracy of the SHT.
- Returns:
Harmonic coefficients \(a_{\ell m}\) of the input map. A set of two arrays (E and B modes) is returned if
spin>0.- Return type:
(array)
- pymaster.utils.alm2map(alm, spin, map_info, alm_info)[source]
Computes the inverse spherical harmonic transform (SHT) for a set of input \(a_{\ell m}\) s. The SHT implementation to be used can be selected via
set_sht_calculator()(see the documentation ofNmtParams).- Parameters:
alm (array) – 2D array with shape
(nmaps, n_lm)wherenmapsis either 1 (for spin-0 fields) or 2 (for spin-s fields), andn_lmis the number of harmonic coefficients.spin (
int) – Field spin.map_info (
NmtMapInfo) – Object describing the pixelization of the output map.alm_info (
NmtAlmInfo) – Object describing the structure of the input spherical harmonic coefficients.
- Returns:
Map reconstructed from the input \(a_{\ell m}\) s. A set of two arrays (e.g. Q and U Stokes parameters) is returned if
spin>0.- Return type:
(array)