SHT

ducc0.sht

Python interface for spherical harmonic transforms and manipulation of spherical harmonic coefficients.

Error conditions are reported by raising exceptions.

ducc0.sht.adjoint_analysis_2d(*, alm: numpy.ndarray, spin: int, lmax: int, geometry: str, ntheta: object = None, nphi: object = None, mmax: object = None, nthreads: int = 1, map: object = None, phi0: float = 0.0, mstart: object = None, lstride: int = 1) numpy.ndarray

Transforms one or two sets of spherical harmonic coefficients to 2D maps. This is the adjoint operation of analysis_2D.

Parameters:
  • alm (numpy.ndarray((ncomp, x), dtype=numpy.complex64 or numpy.complex128)) – the set of spherical harmonic coefficients. The second dimension must be large enough to accommodate all entries, which are stored according to the healpy convention.

  • map (numpy.ndarray((ncomp, ntheta, nphi), dtype=numpy.float of same accuracy as alm)) – storage for the output map. If not supplied, a new array is allocated.

  • ntheta (int > 0) – dimensions of the output map. If not supplied, map must be supplied. If supplied, and map is also supplied, must match with the array dimensions

  • nphi (int > 0) – dimensions of the output map. If not supplied, map must be supplied. If supplied, and map is also supplied, must match with the array dimensions

  • spin (int >= 0) – the spin to use for the transform. If spin==0, ncomp must be 1, otherwise 2

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive).

  • mmax (int >= 0 <= lmax) – the maximum m moment of the transform (inclusive).

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • geometry (one of "CC", "F1", "MW", "MWflip", "GL", "DH", "F2") –

    the distribution of rings over the theta range
    • CC: Clenshaw-Curtis, equidistant, first and last ring on poles

    • F1: Fejer’s first rule, equidistant, first and last ring half a ring width from the poles

    • MW: McEwen & Wiaux scheme, equidistant, first ring half a ring width from the north pole, last ring on the south pole

    • MWflip: flipped McEwen & Wiaux scheme, equidistant, first ring on the north pole, last ring half a ring width from the south pole

    • GL: Gauss-Legendre, non-equidistant

    • DH: Driscoll-Healy, equidistant, first ring on north pole, last ring one ring width from the south pole

    • F2: Fejer’s second rule, equidistant, first and last ring one ring width from the poles.

  • phi0 (float) – the azimuth (in radians) of the first pixel in each ring

  • nthreads (int >= 0) – the number of threads to use for the computation. If 0, use as many threads as there are hardware threads available on the system

Returns:

the computed map. If the map parameter was specified, this is identical with map.

Return type:

numpy.ndarray((ncomp, ntheta, nphi), dtype=numpy.float of same accuracy as alm)

Notes

For limits on lmax and mmax see the documentation of analysis_2d.

ducc0.sht.adjoint_synthesis(*, map: numpy.ndarray, theta: numpy.ndarray, lmax: int, mstart: object = None, nphi: numpy.ndarray, phi0: numpy.ndarray, ringstart: numpy.ndarray, spin: int, lstride: int = 1, pixstride: int = 1, nthreads: int = 1, alm: object = None, mmax: object = None, mode: str = 'STANDARD', theta_interpol: bool = False) numpy.ndarray

Transforms (sets of) one or two maps to spherical harmonic coefficients. This is the adjoint operation of synthesis.

Parameters:
  • alm (None or numpy.ndarray(([ntrans,] nalm, x), dtype=numpy.complex of same precision as map)) – the set of spherical harmonic coefficients. The last dimension must be large enough to accommodate all entries, which are stored according to the parameters lmax, ‘mmax`, mstart, and lstride. if None, a new suitable array is allocated

  • map (numpy.ndarray(([ntrans,] nmaps, x), dtype=numpy.float32 or numpy.float64) – The last dimension must be large enough to accommodate all pixels, which are stored according to the parameters nphi, ringstart, and pixstride.

  • theta (numpy.ndarray((ntheta,), dtype=numpy.float64)) – the colatitudes of the map rings

  • nphi (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – number of pixels in every ring

  • phi0 (numpy.ndarray((ntheta,), dtype=numpy.float64)) – azimuth (in radians) of the first pixel in every ring

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • ringstart (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – the index in the last dimension of map at which the first pixel of every ring is stored

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • pixstride (int) – the index stride in the last dimension of map between two subsequent pixels in a ring

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

  • spin (int >= 0) – the spin to use for the transform.

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive).

  • mmax (int >= 0 <= lmax) – the maximum m moment of the transform (inclusive).

  • mode (str) –

    the transform mode
    ”STANDARD”: standard transform
    ”GRAD_ONLY”: only valid for spin>0, curl a_lm are not computed
    ”DERIV1”: same as “GRAD_ONLY”, but spin is assumed to be 1 and a_lm are multiplied by sqrt(l*(l+1))

  • theta_interpol (bool) – if the input grid is irregularly spaced in theta, try to accelerate the transform by using an intermediate equidistant theta grid and a 1D NUFFT.

Returns:

the set(s) of spherical harmonic coefficients. If alm was supplied, this will be the same object If newly allocated, the smallest possible last dimension will be chosen.

Return type:

numpy.ndarray(([ntrans,] nalm, x), dtype=numpy.complex of same accuracy as map)

Notes

nmaps = 1 if spin == 0 else 2 nalm = 1 if spin == 0 else (2 if mode == “STANDARD” else 1)

ducc0.sht.adjoint_synthesis_2d(*, map: numpy.ndarray, spin: int, lmax: int, geometry: str, mmax: object = None, nthreads: int = 1, alm: object = None, mode: str = 'STANDARD', phi0: float = 0.0, mstart: object = None, lstride: int = 1) numpy.ndarray

Transforms one or two 2D maps to spherical harmonic coefficients. This is the adjoint operation of synthesis_2D.

Parameters:
  • alm (numpy.ndarray((nalm, x), dtype=numpy.complex64 or numpy.complex128)) – storage for the spherical harmonic coefficients. The second dimension must be large enough to accommodate all entries, which are stored according to the healpy convention. If not supplied, a new array is allocated.

  • map (numpy.ndarray((nmaps, ntheta, nphi), dtype=numpy.float of same accuracy as alm)) – The input map.

  • spin (int >= 0) – the spin to use for the transform.

  • lmax (int >= 0) – the maximum l (and m) moment of the transform (inclusive)

  • mmax (int >= 0 <= lmax) – the maximum m moment of the transform (inclusive).

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • geometry (one of "CC", "F1", "MW", "MWflip", "GL", "DH", "F2") –

    the distribution of rings over the theta range
    • CC: Clenshaw-Curtis, equidistant, first and last ring on poles

    • F1: Fejer’s first rule, equidistant, first and last ring half a ring width from the poles

    • MW: McEwen & Wiaux scheme, equidistant, first ring half a ring width from the north pole, last ring on the south pole

    • MWflip: flipped McEwen & Wiaux scheme, equidistant, first ring on the north pole, last ring half a ring width from the south pole

    • GL: Gauss-Legendre, non-equidistant

    • DH: Driscoll-Healy, equidistant, first ring on north pole, last ring one ring width from the south pole

    • F2: Fejer’s second rule, equidistant, first and last ring one ring width from the poles.

  • phi0 (float) – the azimuth (in radians) of the first pixel in each ring

  • nthreads (int >= 0) – the number of threads to use for the computation. If 0, use as many threads as there are hardware threads available on the system

  • mode (str) –

    the transform mode
    ”STANDARD”: standard transform
    ”GRAD_ONLY”: only valid for spin>0, curl a_lm are not computed
    ”DERIV1”: same as “GRAD_ONLY”, but spin is assumed to be 1 and a_lm are multiplied by sqrt(l*(l+1))

Returns:

the computed spherical harmonic coefficients If the alm parameter was specified, this is identical to alm.

Return type:

numpy.ndarray((nalm, x), dtype=numpy.complex64 or numpy.complex128)

Notes

nmaps = 1 if spin == 0 else 2 nalm = 1 if spin == 0 else (2 if mode == “STANDARD” else 1)

ducc0.sht.adjoint_synthesis_general(*, map: numpy.ndarray, spin: int, lmax: int, loc: numpy.ndarray, epsilon: float = 1e-05, mstart: object = None, lstride: int = 1, mmax: object = None, nthreads: int = 1, alm: object = None, sigma_min: float = 1.1, sigma_max: float = 2.6, mode: str = 'STANDARD', verbose: bool = False) numpy.ndarray

This is the adjoint operation of synthesis_general.

Parameters:
  • map (numpy.ndarray((nmaps, npix), dtype=numpy.float32 or numpy.float64) – The pixel values at the locations specified by loc.

  • spin (int >= 0) – the spin to use for the transform.

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive).

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • mmax (int >= 0 and <= lmax) – the maximum m moment of the transform (inclusive). If not supplied, mmax is derived from mstart; if that is not supplied either, it is assumed to be equal to lmax.

  • loc (numpy.array((npix, 2), dtype=numpy.float64)) – the locations on the sphere at which the alm should be evaluated. loc[:, 0] contains colatitude values (range [0;pi]), loc[:, 1] contains longitude values (range [0;2pi])

  • epsilon (float) – desired accuracy for single precision inputs, this must be >1e-6, for double precision it must be >2e-13

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

  • alm (None or numpy.ndarray((nalm, x), dtype=complex, same accuracy as map)) – the set(s) of spherical harmonic coefficients. The last dimension must be large enough to accommodate all entries, which are stored according to the healpy convention. If None, a new suitable array is allocated.

  • sigma_min (float) – minimum and maximum allowed oversampling factors for the NUFFT component 1.2 <= sigma_min < sigma_max <= 2.5

  • sigma_max (float) – minimum and maximum allowed oversampling factors for the NUFFT component 1.2 <= sigma_min < sigma_max <= 2.5

  • mode (str) –

    the transform mode
    ”STANDARD”: standard transform
    ”GRAD_ONLY”: only valid for spin>0, curl a_lm are not computed
    ”DERIV1”: same as “GRAD_ONLY”, but spin is assumed to be 1 and a_lm are multiplied by sqrt(l*(l+1))

Returns:

the computed spherical harmonic coefficients If the alm parameter was specified, this is identical to alm.

Return type:

numpy.ndarray((nalm, x), dtype=complex, same accuracy as map)

Notes

nmaps = 1 if spin == 0 else 2 nalm = 1 if spin == 0 else (2 if mode == “STANDARD” else 1)

ducc0.sht.alm2leg(*, alm: numpy.ndarray, lmax: int, theta: numpy.ndarray, spin: int = 0, mval: object = None, mstart: object = None, lstride: int = 1, nthreads: int = 1, leg: object = None, mode: str = 'STANDARD', theta_interpol: bool = False) numpy.ndarray

Transforms a set of spherical harmonic coefficients to Legendre coefficients dependent on theta and m.

Parameters:
  • alm (numpy.ndarray((nalm, x), dtype=numpy.complex64 or numpy.complex128)) – the set of spherical harmonic coefficients. The second dimension must be large enough to accommodate all entries, which are stored according to the parameters lmax, ‘mval`, mstart, and lstride.

  • leg (None or numpy.ndarray((nleg, ntheta, nm), same dtype as alm)) – output array containing the Legendre coefficients if None, a new suitable array is allocated

  • spin (int >= 0) – the spin to use for the transform

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive)

  • mval (numpy.ndarray((nm,), dtype = numpy.uint64)) – the m moments for which the transform should be carried out entries must be unique and <= lmax

  • mstart (numpy.ndarray((nm,), dtype = numpy.uint64)) – the (hypothetical) index in the second dimension of alm on which the entry with l=0, m=mval[mi] would be stored, for mi in mval

  • lstride (int) – the index stride in the second dimension of alm between the entries for l and l+1, but the same m.

  • theta (numpy.ndarray((ntheta,), dtype=numpy.float64)) – the colatitudes of the map rings

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

  • mode (str) –

    the transform mode
    ”STANDARD”: standard transform
    ”GRAD_ONLY”: only valid for spin>0, curl a_lm are assumed to be zero
    ”DERIV1”: same as “GRAD_ONLY”, but spin is assumed to be 1 and a_lm are multiplied by sqrt(l*(l+1))

  • theta_interpol (bool) – if the input grid is irregularly spaced in theta, try to accelerate the transform by using an intermediate equidistant theta grid and a 1D NUFFT.

Returns:

the Legendre coefficients. If leg was supplied, this will be the same object.

Return type:

numpy.ndarray((nleg, ntheta, nm), same dtype as alm)

Notes

nleg = 1 if spin == 0 else 2 nalm = 1 if spin == 0 else (2 if mode == “STANDARD” else 1)

ducc0.sht.alm2leg_deriv1(*, alm: numpy.ndarray, lmax: int, theta: numpy.ndarray, mval: object = None, mstart: object = None, lstride: int = 1, nthreads: int = 1, leg: object = None, theta_interpol: bool = False) numpy.ndarray

Transforms a set of spin-0 spherical harmonic coefficients to Legendre coefficients of the first derivatives with respect to colatiude and longitude, dependent on theta and m.

Parameters:
  • alm (numpy.ndarray((1, x), dtype=numpy.complex64 or numpy.complex128)) – the set of spherical harmonic coefficients. The second dimension must be large enough to accommodate all entries, which are stored according to the parameters lmax, ‘mval`, mstart, and lstride.

  • leg (None or numpy.ndarray((2, ntheta, nm), same dtype as alm)) – output array containing the Legendre coefficients if None, a new suitable array is allocated

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive)

  • mval (numpy.ndarray((nm,), dtype = numpy.uint64)) – the m moments for which the transform should be carried out entries must be unique and <= lmax

  • mstart (numpy.ndarray((nm,), dtype = numpy.uint64)) – the (hypothetical) index in the second dimension of alm on which the entry with l=0, m=mval[mi] would be stored, for mi in mval

  • lstride (int) – the index stride in the second dimension of alm between the entries for l and l+1, but the same m.

  • theta (numpy.ndarray((ntheta,), dtype=numpy.float64)) – the colatitudes of the map rings

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

  • theta_interpol (bool) – if the input grid is irregularly spaced in theta, try to accelerate the transform by using an intermediate equidistant theta grid and a 1D NUFFT.

Returns:

the Legendre coefficients. If leg was supplied, this will be the same object. The first component contains coefficients representing a map of df/dtheta; the second component those of 1./sin(theta) df/dphi.

Return type:

numpy.ndarray((2, ntheta, nm), same dtype as alm)

ducc0.sht.analysis_2d(*, map: numpy.ndarray, spin: int, lmax: int, geometry: str, mmax: object = None, nthreads: int = 1, alm: object = None, phi0: float = 0.0, mstart: object = None, lstride: int = 1) numpy.ndarray

Transforms one or two 2D maps to spherical harmonic coefficients. This is the inverse operation of synthesis_2D.

Parameters:
  • alm (numpy.ndarray((ncomp, x), dtype=numpy.complex64 or numpy.complex128)) – storage for the spherical harmonic coefficients. The second dimension must be large enough to accommodate all entries, which are stored according to the healpy convention If not supplied, a new array is allocated.

  • map (numpy.ndarray((ncomp, ntheta, nphi), dtype=numpy.float of same accuracy as alm)) – The input map.

  • spin (int >= 0) – the spin to use for the transform. If spin==0, ncomp must be 1, otherwise 2

  • lmax (int >= 0) – the maximum l (and m) moment of the transform (inclusive)

  • mmax (int >= 0 <= lmax) – the maximum m moment of the transform (inclusive).

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • geometry (one of "CC", "F1", "MW", "MWflip", "GL", "DH", "F2") –

    the distribution of rings over the theta range
    • CC: Clenshaw-Curtis, equidistant, first and last ring on poles

    • F1: Fejer’s first rule, equidistant, first and last ring half a ring width from the poles

    • MW: McEwen & Wiaux scheme, equidistant, first ring half a ring width from the north pole, last ring on the south pole

    • MWflip: flipped McEwen & Wiaux scheme, equidistant, first ring on the north pole, last ring half a ring width from the south pole

    • GL: Gauss-Legendre, non-equidistant

    • DH: Driscoll-Healy, equidistant, first ring on north pole, last ring one ring width from the south pole

    • F2: Fejer’s second rule, equidistant, first and last ring one ring width from the poles.

  • phi0 (float) – the azimuth (in radians) of the first pixel in each ring

  • nthreads (int >= 0) – the number of threads to use for the computation. If 0, use as many threads as there are hardware threads available on the system

Returns:

the computed spherical harmonic coefficients If the alm parameter was specified, this is identical to alm.

Return type:

numpy.ndarray((ncomp, x), dtype=numpy.complex64 or numpy.complex128)

Notes

The maximum m moment to which this function can analyze its input map is min(lmax, (nphi-1)//2).

The maximum l moment to which this function can analyze its input map depends on the geometry, and is

  • ntheta-2 for CC

  • ntheta-1 for F1, MW, MWflip, and GL

  • (ntheta-2)//2 for DH

  • (ntheta-1)//2 for F2

For the CC and F1 geometries this limit is considerably higher than the one obtainable by simply applying quadrature weights. This improvement is achieved by temporary upsampling along meridians to apply the weights at a higher resolution.

ducc0.sht.get_gridweights(geometry: str, ntheta: int) numpy.ndarray

Returns the quadrature weights for a given grid geometry and number of rings.

Parameters:
  • geometry (one of "CC", "F1", "MW", "MWflip", "GL", "DH", "F2") –

    the distribution of rings over the theta range
    • CC: Clenshaw-Curtis, equidistant, first and last ring on poles

    • F1: Fejer’s first rule, equidistant, first and last ring half a ring width from the poles

    • MW: McEwen & Wiaux scheme, equidistant, first ring half a ring width from the north pole, last ring on the south pole

    • MWflip: flipped McEwen & Wiaux scheme, equidistant, first ring on the north pole, last ring half a ring width from the south pole

    • GL: Gauss-Legendre, non-equidistant

    • DH: Driscoll-Healy, equidistant, first ring on north pole, last ring one ring width from the south pole

    • F2: Fejer’s second rule, equidistant, first and last ring one ring width from the poles.

  • ntheta (int > 0) – number of rings in the grid.

Returns:

the quadrature weights for the individual rings. Please note that these weights need to be divided by the number of pixels per ring to obtain the actual quadrature weights for a particular map.

Return type:

numpy.ndarray((ntheta,), dtype=numpy.float64)

ducc0.sht.leg2alm(*, leg: numpy.ndarray, lmax: int, theta: numpy.ndarray, spin: int = 0, mval: object = None, mstart: object = None, lstride: int = 1, nthreads: int = 1, alm: object = None, mode: str = 'STANDARD', theta_interpol: bool = False) numpy.ndarray

Transforms a set of Legendre coefficients to spherical harmonic coefficients

Parameters:
  • leg (numpy.ndarray((nleg, ntheta, nm), dtype=numpy.complex64 or numpy.complex128)) –

  • alm (None or numpy.ndarray((nalm, x), same dtype as leg)) – the set of spherical harmonic coefficients. The second dimension must be large enough to accommodate all entries, which are stored according to the parameters lmax, ‘mval`, mstart, and lstride. if None, a new suitable array is allocated

  • spin (int >= 0) – the spin to use for the transform

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive)

  • mval (numpy.ndarray((nm,), dtype = numpy.uint64)) – the m moments for which the transform should be carried out entries must be unique and <= lmax

  • mstart (numpy.ndarray((nm,), dtype = numpy.uint64)) – the (hypothetical) index in the second dimension of alm on which the entry with l=0, m=mval[mi] would be stored, for mi in mval

  • lstride (int) – the index stride in the second dimension of alm between the entries for l and l+1, but the same m.

  • theta (numpy.ndarray((ntheta,), dtype=numpy.float64)) – the colatitudes of the map rings

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

  • mode (str) –

    the transform mode
    ”STANDARD”: standard transform
    ”GRAD_ONLY”: only valid for spin>0, curl a_lm are not computed
    ”DERIV1”: same as “GRAD_ONLY”, but spin is assumed to be 1 and a_lm are multiplied by sqrt(l*(l+1))

  • theta_interpol (bool) – if the input grid is irregularly spaced in theta, try to accelerate the transform by using an intermediate equidistant theta grid and a 1D NUFFT.

Returns:

the Legendre coefficients. if alm was supplied, this will be the same object If newly allocated, the smallest possible second dimensions will be chosen.

Return type:

numpy.ndarray((nalm, *), same dtype as leg)

Notes

nleg = 1 if spin == 0 else 2 nalm = 1 if spin == 0 else (2 if mode == “STANDARD” else 1)

ducc0.sht.leg2map(*, leg: numpy.ndarray, nphi: numpy.ndarray, phi0: numpy.ndarray, ringstart: numpy.ndarray, pixstride: int = 1, nthreads: int = 1, map: object = None) numpy.ndarray

Transforms one or more sets of Legendre coefficients to maps.

Parameters:
  • leg (numpy.ndarray((ncomp, ntheta, mmax+1), numppy.complex64 or numpy.complex128)) – input array containing the Legendre coefficients. The entries in leg[:,:,m] correspond to quantum number m, i.e. the m values must be stored in ascending order, and complete.

  • map (None or numpy.ndarray((ncomp, x), dtype=numpy.float of same accuracy as leg) – the map pixel data. The second dimension must be large enough to accommodate all pixels, which are stored according to the parameters nphi, ‘ringstart`, and pixstride. if None, a new suitable array is allocated

  • nphi (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – number of pixels in every ring

  • phi0 (numpy.ndarray((ntheta,), dtype=numpy.float64)) – azimuth (in radians) of the first pixel in every ring

  • ringstart (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – the index in the second dimension of map at which the first pixel of every ring is stored

  • pixstride (int) – the index stride in the second dimension of map between two subsequent pixels in a ring

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

Returns:

the map pixel data. If map was supplied, this will be the same object If newly allocated, the smallest possible second dimensions will be chosen.

Return type:

numpy.ndarray((ncomp, x), dtype=numpy.float of same accuracy as leg)

Notes

In contrast to leg2alm and alm2leg the m values are assumed to form a range from 0 to mmax, inclusively.

ducc0.sht.map2leg(*, map: numpy.ndarray, nphi: numpy.ndarray, phi0: numpy.ndarray, ringstart: numpy.ndarray, mmax: int, pixstride: int = 1, nthreads: int = 1, leg: object = None) numpy.ndarray

Transforms a map or several maps to Legendre coefficients dependent on theta and m.

Parameters:
  • map (numpy.ndarray((ncomp, x), dtype=numpy.float32 or numpy.float64)) – the map pixel data. The second dimension must be large enough to accommodate all pixels, which are stored according to the parameters nphi, ‘ringstart`, and pixstride.

  • leg (None or numpy.ndarray((ncomp, ntheta, mmax+1), dtype=numpy.complex of same accuracy as map)) – output array containing the Legendre coefficients if None, a new suitable array is allocated The entries in leg[:,:,m] correspond to quantum number m, i.e. the m values will be stored in ascending order.

  • nphi (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – number of pixels in every ring

  • phi0 (numpy.ndarray((ntheta,), dtype=numpy.float64)) – azimuth (in radians) of the first pixel in every ring

  • ringstart (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – the index in the second dimension of map at which the first pixel of every ring is stored

  • pixstride (int) – the index stride in the second dimension of map between two subsequent pixels in a ring

  • mmax (int) – the maximum m moment to compute in this transform. If leg is provided, mmax must be equal to leg.shape[2]=1.

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

Returns:

the Legendre coefficients if leg was supplied, this will be the same object

Return type:

numpy.ndarray((ncomp, ntheta, nm), dtype=numpy.complex of same accuracy as map)

Notes

In contrast to leg2alm and alm2leg the m values are assumed to form a range from 0 to mmax, inclusively.

ducc0.sht.maxium_safe_l(geometry: str, ntheta: int) int

Returns the maximum l moment that can be safely stored (i.e. is guaranteed to be recoverable using analysis_2d) in a map with the specified geometry and number of rings.

Parameters:
  • geometry (one of "CC", "F1", "MW", "MWflip", "GL", "DH", "F2") –

    the distribution of rings over the theta range
    • CC: Clenshaw-Curtis, equidistant, first and last ring on poles

    • F1: Fejer’s first rule, equidistant, first and last ring half a ring width from the poles

    • MW: McEwen & Wiaux scheme, equidistant, first ring half a ring width from the north pole, last ring on the south pole

    • MWflip: flipped McEwen & Wiaux scheme, equidistant, first ring on the north pole, last ring half a ring width from the south pole

    • GL: Gauss-Legendre, non-equidistant

    • DH: Driscoll-Healy, equidistant, first ring on north pole, last ring one ring width from the south pole

    • F2: Fejer’s second rule, equidistant, first and last ring one ring width from the poles.

  • ntheta (int > 0 (> 1 for geometries CC, DH, and F1)) – number of rings in the grid.

Returns:

lmax – The maximum l moment that can be safely stored on the specified grid.

Return type:

int>=0

ducc0.sht.pseudo_analysis(*, map: numpy.ndarray, theta: numpy.ndarray, lmax: int, mstart: object = None, nphi: numpy.ndarray, phi0: numpy.ndarray, ringstart: numpy.ndarray, spin: int, lstride: int = 1, pixstride: int = 1, nthreads: int = 1, alm: object = None, maxiter: int, epsilon: float, mmax: object = None, theta_interpol: bool = False) object

Tries to extract spherical harmonic coefficients from (sets of) one or two maps by using the iterative LSMR algorithm.

Parameters:
  • alm (None or numpy.ndarray(([ntrans,] ncomp, x), dtype=numpy.complex of same precision as map)) – the set of spherical harmonic coefficients. The last dimension must be large enough to accommodate all entries, which are stored according to the parameters lmax, ‘mmax`, mstart, and lstride. if None, a new suitable array is allocated

  • map (numpy.ndarray(([ntrans,] ncomp, x), dtype=numpy.float32 or numpy.float64) – The last dimension must be large enough to accommodate all pixels, which are stored according to the parameters nphi, ‘ringstart`, and pixstride.

  • theta (numpy.ndarray((ntheta,), dtype=numpy.float64)) – the colatitudes of the map rings

  • nphi (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – number of pixels in every ring

  • phi0 (numpy.ndarray((ntheta,), dtype=numpy.float64)) – azimuth (in radians) of the first pixel in every ring

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • ringstart (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – the index in the last dimension of map at which the first pixel of every ring is stored

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • pixstride (int) – the index stride in the last dimension of map between two subsequent pixels in a ring

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

  • spin (int >= 0) – the spin to use for the transform. If spin==0, ncomp must be 1, otherwise 2

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive).

  • maxiter (int >= 0) – the maximum number of iterations before stopping the algorithm

  • epsilon (float >= 0) – the relative tolerance used as a stopping criterion

  • mmax (int >= 0 <= lmax) – the maximum m moment of the transform (inclusive).

  • theta_interpol (bool) – if the input grid is irregularly spaced in theta, try to accelerate the transform by using an intermediate equidistant theta grid and a 1D NUFFT.

Returns:

  • numpy.ndarray(([ntrans,] ncomp, x), dtype=numpy.complex of same accuracy as map) – the set(s) of spherical harmonic coefficients. If alm was supplied, this will be the same object If newly allocated, the smallest possible last dimension will be chosen.

  • int or list(int) – the reason for stopping the iteration 1: approximate solution to the equation system found 2: approximate least-squares solution found 3: condition number of the equation system too large 7: maximum number of iterations reached

  • int or list(int) – the iteration count(s)

  • float or list(float) – the residual norm, divided by the norm of map

  • float or list(float) – the quality of the least-squares solution

ducc0.sht.pseudo_analysis_general(*, lmax: int, map: numpy.ndarray, loc: numpy.ndarray, spin: int, nthreads: int, maxiter: int, epsilon: float = 1e-05, sigma_min: float = 1.1, sigma_max: float = 2.6, mstart: object = None, lstride: int = 1, mmax: object = None, alm: object = None) tuple

Tries to extract spherical harmonic coefficients from one or two maps by using the iterative LSMR algorithm.

Parameters:
  • map (numpy.ndarray((ncomp, npix), dtype=numpy.float32 or numpy.float64) – The pixel values at the locations specified by loc.

  • spin (int >= 0) – the spin to use for the transform. If spin==0, ncomp must be 1, otherwise 2

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive).

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • mmax (int >= 0 and <= lmax) – the maximum m moment of the transform (inclusive). If not supplied, mmax is derived from mstart; if that is not supplied either, it is assumed to be equal to lmax.

  • loc (numpy.array((npix, 2), dtype=numpy.float64)) – the locations on the sphere at which the alm should be evaluated. loc[:, 0] contains colatitude values (range [0;pi]), loc[:, 1] contains longitude values (range [0;2pi])

  • epsilon (float >= 0) – the relative tolerance used as a stopping criterion NOTE: for the “epsilon” paraeter of the underlyig NUFFT calls, 0.1*epsilon will be used.

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

  • alm (None or numpy.ndarray((ncomp, x), dtype=complex, same accuracy as map)) – the set(s) of spherical harmonic coefficients. The last dimension must be large enough to accommodate all entries, which are stored according to the healpy convention. If None, a new suitable array is allocated.

  • sigma_min (float) – minimum and maximum allowed oversampling factors for the NUFFT component 1.2 <= sigma_min < sigma_max <= 2.5

  • sigma_max (float) – minimum and maximum allowed oversampling factors for the NUFFT component 1.2 <= sigma_min < sigma_max <= 2.5

  • maxiter (int >= 0) – the maximum number of iterations before stopping the algorithm

Returns:

  • numpy.ndarray((ncomp, x), dtype=complex, same accuracy as map) – the computed spherical harmonic coefficients If the alm parameter was specified, this is identical to alm. If newly allocated, the smallest possible last dimension will be chosen.

  • int – the reason for stopping the iteration 1: approximate solution to the equation system found 2: approximate least-squares solution found 3: condition number of the equation system too large 7: maximum number of iterations reached

  • int – the iteration count(s)

  • float – the residual norm, divided by the norm of map

  • float – the quality of the least-squares solution

ducc0.sht.rotate_alm(alm: numpy.ndarray, lmax: int, psi: float, theta: float, phi: float, nthreads: int = 1) numpy.ndarray

Rotates a set of spherical harmonic coefficients according to the given Euler angles.

Parameters:
  • alm (numpy.ndarray(((lmax+1)*(lmax=2)/2,), dtype=numpy complex64 or numpy.complex128)) – the spherical harmonic coefficients, in the order (0,0), (1,0), (2,0), … (lmax,0), (1,1), (2,1), …, (lmax, lmax)

  • lmax (int >= 0) – Maximum multipole order l of the data set.

  • psi (float) – First rotation angle about the z-axis. All angles are in radians, the rotations are active and the referential system is assumed to be right handed.

  • theta (float) – Second rotation angl about the original (unrotated) y-axis

  • phi (float) – Third rotation angle about the original (unrotated) z-axis.

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

Return type:

numpy.ndarray(same shape and dtype as alm)

class ducc0.sht.sharpjob_d

Interface class to some of libsharp2’s functionality.

Notes

This class is considered obsolescent and will be removed in the future.

alm2map(self: ducc0.sht.sharpjob_d, alm: numpy.ndarray[numpy.complex128]) numpy.ndarray
alm2map_adjoint(self: ducc0.sht.sharpjob_d, map: numpy.ndarray[numpy.float64]) numpy.ndarray
alm2map_spin(self: ducc0.sht.sharpjob_d, alm: numpy.ndarray[numpy.complex128], spin: int) numpy.ndarray
map2alm(self: ducc0.sht.sharpjob_d, map: numpy.ndarray[numpy.float64]) numpy.ndarray
map2alm_spin(self: ducc0.sht.sharpjob_d, map: numpy.ndarray[numpy.float64], spin: int) numpy.ndarray
n_alm(self: ducc0.sht.sharpjob_d) int
set_cc_geometry(self: ducc0.sht.sharpjob_d, ntheta: int, nphi: int) None
set_dh_geometry(self: ducc0.sht.sharpjob_d, ntheta: int, nphi: int) None
set_fejer1_geometry(self: ducc0.sht.sharpjob_d, ntheta: int, nphi: int) None
set_fejer2_geometry(self: ducc0.sht.sharpjob_d, ntheta: int, nphi: int) None
set_gauss_geometry(self: ducc0.sht.sharpjob_d, ntheta: int, nphi: int) None
set_healpix_geometry(self: ducc0.sht.sharpjob_d, nside: int) None
set_mw_geometry(self: ducc0.sht.sharpjob_d, ntheta: int, nphi: int) None
set_nthreads(self: ducc0.sht.sharpjob_d, nthreads: int) None
set_triangular_alm_info(self: ducc0.sht.sharpjob_d, lmax: int, mmax: int) None
ducc0.sht.synthesis(*, alm: numpy.ndarray, theta: numpy.ndarray, lmax: int, mstart: object = None, nphi: numpy.ndarray, phi0: numpy.ndarray, ringstart: numpy.ndarray, spin: int, lstride: int = 1, pixstride: int = 1, nthreads: int = 1, map: object = None, mmax: object = None, mode: str = 'STANDARD', theta_interpol: bool = False) numpy.ndarray

Transforms (sets of) one or two sets of spherical harmonic coefficients to maps on the sphere.

Parameters:
  • alm (numpy.ndarray(([ntrans,] nalm, x), dtype=numpy.complex64 or numpy.complex128)) – the set(s) of spherical harmonic coefficients. The last dimension must be large enough to accommodate all entries, which are stored according to the parameters lmax, ‘mmax`, mstart, and lstride.

  • map (None or numpy.ndarray(([ntrans,] nmaps, x), dtype=numpy.float of same accuracy as alm) – the map pixel data. The last dimension must be large enough to accommodate all pixels, which are stored according to the parameters nphi, ‘ringstart`, and pixstride. if None, a new suitable array is allocated

  • theta (numpy.ndarray((ntheta,), dtype=numpy.float64)) – the colatitudes of the map rings

  • nphi (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – number of pixels in every ring

  • phi0 (numpy.ndarray((ntheta,), dtype=numpy.float64)) – azimuth (in radians) of the first pixel in every ring

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • ringstart (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – the index in the last dimension of map at which the first pixel of every ring is stored

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • pixstride (int) – the index stride in the last dimension of map between two subsequent pixels in a ring

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

  • spin (int >= 0) – the spin to use for the transform.

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive).

  • mmax (int >= 0 <= lmax) – the maximum m moment of the transform (inclusive).

  • mode (str) –

    the transform mode
    ”STANDARD”: standard transform
    ”GRAD_ONLY”: only valid for spin>0, curl a_lm are assumed to be zero
    ”DERIV1”: same as “GRAD_ONLY”, but spin is assumed to be 1 and a_lm are multiplied by sqrt(l*(l+1))

  • theta_interpol (bool) – if the input grid is irregularly spaced in theta, try to accelerate the transform by using an intermediate equidistant theta grid and a 1D NUFFT.

Returns:

the map pixel data. If map was supplied, this will be the same object If newly allocated, the smallest possible last dimension will be chosen.

Return type:

numpy.ndarray(([ntrans,] nmaps, x), dtype=numpy.float of same accuracy as alm)

Notes

nmaps = 1 if spin == 0 else 2 nalm = 1 if spin == 0 else (2 if mode == “STANDARD” else 1)

ducc0.sht.synthesis_2d(*, alm: numpy.ndarray, spin: int, lmax: int, geometry: str, ntheta: object = None, nphi: object = None, mmax: object = None, nthreads: int = 1, map: object = None, mode: str = 'STANDARD', phi0: float = 0.0, mstart: object = None, lstride: int = 1) numpy.ndarray

Transforms one or two sets of spherical harmonic coefficients to 2D maps.

Parameters:
  • alm (numpy.ndarray((nalm, x), dtype=numpy.complex64 or numpy.complex128)) – the set of spherical harmonic coefficients. The second dimension must be large enough to accommodate all entries, which are stored according to the healpy convention.

  • map (numpy.ndarray((nmaps, ntheta, nphi), dtype=numpy.float of same accuracy as alm)) – storage for the output map. If not supplied, a new array is allocated.

  • ntheta (int > 0) – dimensions of the output map. If not supplied, map must be supplied. If supplied, and map is also supplied, must match with the array dimensions

  • nphi (int > 0) – dimensions of the output map. If not supplied, map must be supplied. If supplied, and map is also supplied, must match with the array dimensions

  • spin (int >= 0) – the spin to use for the transform.

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive).

  • mmax (int >= 0 <= lmax) – the maximum m moment of the transform (inclusive).

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • geometry (one of "CC", "F1", "MW", "MWflip", "GL", "DH", "F2") –

    the distribution of rings over the theta range
    • CC: Clenshaw-Curtis, equidistant, first and last ring on poles

    • F1: Fejer’s first rule, equidistant, first and last ring half a ring width from the poles

    • MW: McEwen & Wiaux scheme, equidistant, first ring half a ring width from the north pole, last ring on the south pole

    • MWflip: flipped McEwen & Wiaux scheme, equidistant, first ring on the north pole, last ring half a ring width from the south pole

    • GL: Gauss-Legendre, non-equidistant

    • DH: Driscoll-Healy, equidistant, first ring on north pole, last ring one ring width from the south pole

    • F2: Fejer’s second rule, equidistant, first and last ring one ring width from the poles.

  • phi0 (float) – the azimuth (in radians) of the first pixel in each ring

  • nthreads (int >= 0) – the number of threads to use for the computation. If 0, use as many threads as there are hardware threads available on the system

  • mode (str) –

    the transform mode
    ”STANDARD”: standard transform
    ”GRAD_ONLY”: only valid for spin>0, curl a_lm are assumed to be zero
    ”DERIV1”: same as “GRAD_ONLY”, but spin is assumed to be 1 and a_lm are multiplied by sqrt(l*(l+1))

Returns:

the computed map. If the map parameter was specified, this is identical with map.

Return type:

numpy.ndarray((nmaps, ntheta, nphi), dtype=numpy.float of same accuracy as alm)

Notes

nmaps = 1 if spin == 0 else 2 nalm = 1 if spin == 0 else (2 if mode == “STANDARD” else 1)

ducc0.sht.synthesis_2d_deriv1(*, alm: numpy.ndarray, lmax: int, geometry: str, ntheta: object = None, nphi: object = None, mmax: object = None, nthreads: int = 1, map: object = None, phi0: float = 0.0, mstart: object = None, lstride: int = 1) numpy.ndarray

Transforms a set of spherical harmonic coefficients to two 2D maps containing the derivatives with respect to theta and phi.

Parameters:
  • alm (numpy.ndarray((1, x), dtype=numpy.complex64 or numpy.complex128)) – the set of spherical harmonic coefficients. The second dimension must be large enough to accommodate all entries, which are stored according to the healpy convention.

  • map (numpy.ndarray((2, ntheta, nphi), dtype=numpy.float of same accuracy as alm)) – storage for the output map. If not supplied, a new array is allocated.

  • ntheta (int > 0) – dimensions of the output map. If not supplied, map must be supplied. If supplied, and map is also supplied, must match with the array dimensions

  • nphi (int > 0) – dimensions of the output map. If not supplied, map must be supplied. If supplied, and map is also supplied, must match with the array dimensions

  • lmax (int >= 0) – the maximum l (and m) moment of the transform (inclusive)

  • mmax (int >= 0 <= lmax) – the maximum m moment of the transform (inclusive).

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • geometry (one of "CC", "F1", "MW", "MWflip", "GL", "DH", "F2") –

    the distribution of rings over the theta range
    • CC: Clenshaw-Curtis, equidistant, first and last ring on poles

    • F1: Fejer’s first rule, equidistant, first and last ring half a ring width from the poles

    • MW: McEwen & Wiaux scheme, equidistant, first ring half a ring width from the north pole, last ring on the south pole

    • MWflip: flipped McEwen & Wiaux scheme, equidistant, first ring on the north pole, last ring half a ring width from the south pole

    • GL: Gauss-Legendre, non-equidistant

    • DH: Driscoll-Healy, equidistant, first ring on north pole, last ring one ring width from the south pole

    • F2: Fejer’s second rule, equidistant, first and last ring one ring width from the poles.

  • phi0 (float) – the azimuth (in radians) of the first pixel in each ring

  • nthreads (int >= 0) – the number of threads to use for the computation. If 0, use as many threads as there are hardware threads available on the system

Returns:

the maps containing the derivatives with respect to theta and phi. If the map parameter was specified, this is identical with map.

Return type:

numpy.ndarray((2, ntheta, nphi), dtype=numpy.float of same accuracy as alm)

ducc0.sht.synthesis_deriv1(*, alm: numpy.ndarray, theta: numpy.ndarray, lmax: int, mstart: object = None, nphi: numpy.ndarray, phi0: numpy.ndarray, ringstart: numpy.ndarray, lstride: int = 1, pixstride: int = 1, nthreads: int = 1, map: object = None, mmax: object = None, theta_interpol: bool = False) numpy.ndarray

Transforms a set (or sets) of spherical harmonic coefficients to two maps containing the derivatives with respect to theta and phi.

Parameters:
  • alm (numpy.ndarray(([ntrans,] 1, x), dtype=numpy.complex64 or numpy.complex128)) – the set(s) of spherical harmonic coefficients. The last dimension must be large enough to accommodate all entries, which are stored according to the parameters lmax, ‘mmax`, mstart, and lstride.

  • map (None or numpy.ndarray(([ntrans,] 2, x), dtype=numpy.float of same accuracy as alm) – the map pixel data. The last dimension must be large enough to accommodate all pixels, which are stored according to the parameters nphi, ‘ringstart`, and pixstride. if None, a new suitable array is allocated

  • theta (numpy.ndarray((ntheta,), dtype=numpy.float64)) – the colatitudes of the map rings

  • nphi (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – number of pixels in every ring

  • phi0 (numpy.ndarray((ntheta,), dtype=numpy.float64)) – azimuth (in radians) of the first pixel in every ring

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • ringstart (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – the index in the last dimension of map at which the first pixel of every ring is stored

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • pixstride (int) – the index stride in the last dimension of map between two subsequent pixels in a ring

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive).

  • mmax (int >= 0 <= lmax) – the maximum m moment of the transform (inclusive).

  • theta_interpol (bool) – if the input grid is irregularly spaced in theta, try to accelerate the transform by using an intermediate equidistant theta grid and a 1D NUFFT.

Returns:

the map pixel data. If map was supplied, this will be the same object If newly allocated, the smallest possible last dimension will be chosen.

Return type:

numpy.ndarray((2, x), dtype=numpy.float of same accuracy as alm)

ducc0.sht.synthesis_general(*, alm: numpy.ndarray, spin: int, lmax: int, loc: numpy.ndarray, epsilon: float = 1e-05, mstart: object = None, lstride: int = 1, mmax: object = None, nthreads: int = 1, map: object = None, sigma_min: float = 1.1, sigma_max: float = 2.6, mode: str = 'STANDARD', verbose: bool = False) numpy.ndarray

Evaluate a_lm at arbitrary positions on the sphere

Parameters:
  • alm (numpy.ndarray((nalm, x), dtype=numpy.complex64 or numpy.complex128)) – the set(s) of spherical harmonic coefficients. The last dimension must be large enough to accommodate all entries, which are stored according to the healpy convention.

  • spin (int >= 0) – the spin to use for the transform.

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive).

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • mmax (int >= 0 and <= lmax) – the maximum m moment of the transform (inclusive). If not supplied, mmax is derived from mstart; if that is not supplied either, it is assumed to be equal to lmax.

  • loc (numpy.array((npix, 2), dtype=numpy.float64)) – the locations on the sphere at which the alm should be evaluated. loc[:, 0] contains colatitude values (range [0;pi]), loc[:, 1] contains longitude values (range [0;2pi])

  • epsilon (float) – desired accuracy for single precision inputs, this must be >1e-6, for double precision it must be >2e-13

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

  • map (None or numpy.ndarray((nmaps, npix), dtype=numpy.float of same accuracy as alm) – the map pixel data. If None, a new suitable array is allocated.

  • sigma_min (float) – minimum and maximum allowed oversampling factors for the NUFFT component 1.2 <= sigma_min < sigma_max <= 2.5

  • sigma_max (float) – minimum and maximum allowed oversampling factors for the NUFFT component 1.2 <= sigma_min < sigma_max <= 2.5

  • mode (str) –

    the transform mode
    ”STANDARD”: standard transform
    ”GRAD_ONLY”: only valid for spin>0, curl a_lm are assumed to be zero
    ”DERIV1”: same as “GRAD_ONLY”, but spin is assumed to be 1 and a_lm are multiplied by sqrt(l*(l+1))

Returns:

the pixel values at the locations specified by loc. If the map parameter was specified, this is identical with map.

Return type:

numpy.ndarray((nmaps, npix), dtype=numpy.float of same accuracy as alm

Notes

nmaps = 1 if spin == 0 else 2 nalm = 1 if spin == 0 else (2 if mode == “STANDARD” else 1)

ducc0.sht.experimental

Experimental interface to the SHT functionality.

Notes

The functionality in this module is not considered to have a stable interface and also may be moved to other modules in the future. If you use it, be prepared to adjust your code at some point ion the future!

ducc0.sht.experimental.adjoint_analysis_2d(*, alm: numpy.ndarray, spin: int, lmax: int, geometry: str, ntheta: object = None, nphi: object = None, mmax: object = None, nthreads: int = 1, map: object = None, phi0: float = 0.0, mstart: object = None, lstride: int = 1) numpy.ndarray

Transforms one or two sets of spherical harmonic coefficients to 2D maps. This is the adjoint operation of analysis_2D.

Parameters:
  • alm (numpy.ndarray((ncomp, x), dtype=numpy.complex64 or numpy.complex128)) – the set of spherical harmonic coefficients. The second dimension must be large enough to accommodate all entries, which are stored according to the healpy convention.

  • map (numpy.ndarray((ncomp, ntheta, nphi), dtype=numpy.float of same accuracy as alm)) – storage for the output map. If not supplied, a new array is allocated.

  • ntheta (int > 0) – dimensions of the output map. If not supplied, map must be supplied. If supplied, and map is also supplied, must match with the array dimensions

  • nphi (int > 0) – dimensions of the output map. If not supplied, map must be supplied. If supplied, and map is also supplied, must match with the array dimensions

  • spin (int >= 0) – the spin to use for the transform. If spin==0, ncomp must be 1, otherwise 2

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive).

  • mmax (int >= 0 <= lmax) – the maximum m moment of the transform (inclusive).

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • geometry (one of "CC", "F1", "MW", "MWflip", "GL", "DH", "F2") –

    the distribution of rings over the theta range
    • CC: Clenshaw-Curtis, equidistant, first and last ring on poles

    • F1: Fejer’s first rule, equidistant, first and last ring half a ring width from the poles

    • MW: McEwen & Wiaux scheme, equidistant, first ring half a ring width from the north pole, last ring on the south pole

    • MWflip: flipped McEwen & Wiaux scheme, equidistant, first ring on the north pole, last ring half a ring width from the south pole

    • GL: Gauss-Legendre, non-equidistant

    • DH: Driscoll-Healy, equidistant, first ring on north pole, last ring one ring width from the south pole

    • F2: Fejer’s second rule, equidistant, first and last ring one ring width from the poles.

  • phi0 (float) – the azimuth (in radians) of the first pixel in each ring

  • nthreads (int >= 0) – the number of threads to use for the computation. If 0, use as many threads as there are hardware threads available on the system

Returns:

the computed map. If the map parameter was specified, this is identical with map.

Return type:

numpy.ndarray((ncomp, ntheta, nphi), dtype=numpy.float of same accuracy as alm)

Notes

For limits on lmax and mmax see the documentation of analysis_2d.

ducc0.sht.experimental.adjoint_synthesis(*, map: numpy.ndarray, theta: numpy.ndarray, lmax: int, mstart: object = None, nphi: numpy.ndarray, phi0: numpy.ndarray, ringstart: numpy.ndarray, spin: int, lstride: int = 1, pixstride: int = 1, nthreads: int = 1, alm: object = None, mmax: object = None, mode: str = 'STANDARD', theta_interpol: bool = False) numpy.ndarray

Transforms (sets of) one or two maps to spherical harmonic coefficients. This is the adjoint operation of synthesis.

Parameters:
  • alm (None or numpy.ndarray(([ntrans,] nalm, x), dtype=numpy.complex of same precision as map)) – the set of spherical harmonic coefficients. The last dimension must be large enough to accommodate all entries, which are stored according to the parameters lmax, ‘mmax`, mstart, and lstride. if None, a new suitable array is allocated

  • map (numpy.ndarray(([ntrans,] nmaps, x), dtype=numpy.float32 or numpy.float64) – The last dimension must be large enough to accommodate all pixels, which are stored according to the parameters nphi, ringstart, and pixstride.

  • theta (numpy.ndarray((ntheta,), dtype=numpy.float64)) – the colatitudes of the map rings

  • nphi (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – number of pixels in every ring

  • phi0 (numpy.ndarray((ntheta,), dtype=numpy.float64)) – azimuth (in radians) of the first pixel in every ring

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • ringstart (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – the index in the last dimension of map at which the first pixel of every ring is stored

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • pixstride (int) – the index stride in the last dimension of map between two subsequent pixels in a ring

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

  • spin (int >= 0) – the spin to use for the transform.

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive).

  • mmax (int >= 0 <= lmax) – the maximum m moment of the transform (inclusive).

  • mode (str) –

    the transform mode
    ”STANDARD”: standard transform
    ”GRAD_ONLY”: only valid for spin>0, curl a_lm are not computed
    ”DERIV1”: same as “GRAD_ONLY”, but spin is assumed to be 1 and a_lm are multiplied by sqrt(l*(l+1))

  • theta_interpol (bool) – if the input grid is irregularly spaced in theta, try to accelerate the transform by using an intermediate equidistant theta grid and a 1D NUFFT.

Returns:

the set(s) of spherical harmonic coefficients. If alm was supplied, this will be the same object If newly allocated, the smallest possible last dimension will be chosen.

Return type:

numpy.ndarray(([ntrans,] nalm, x), dtype=numpy.complex of same accuracy as map)

Notes

nmaps = 1 if spin == 0 else 2 nalm = 1 if spin == 0 else (2 if mode == “STANDARD” else 1)

ducc0.sht.experimental.adjoint_synthesis_2d(*, map: numpy.ndarray, spin: int, lmax: int, geometry: str, mmax: object = None, nthreads: int = 1, alm: object = None, mode: str = 'STANDARD', phi0: float = 0.0, mstart: object = None, lstride: int = 1) numpy.ndarray

Transforms one or two 2D maps to spherical harmonic coefficients. This is the adjoint operation of synthesis_2D.

Parameters:
  • alm (numpy.ndarray((nalm, x), dtype=numpy.complex64 or numpy.complex128)) – storage for the spherical harmonic coefficients. The second dimension must be large enough to accommodate all entries, which are stored according to the healpy convention. If not supplied, a new array is allocated.

  • map (numpy.ndarray((nmaps, ntheta, nphi), dtype=numpy.float of same accuracy as alm)) – The input map.

  • spin (int >= 0) – the spin to use for the transform.

  • lmax (int >= 0) – the maximum l (and m) moment of the transform (inclusive)

  • mmax (int >= 0 <= lmax) – the maximum m moment of the transform (inclusive).

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • geometry (one of "CC", "F1", "MW", "MWflip", "GL", "DH", "F2") –

    the distribution of rings over the theta range
    • CC: Clenshaw-Curtis, equidistant, first and last ring on poles

    • F1: Fejer’s first rule, equidistant, first and last ring half a ring width from the poles

    • MW: McEwen & Wiaux scheme, equidistant, first ring half a ring width from the north pole, last ring on the south pole

    • MWflip: flipped McEwen & Wiaux scheme, equidistant, first ring on the north pole, last ring half a ring width from the south pole

    • GL: Gauss-Legendre, non-equidistant

    • DH: Driscoll-Healy, equidistant, first ring on north pole, last ring one ring width from the south pole

    • F2: Fejer’s second rule, equidistant, first and last ring one ring width from the poles.

  • phi0 (float) – the azimuth (in radians) of the first pixel in each ring

  • nthreads (int >= 0) – the number of threads to use for the computation. If 0, use as many threads as there are hardware threads available on the system

  • mode (str) –

    the transform mode
    ”STANDARD”: standard transform
    ”GRAD_ONLY”: only valid for spin>0, curl a_lm are not computed
    ”DERIV1”: same as “GRAD_ONLY”, but spin is assumed to be 1 and a_lm are multiplied by sqrt(l*(l+1))

Returns:

the computed spherical harmonic coefficients If the alm parameter was specified, this is identical to alm.

Return type:

numpy.ndarray((nalm, x), dtype=numpy.complex64 or numpy.complex128)

Notes

nmaps = 1 if spin == 0 else 2 nalm = 1 if spin == 0 else (2 if mode == “STANDARD” else 1)

ducc0.sht.experimental.adjoint_synthesis_general(*, map: numpy.ndarray, spin: int, lmax: int, loc: numpy.ndarray, epsilon: float = 1e-05, mstart: object = None, lstride: int = 1, mmax: object = None, nthreads: int = 1, alm: object = None, sigma_min: float = 1.1, sigma_max: float = 2.6, mode: str = 'STANDARD', verbose: bool = False) numpy.ndarray

This is the adjoint operation of synthesis_general.

Parameters:
  • map (numpy.ndarray((nmaps, npix), dtype=numpy.float32 or numpy.float64) – The pixel values at the locations specified by loc.

  • spin (int >= 0) – the spin to use for the transform.

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive).

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • mmax (int >= 0 and <= lmax) – the maximum m moment of the transform (inclusive). If not supplied, mmax is derived from mstart; if that is not supplied either, it is assumed to be equal to lmax.

  • loc (numpy.array((npix, 2), dtype=numpy.float64)) – the locations on the sphere at which the alm should be evaluated. loc[:, 0] contains colatitude values (range [0;pi]), loc[:, 1] contains longitude values (range [0;2pi])

  • epsilon (float) – desired accuracy for single precision inputs, this must be >1e-6, for double precision it must be >2e-13

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

  • alm (None or numpy.ndarray((nalm, x), dtype=complex, same accuracy as map)) – the set(s) of spherical harmonic coefficients. The last dimension must be large enough to accommodate all entries, which are stored according to the healpy convention. If None, a new suitable array is allocated.

  • sigma_min (float) – minimum and maximum allowed oversampling factors for the NUFFT component 1.2 <= sigma_min < sigma_max <= 2.5

  • sigma_max (float) – minimum and maximum allowed oversampling factors for the NUFFT component 1.2 <= sigma_min < sigma_max <= 2.5

  • mode (str) –

    the transform mode
    ”STANDARD”: standard transform
    ”GRAD_ONLY”: only valid for spin>0, curl a_lm are not computed
    ”DERIV1”: same as “GRAD_ONLY”, but spin is assumed to be 1 and a_lm are multiplied by sqrt(l*(l+1))

Returns:

the computed spherical harmonic coefficients If the alm parameter was specified, this is identical to alm.

Return type:

numpy.ndarray((nalm, x), dtype=complex, same accuracy as map)

Notes

nmaps = 1 if spin == 0 else 2 nalm = 1 if spin == 0 else (2 if mode == “STANDARD” else 1)

ducc0.sht.experimental.alm2flm(alm: numpy.ndarray, spin: int, flm: object = None) numpy.ndarray
ducc0.sht.experimental.alm2leg(*, alm: numpy.ndarray, lmax: int, theta: numpy.ndarray, spin: int = 0, mval: object = None, mstart: object = None, lstride: int = 1, nthreads: int = 1, leg: object = None, mode: str = 'STANDARD', theta_interpol: bool = False) numpy.ndarray

Transforms a set of spherical harmonic coefficients to Legendre coefficients dependent on theta and m.

Parameters:
  • alm (numpy.ndarray((nalm, x), dtype=numpy.complex64 or numpy.complex128)) – the set of spherical harmonic coefficients. The second dimension must be large enough to accommodate all entries, which are stored according to the parameters lmax, ‘mval`, mstart, and lstride.

  • leg (None or numpy.ndarray((nleg, ntheta, nm), same dtype as alm)) – output array containing the Legendre coefficients if None, a new suitable array is allocated

  • spin (int >= 0) – the spin to use for the transform

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive)

  • mval (numpy.ndarray((nm,), dtype = numpy.uint64)) – the m moments for which the transform should be carried out entries must be unique and <= lmax

  • mstart (numpy.ndarray((nm,), dtype = numpy.uint64)) – the (hypothetical) index in the second dimension of alm on which the entry with l=0, m=mval[mi] would be stored, for mi in mval

  • lstride (int) – the index stride in the second dimension of alm between the entries for l and l+1, but the same m.

  • theta (numpy.ndarray((ntheta,), dtype=numpy.float64)) – the colatitudes of the map rings

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

  • mode (str) –

    the transform mode
    ”STANDARD”: standard transform
    ”GRAD_ONLY”: only valid for spin>0, curl a_lm are assumed to be zero
    ”DERIV1”: same as “GRAD_ONLY”, but spin is assumed to be 1 and a_lm are multiplied by sqrt(l*(l+1))

  • theta_interpol (bool) – if the input grid is irregularly spaced in theta, try to accelerate the transform by using an intermediate equidistant theta grid and a 1D NUFFT.

Returns:

the Legendre coefficients. If leg was supplied, this will be the same object.

Return type:

numpy.ndarray((nleg, ntheta, nm), same dtype as alm)

Notes

nleg = 1 if spin == 0 else 2 nalm = 1 if spin == 0 else (2 if mode == “STANDARD” else 1)

ducc0.sht.experimental.alm2leg_deriv1(*, alm: numpy.ndarray, lmax: int, theta: numpy.ndarray, mval: object = None, mstart: object = None, lstride: int = 1, nthreads: int = 1, leg: object = None, theta_interpol: bool = False) numpy.ndarray

Transforms a set of spin-0 spherical harmonic coefficients to Legendre coefficients of the first derivatives with respect to colatiude and longitude, dependent on theta and m.

Parameters:
  • alm (numpy.ndarray((1, x), dtype=numpy.complex64 or numpy.complex128)) – the set of spherical harmonic coefficients. The second dimension must be large enough to accommodate all entries, which are stored according to the parameters lmax, ‘mval`, mstart, and lstride.

  • leg (None or numpy.ndarray((2, ntheta, nm), same dtype as alm)) – output array containing the Legendre coefficients if None, a new suitable array is allocated

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive)

  • mval (numpy.ndarray((nm,), dtype = numpy.uint64)) – the m moments for which the transform should be carried out entries must be unique and <= lmax

  • mstart (numpy.ndarray((nm,), dtype = numpy.uint64)) – the (hypothetical) index in the second dimension of alm on which the entry with l=0, m=mval[mi] would be stored, for mi in mval

  • lstride (int) – the index stride in the second dimension of alm between the entries for l and l+1, but the same m.

  • theta (numpy.ndarray((ntheta,), dtype=numpy.float64)) – the colatitudes of the map rings

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

  • theta_interpol (bool) – if the input grid is irregularly spaced in theta, try to accelerate the transform by using an intermediate equidistant theta grid and a 1D NUFFT.

Returns:

the Legendre coefficients. If leg was supplied, this will be the same object. The first component contains coefficients representing a map of df/dtheta; the second component those of 1./sin(theta) df/dphi.

Return type:

numpy.ndarray((2, ntheta, nm), same dtype as alm)

ducc0.sht.experimental.analysis_2d(*, map: numpy.ndarray, spin: int, lmax: int, geometry: str, mmax: object = None, nthreads: int = 1, alm: object = None, phi0: float = 0.0, mstart: object = None, lstride: int = 1) numpy.ndarray

Transforms one or two 2D maps to spherical harmonic coefficients. This is the inverse operation of synthesis_2D.

Parameters:
  • alm (numpy.ndarray((ncomp, x), dtype=numpy.complex64 or numpy.complex128)) – storage for the spherical harmonic coefficients. The second dimension must be large enough to accommodate all entries, which are stored according to the healpy convention If not supplied, a new array is allocated.

  • map (numpy.ndarray((ncomp, ntheta, nphi), dtype=numpy.float of same accuracy as alm)) – The input map.

  • spin (int >= 0) – the spin to use for the transform. If spin==0, ncomp must be 1, otherwise 2

  • lmax (int >= 0) – the maximum l (and m) moment of the transform (inclusive)

  • mmax (int >= 0 <= lmax) – the maximum m moment of the transform (inclusive).

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • geometry (one of "CC", "F1", "MW", "MWflip", "GL", "DH", "F2") –

    the distribution of rings over the theta range
    • CC: Clenshaw-Curtis, equidistant, first and last ring on poles

    • F1: Fejer’s first rule, equidistant, first and last ring half a ring width from the poles

    • MW: McEwen & Wiaux scheme, equidistant, first ring half a ring width from the north pole, last ring on the south pole

    • MWflip: flipped McEwen & Wiaux scheme, equidistant, first ring on the north pole, last ring half a ring width from the south pole

    • GL: Gauss-Legendre, non-equidistant

    • DH: Driscoll-Healy, equidistant, first ring on north pole, last ring one ring width from the south pole

    • F2: Fejer’s second rule, equidistant, first and last ring one ring width from the poles.

  • phi0 (float) – the azimuth (in radians) of the first pixel in each ring

  • nthreads (int >= 0) – the number of threads to use for the computation. If 0, use as many threads as there are hardware threads available on the system

Returns:

the computed spherical harmonic coefficients If the alm parameter was specified, this is identical to alm.

Return type:

numpy.ndarray((ncomp, x), dtype=numpy.complex64 or numpy.complex128)

Notes

The maximum m moment to which this function can analyze its input map is min(lmax, (nphi-1)//2).

The maximum l moment to which this function can analyze its input map depends on the geometry, and is

  • ntheta-2 for CC

  • ntheta-1 for F1, MW, MWflip, and GL

  • (ntheta-2)//2 for DH

  • (ntheta-1)//2 for F2

For the CC and F1 geometries this limit is considerably higher than the one obtainable by simply applying quadrature weights. This improvement is achieved by temporary upsampling along meridians to apply the weights at a higher resolution.

ducc0.sht.experimental.flm2alm(flm: numpy.ndarray, spin: int, alm: object = None, real: bool = False) numpy.ndarray
ducc0.sht.experimental.get_gridweights(geometry: str, ntheta: int) numpy.ndarray

Returns the quadrature weights for a given grid geometry and number of rings.

Parameters:
  • geometry (one of "CC", "F1", "MW", "MWflip", "GL", "DH", "F2") –

    the distribution of rings over the theta range
    • CC: Clenshaw-Curtis, equidistant, first and last ring on poles

    • F1: Fejer’s first rule, equidistant, first and last ring half a ring width from the poles

    • MW: McEwen & Wiaux scheme, equidistant, first ring half a ring width from the north pole, last ring on the south pole

    • MWflip: flipped McEwen & Wiaux scheme, equidistant, first ring on the north pole, last ring half a ring width from the south pole

    • GL: Gauss-Legendre, non-equidistant

    • DH: Driscoll-Healy, equidistant, first ring on north pole, last ring one ring width from the south pole

    • F2: Fejer’s second rule, equidistant, first and last ring one ring width from the poles.

  • ntheta (int > 0) – number of rings in the grid.

Returns:

the quadrature weights for the individual rings. Please note that these weights need to be divided by the number of pixels per ring to obtain the actual quadrature weights for a particular map.

Return type:

numpy.ndarray((ntheta,), dtype=numpy.float64)

ducc0.sht.experimental.leg2alm(*, leg: numpy.ndarray, lmax: int, theta: numpy.ndarray, spin: int = 0, mval: object = None, mstart: object = None, lstride: int = 1, nthreads: int = 1, alm: object = None, mode: str = 'STANDARD', theta_interpol: bool = False) numpy.ndarray

Transforms a set of Legendre coefficients to spherical harmonic coefficients

Parameters:
  • leg (numpy.ndarray((nleg, ntheta, nm), dtype=numpy.complex64 or numpy.complex128)) –

  • alm (None or numpy.ndarray((nalm, x), same dtype as leg)) – the set of spherical harmonic coefficients. The second dimension must be large enough to accommodate all entries, which are stored according to the parameters lmax, ‘mval`, mstart, and lstride. if None, a new suitable array is allocated

  • spin (int >= 0) – the spin to use for the transform

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive)

  • mval (numpy.ndarray((nm,), dtype = numpy.uint64)) – the m moments for which the transform should be carried out entries must be unique and <= lmax

  • mstart (numpy.ndarray((nm,), dtype = numpy.uint64)) – the (hypothetical) index in the second dimension of alm on which the entry with l=0, m=mval[mi] would be stored, for mi in mval

  • lstride (int) – the index stride in the second dimension of alm between the entries for l and l+1, but the same m.

  • theta (numpy.ndarray((ntheta,), dtype=numpy.float64)) – the colatitudes of the map rings

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

  • mode (str) –

    the transform mode
    ”STANDARD”: standard transform
    ”GRAD_ONLY”: only valid for spin>0, curl a_lm are not computed
    ”DERIV1”: same as “GRAD_ONLY”, but spin is assumed to be 1 and a_lm are multiplied by sqrt(l*(l+1))

  • theta_interpol (bool) – if the input grid is irregularly spaced in theta, try to accelerate the transform by using an intermediate equidistant theta grid and a 1D NUFFT.

Returns:

the Legendre coefficients. if alm was supplied, this will be the same object If newly allocated, the smallest possible second dimensions will be chosen.

Return type:

numpy.ndarray((nalm, *), same dtype as leg)

Notes

nleg = 1 if spin == 0 else 2 nalm = 1 if spin == 0 else (2 if mode == “STANDARD” else 1)

ducc0.sht.experimental.leg2map(*, leg: numpy.ndarray, nphi: numpy.ndarray, phi0: numpy.ndarray, ringstart: numpy.ndarray, pixstride: int = 1, nthreads: int = 1, map: object = None) numpy.ndarray

Transforms one or more sets of Legendre coefficients to maps.

Parameters:
  • leg (numpy.ndarray((ncomp, ntheta, mmax+1), numppy.complex64 or numpy.complex128)) – input array containing the Legendre coefficients. The entries in leg[:,:,m] correspond to quantum number m, i.e. the m values must be stored in ascending order, and complete.

  • map (None or numpy.ndarray((ncomp, x), dtype=numpy.float of same accuracy as leg) – the map pixel data. The second dimension must be large enough to accommodate all pixels, which are stored according to the parameters nphi, ‘ringstart`, and pixstride. if None, a new suitable array is allocated

  • nphi (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – number of pixels in every ring

  • phi0 (numpy.ndarray((ntheta,), dtype=numpy.float64)) – azimuth (in radians) of the first pixel in every ring

  • ringstart (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – the index in the second dimension of map at which the first pixel of every ring is stored

  • pixstride (int) – the index stride in the second dimension of map between two subsequent pixels in a ring

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

Returns:

the map pixel data. If map was supplied, this will be the same object If newly allocated, the smallest possible second dimensions will be chosen.

Return type:

numpy.ndarray((ncomp, x), dtype=numpy.float of same accuracy as leg)

Notes

In contrast to leg2alm and alm2leg the m values are assumed to form a range from 0 to mmax, inclusively.

ducc0.sht.experimental.map2leg(*, map: numpy.ndarray, nphi: numpy.ndarray, phi0: numpy.ndarray, ringstart: numpy.ndarray, mmax: int, pixstride: int = 1, nthreads: int = 1, leg: object = None) numpy.ndarray

Transforms a map or several maps to Legendre coefficients dependent on theta and m.

Parameters:
  • map (numpy.ndarray((ncomp, x), dtype=numpy.float32 or numpy.float64)) – the map pixel data. The second dimension must be large enough to accommodate all pixels, which are stored according to the parameters nphi, ‘ringstart`, and pixstride.

  • leg (None or numpy.ndarray((ncomp, ntheta, mmax+1), dtype=numpy.complex of same accuracy as map)) – output array containing the Legendre coefficients if None, a new suitable array is allocated The entries in leg[:,:,m] correspond to quantum number m, i.e. the m values will be stored in ascending order.

  • nphi (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – number of pixels in every ring

  • phi0 (numpy.ndarray((ntheta,), dtype=numpy.float64)) – azimuth (in radians) of the first pixel in every ring

  • ringstart (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – the index in the second dimension of map at which the first pixel of every ring is stored

  • pixstride (int) – the index stride in the second dimension of map between two subsequent pixels in a ring

  • mmax (int) – the maximum m moment to compute in this transform. If leg is provided, mmax must be equal to leg.shape[2]=1.

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

Returns:

the Legendre coefficients if leg was supplied, this will be the same object

Return type:

numpy.ndarray((ncomp, ntheta, nm), dtype=numpy.complex of same accuracy as map)

Notes

In contrast to leg2alm and alm2leg the m values are assumed to form a range from 0 to mmax, inclusively.

ducc0.sht.experimental.maxium_safe_l(geometry: str, ntheta: int) int

Returns the maximum l moment that can be safely stored (i.e. is guaranteed to be recoverable using analysis_2d) in a map with the specified geometry and number of rings.

Parameters:
  • geometry (one of "CC", "F1", "MW", "MWflip", "GL", "DH", "F2") –

    the distribution of rings over the theta range
    • CC: Clenshaw-Curtis, equidistant, first and last ring on poles

    • F1: Fejer’s first rule, equidistant, first and last ring half a ring width from the poles

    • MW: McEwen & Wiaux scheme, equidistant, first ring half a ring width from the north pole, last ring on the south pole

    • MWflip: flipped McEwen & Wiaux scheme, equidistant, first ring on the north pole, last ring half a ring width from the south pole

    • GL: Gauss-Legendre, non-equidistant

    • DH: Driscoll-Healy, equidistant, first ring on north pole, last ring one ring width from the south pole

    • F2: Fejer’s second rule, equidistant, first and last ring one ring width from the poles.

  • ntheta (int > 0 (> 1 for geometries CC, DH, and F1)) – number of rings in the grid.

Returns:

lmax – The maximum l moment that can be safely stored on the specified grid.

Return type:

int>=0

ducc0.sht.experimental.pseudo_analysis(*, map: numpy.ndarray, theta: numpy.ndarray, lmax: int, mstart: object = None, nphi: numpy.ndarray, phi0: numpy.ndarray, ringstart: numpy.ndarray, spin: int, lstride: int = 1, pixstride: int = 1, nthreads: int = 1, alm: object = None, maxiter: int, epsilon: float, mmax: object = None, theta_interpol: bool = False) object

Tries to extract spherical harmonic coefficients from (sets of) one or two maps by using the iterative LSMR algorithm.

Parameters:
  • alm (None or numpy.ndarray(([ntrans,] ncomp, x), dtype=numpy.complex of same precision as map)) – the set of spherical harmonic coefficients. The last dimension must be large enough to accommodate all entries, which are stored according to the parameters lmax, ‘mmax`, mstart, and lstride. if None, a new suitable array is allocated

  • map (numpy.ndarray(([ntrans,] ncomp, x), dtype=numpy.float32 or numpy.float64) – The last dimension must be large enough to accommodate all pixels, which are stored according to the parameters nphi, ‘ringstart`, and pixstride.

  • theta (numpy.ndarray((ntheta,), dtype=numpy.float64)) – the colatitudes of the map rings

  • nphi (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – number of pixels in every ring

  • phi0 (numpy.ndarray((ntheta,), dtype=numpy.float64)) – azimuth (in radians) of the first pixel in every ring

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • ringstart (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – the index in the last dimension of map at which the first pixel of every ring is stored

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • pixstride (int) – the index stride in the last dimension of map between two subsequent pixels in a ring

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

  • spin (int >= 0) – the spin to use for the transform. If spin==0, ncomp must be 1, otherwise 2

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive).

  • maxiter (int >= 0) – the maximum number of iterations before stopping the algorithm

  • epsilon (float >= 0) – the relative tolerance used as a stopping criterion

  • mmax (int >= 0 <= lmax) – the maximum m moment of the transform (inclusive).

  • theta_interpol (bool) – if the input grid is irregularly spaced in theta, try to accelerate the transform by using an intermediate equidistant theta grid and a 1D NUFFT.

Returns:

  • numpy.ndarray(([ntrans,] ncomp, x), dtype=numpy.complex of same accuracy as map) – the set(s) of spherical harmonic coefficients. If alm was supplied, this will be the same object If newly allocated, the smallest possible last dimension will be chosen.

  • int or list(int) – the reason for stopping the iteration 1: approximate solution to the equation system found 2: approximate least-squares solution found 3: condition number of the equation system too large 7: maximum number of iterations reached

  • int or list(int) – the iteration count(s)

  • float or list(float) – the residual norm, divided by the norm of map

  • float or list(float) – the quality of the least-squares solution

ducc0.sht.experimental.pseudo_analysis_general(*, lmax: int, map: numpy.ndarray, loc: numpy.ndarray, spin: int, nthreads: int, maxiter: int, epsilon: float = 1e-05, sigma_min: float = 1.1, sigma_max: float = 2.6, mstart: object = None, lstride: int = 1, mmax: object = None, alm: object = None) tuple

Tries to extract spherical harmonic coefficients from one or two maps by using the iterative LSMR algorithm.

Parameters:
  • map (numpy.ndarray((ncomp, npix), dtype=numpy.float32 or numpy.float64) – The pixel values at the locations specified by loc.

  • spin (int >= 0) – the spin to use for the transform. If spin==0, ncomp must be 1, otherwise 2

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive).

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • mmax (int >= 0 and <= lmax) – the maximum m moment of the transform (inclusive). If not supplied, mmax is derived from mstart; if that is not supplied either, it is assumed to be equal to lmax.

  • loc (numpy.array((npix, 2), dtype=numpy.float64)) – the locations on the sphere at which the alm should be evaluated. loc[:, 0] contains colatitude values (range [0;pi]), loc[:, 1] contains longitude values (range [0;2pi])

  • epsilon (float >= 0) – the relative tolerance used as a stopping criterion NOTE: for the “epsilon” paraeter of the underlyig NUFFT calls, 0.1*epsilon will be used.

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

  • alm (None or numpy.ndarray((ncomp, x), dtype=complex, same accuracy as map)) – the set(s) of spherical harmonic coefficients. The last dimension must be large enough to accommodate all entries, which are stored according to the healpy convention. If None, a new suitable array is allocated.

  • sigma_min (float) – minimum and maximum allowed oversampling factors for the NUFFT component 1.2 <= sigma_min < sigma_max <= 2.5

  • sigma_max (float) – minimum and maximum allowed oversampling factors for the NUFFT component 1.2 <= sigma_min < sigma_max <= 2.5

  • maxiter (int >= 0) – the maximum number of iterations before stopping the algorithm

Returns:

  • numpy.ndarray((ncomp, x), dtype=complex, same accuracy as map) – the computed spherical harmonic coefficients If the alm parameter was specified, this is identical to alm. If newly allocated, the smallest possible last dimension will be chosen.

  • int – the reason for stopping the iteration 1: approximate solution to the equation system found 2: approximate least-squares solution found 3: condition number of the equation system too large 7: maximum number of iterations reached

  • int – the iteration count(s)

  • float – the residual norm, divided by the norm of map

  • float – the quality of the least-squares solution

ducc0.sht.experimental.synthesis(*, alm: numpy.ndarray, theta: numpy.ndarray, lmax: int, mstart: object = None, nphi: numpy.ndarray, phi0: numpy.ndarray, ringstart: numpy.ndarray, spin: int, lstride: int = 1, pixstride: int = 1, nthreads: int = 1, map: object = None, mmax: object = None, mode: str = 'STANDARD', theta_interpol: bool = False) numpy.ndarray

Transforms (sets of) one or two sets of spherical harmonic coefficients to maps on the sphere.

Parameters:
  • alm (numpy.ndarray(([ntrans,] nalm, x), dtype=numpy.complex64 or numpy.complex128)) – the set(s) of spherical harmonic coefficients. The last dimension must be large enough to accommodate all entries, which are stored according to the parameters lmax, ‘mmax`, mstart, and lstride.

  • map (None or numpy.ndarray(([ntrans,] nmaps, x), dtype=numpy.float of same accuracy as alm) – the map pixel data. The last dimension must be large enough to accommodate all pixels, which are stored according to the parameters nphi, ‘ringstart`, and pixstride. if None, a new suitable array is allocated

  • theta (numpy.ndarray((ntheta,), dtype=numpy.float64)) – the colatitudes of the map rings

  • nphi (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – number of pixels in every ring

  • phi0 (numpy.ndarray((ntheta,), dtype=numpy.float64)) – azimuth (in radians) of the first pixel in every ring

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • ringstart (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – the index in the last dimension of map at which the first pixel of every ring is stored

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • pixstride (int) – the index stride in the last dimension of map between two subsequent pixels in a ring

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

  • spin (int >= 0) – the spin to use for the transform.

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive).

  • mmax (int >= 0 <= lmax) – the maximum m moment of the transform (inclusive).

  • mode (str) –

    the transform mode
    ”STANDARD”: standard transform
    ”GRAD_ONLY”: only valid for spin>0, curl a_lm are assumed to be zero
    ”DERIV1”: same as “GRAD_ONLY”, but spin is assumed to be 1 and a_lm are multiplied by sqrt(l*(l+1))

  • theta_interpol (bool) – if the input grid is irregularly spaced in theta, try to accelerate the transform by using an intermediate equidistant theta grid and a 1D NUFFT.

Returns:

the map pixel data. If map was supplied, this will be the same object If newly allocated, the smallest possible last dimension will be chosen.

Return type:

numpy.ndarray(([ntrans,] nmaps, x), dtype=numpy.float of same accuracy as alm)

Notes

nmaps = 1 if spin == 0 else 2 nalm = 1 if spin == 0 else (2 if mode == “STANDARD” else 1)

ducc0.sht.experimental.synthesis_2d(*, alm: numpy.ndarray, spin: int, lmax: int, geometry: str, ntheta: object = None, nphi: object = None, mmax: object = None, nthreads: int = 1, map: object = None, mode: str = 'STANDARD', phi0: float = 0.0, mstart: object = None, lstride: int = 1) numpy.ndarray

Transforms one or two sets of spherical harmonic coefficients to 2D maps.

Parameters:
  • alm (numpy.ndarray((nalm, x), dtype=numpy.complex64 or numpy.complex128)) – the set of spherical harmonic coefficients. The second dimension must be large enough to accommodate all entries, which are stored according to the healpy convention.

  • map (numpy.ndarray((nmaps, ntheta, nphi), dtype=numpy.float of same accuracy as alm)) – storage for the output map. If not supplied, a new array is allocated.

  • ntheta (int > 0) – dimensions of the output map. If not supplied, map must be supplied. If supplied, and map is also supplied, must match with the array dimensions

  • nphi (int > 0) – dimensions of the output map. If not supplied, map must be supplied. If supplied, and map is also supplied, must match with the array dimensions

  • spin (int >= 0) – the spin to use for the transform.

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive).

  • mmax (int >= 0 <= lmax) – the maximum m moment of the transform (inclusive).

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • geometry (one of "CC", "F1", "MW", "MWflip", "GL", "DH", "F2") –

    the distribution of rings over the theta range
    • CC: Clenshaw-Curtis, equidistant, first and last ring on poles

    • F1: Fejer’s first rule, equidistant, first and last ring half a ring width from the poles

    • MW: McEwen & Wiaux scheme, equidistant, first ring half a ring width from the north pole, last ring on the south pole

    • MWflip: flipped McEwen & Wiaux scheme, equidistant, first ring on the north pole, last ring half a ring width from the south pole

    • GL: Gauss-Legendre, non-equidistant

    • DH: Driscoll-Healy, equidistant, first ring on north pole, last ring one ring width from the south pole

    • F2: Fejer’s second rule, equidistant, first and last ring one ring width from the poles.

  • phi0 (float) – the azimuth (in radians) of the first pixel in each ring

  • nthreads (int >= 0) – the number of threads to use for the computation. If 0, use as many threads as there are hardware threads available on the system

  • mode (str) –

    the transform mode
    ”STANDARD”: standard transform
    ”GRAD_ONLY”: only valid for spin>0, curl a_lm are assumed to be zero
    ”DERIV1”: same as “GRAD_ONLY”, but spin is assumed to be 1 and a_lm are multiplied by sqrt(l*(l+1))

Returns:

the computed map. If the map parameter was specified, this is identical with map.

Return type:

numpy.ndarray((nmaps, ntheta, nphi), dtype=numpy.float of same accuracy as alm)

Notes

nmaps = 1 if spin == 0 else 2 nalm = 1 if spin == 0 else (2 if mode == “STANDARD” else 1)

ducc0.sht.experimental.synthesis_2d_deriv1(*, alm: numpy.ndarray, lmax: int, geometry: str, ntheta: object = None, nphi: object = None, mmax: object = None, nthreads: int = 1, map: object = None, phi0: float = 0.0, mstart: object = None, lstride: int = 1) numpy.ndarray

Transforms a set of spherical harmonic coefficients to two 2D maps containing the derivatives with respect to theta and phi.

Parameters:
  • alm (numpy.ndarray((1, x), dtype=numpy.complex64 or numpy.complex128)) – the set of spherical harmonic coefficients. The second dimension must be large enough to accommodate all entries, which are stored according to the healpy convention.

  • map (numpy.ndarray((2, ntheta, nphi), dtype=numpy.float of same accuracy as alm)) – storage for the output map. If not supplied, a new array is allocated.

  • ntheta (int > 0) – dimensions of the output map. If not supplied, map must be supplied. If supplied, and map is also supplied, must match with the array dimensions

  • nphi (int > 0) – dimensions of the output map. If not supplied, map must be supplied. If supplied, and map is also supplied, must match with the array dimensions

  • lmax (int >= 0) – the maximum l (and m) moment of the transform (inclusive)

  • mmax (int >= 0 <= lmax) – the maximum m moment of the transform (inclusive).

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • geometry (one of "CC", "F1", "MW", "MWflip", "GL", "DH", "F2") –

    the distribution of rings over the theta range
    • CC: Clenshaw-Curtis, equidistant, first and last ring on poles

    • F1: Fejer’s first rule, equidistant, first and last ring half a ring width from the poles

    • MW: McEwen & Wiaux scheme, equidistant, first ring half a ring width from the north pole, last ring on the south pole

    • MWflip: flipped McEwen & Wiaux scheme, equidistant, first ring on the north pole, last ring half a ring width from the south pole

    • GL: Gauss-Legendre, non-equidistant

    • DH: Driscoll-Healy, equidistant, first ring on north pole, last ring one ring width from the south pole

    • F2: Fejer’s second rule, equidistant, first and last ring one ring width from the poles.

  • phi0 (float) – the azimuth (in radians) of the first pixel in each ring

  • nthreads (int >= 0) – the number of threads to use for the computation. If 0, use as many threads as there are hardware threads available on the system

Returns:

the maps containing the derivatives with respect to theta and phi. If the map parameter was specified, this is identical with map.

Return type:

numpy.ndarray((2, ntheta, nphi), dtype=numpy.float of same accuracy as alm)

ducc0.sht.experimental.synthesis_deriv1(*, alm: numpy.ndarray, theta: numpy.ndarray, lmax: int, mstart: object = None, nphi: numpy.ndarray, phi0: numpy.ndarray, ringstart: numpy.ndarray, lstride: int = 1, pixstride: int = 1, nthreads: int = 1, map: object = None, mmax: object = None, theta_interpol: bool = False) numpy.ndarray

Transforms a set (or sets) of spherical harmonic coefficients to two maps containing the derivatives with respect to theta and phi.

Parameters:
  • alm (numpy.ndarray(([ntrans,] 1, x), dtype=numpy.complex64 or numpy.complex128)) – the set(s) of spherical harmonic coefficients. The last dimension must be large enough to accommodate all entries, which are stored according to the parameters lmax, ‘mmax`, mstart, and lstride.

  • map (None or numpy.ndarray(([ntrans,] 2, x), dtype=numpy.float of same accuracy as alm) – the map pixel data. The last dimension must be large enough to accommodate all pixels, which are stored according to the parameters nphi, ‘ringstart`, and pixstride. if None, a new suitable array is allocated

  • theta (numpy.ndarray((ntheta,), dtype=numpy.float64)) – the colatitudes of the map rings

  • nphi (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – number of pixels in every ring

  • phi0 (numpy.ndarray((ntheta,), dtype=numpy.float64)) – azimuth (in radians) of the first pixel in every ring

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • ringstart (numpy.ndarray((ntheta,), dtype=numpy.uint64)) – the index in the last dimension of map at which the first pixel of every ring is stored

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • pixstride (int) – the index stride in the last dimension of map between two subsequent pixels in a ring

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive).

  • mmax (int >= 0 <= lmax) – the maximum m moment of the transform (inclusive).

  • theta_interpol (bool) – if the input grid is irregularly spaced in theta, try to accelerate the transform by using an intermediate equidistant theta grid and a 1D NUFFT.

Returns:

the map pixel data. If map was supplied, this will be the same object If newly allocated, the smallest possible last dimension will be chosen.

Return type:

numpy.ndarray((2, x), dtype=numpy.float of same accuracy as alm)

ducc0.sht.experimental.synthesis_general(*, alm: numpy.ndarray, spin: int, lmax: int, loc: numpy.ndarray, epsilon: float = 1e-05, mstart: object = None, lstride: int = 1, mmax: object = None, nthreads: int = 1, map: object = None, sigma_min: float = 1.1, sigma_max: float = 2.6, mode: str = 'STANDARD', verbose: bool = False) numpy.ndarray

Evaluate a_lm at arbitrary positions on the sphere

Parameters:
  • alm (numpy.ndarray((nalm, x), dtype=numpy.complex64 or numpy.complex128)) – the set(s) of spherical harmonic coefficients. The last dimension must be large enough to accommodate all entries, which are stored according to the healpy convention.

  • spin (int >= 0) – the spin to use for the transform.

  • lmax (int >= 0) – the maximum l moment of the transform (inclusive).

  • mstart (numpy.ndarray((mmax+1,), dtype = numpy.uint64)) – the (hypothetical) index in the last dimension of alm on which the entry with (l=0, m) would be stored. If not supplied, a contiguous storage scheme in the order m=0,1,2,… is assumed.

  • lstride (int) – the index stride in the last dimension of alm between the entries for l and l+1, but the same m.

  • mmax (int >= 0 and <= lmax) – the maximum m moment of the transform (inclusive). If not supplied, mmax is derived from mstart; if that is not supplied either, it is assumed to be equal to lmax.

  • loc (numpy.array((npix, 2), dtype=numpy.float64)) – the locations on the sphere at which the alm should be evaluated. loc[:, 0] contains colatitude values (range [0;pi]), loc[:, 1] contains longitude values (range [0;2pi])

  • epsilon (float) – desired accuracy for single precision inputs, this must be >1e-6, for double precision it must be >2e-13

  • nthreads (int >= 0) – the number of threads to use for the computation if 0, use as many threads as there are hardware threads available on the system

  • map (None or numpy.ndarray((nmaps, npix), dtype=numpy.float of same accuracy as alm) – the map pixel data. If None, a new suitable array is allocated.

  • sigma_min (float) – minimum and maximum allowed oversampling factors for the NUFFT component 1.2 <= sigma_min < sigma_max <= 2.5

  • sigma_max (float) – minimum and maximum allowed oversampling factors for the NUFFT component 1.2 <= sigma_min < sigma_max <= 2.5

  • mode (str) –

    the transform mode
    ”STANDARD”: standard transform
    ”GRAD_ONLY”: only valid for spin>0, curl a_lm are assumed to be zero
    ”DERIV1”: same as “GRAD_ONLY”, but spin is assumed to be 1 and a_lm are multiplied by sqrt(l*(l+1))

Returns:

the pixel values at the locations specified by loc. If the map parameter was specified, this is identical with map.

Return type:

numpy.ndarray((nmaps, npix), dtype=numpy.float of same accuracy as alm

Notes

nmaps = 1 if spin == 0 else 2 nalm = 1 if spin == 0 else (2 if mode == “STANDARD” else 1)