healpy.sphtfunc.harmonic_ud_grade#
- healpy.sphtfunc.harmonic_ud_grade(map_in, nside_out, nside_in=None, lmax=None, mmax=None, iter=None, pol=True, pixwin=True, fwhm_in=0, fwhm_out=0, beam_window_in=None, beam_window_out=None, use_weights=False, datapath=None, use_pixel_weights=True, dtype=None, input_type='map')#
Change map NSIDE using spherical-harmonic transforms.
The input map is analyzed into \(a_{\ell m}\) coefficients, optionally corrected for pixel-window and beam transfer functions following the effective beam transfer ratio (Planck Collaboration 2015 X, arXiv:1502.01588), and synthesized at
nside_out:\[a^{\rm out}_{\ell m} = \frac{p^{\rm out}_\ell}{p^{\rm in}_\ell} \frac{b^{\rm out}_\ell}{b^{\rm in}_\ell} \, a^{\rm in}_{\ell m}\]where \(p_\ell\) is the HEALPix pixel window function and \(b_\ell\) is a beam transfer function (Gaussian or custom).
The transfer is applied as a single-step ratio, so there is no intermediate deconvolution that would amplify high-\(\ell\) noise.
- Parameters:
- map_inarray-like, shape (Npix,) or (n, Npix)
Input map(s), in RING ordering. When
input_type='alm', this is interpreted as \(a_{\ell m}\) coefficients instead (seeinput_typebelow).- nside_outint
Desired NSIDE of the output map(s).
- nside_inint, optional
NSIDE of the input. Required when
input_type='alm'(cannot be inferred from \(a_{\ell m}\)). Ignored wheninput_type='map'(inferred from the input map).- lmaxint, optional
Maximum multipole to retain. If None, defaults to
min(3 * nside_out - 1, 3 * nside_in - 1).- mmaxint, optional
Maximum m of the alm. Default:
lmax. Wheninput_type='alm', the input \(a_{\ell m}\) are assumed to follow the standard healpy conventionmmax_in == lmax_in(as produced bysynalm,map2alm, and friends).- iterint, optional
Number of map2alm iterations. If None, defaults to 0 when using pixel weights with
lmax <= 1.5 * nside_inand 3 otherwise. This follows the HEALPix guidance that pixel weights alone are sufficient without iteration only in that lower-bandlimit regime.- polbool, optional
If True, treat 1- or 3-component input as polarized transforms (TQU/TEB conventions). If False, transform each map independently as spin-0.
- pixwinbool, optional
If True (default), deconvolve the input pixel window \(p^{\rm in}_\ell\) and apply the output pixel window \(p^{\rm out}_\ell\). For polarized maps the temperature and polarization windows are handled separately.
- fwhm_infloat, optional
FWHM in radians of a Gaussian beam to deconvolve from the input \(a_{\ell m}\). Default:
0(no input beam).- fwhm_outfloat or None, optional
FWHM in radians of a Gaussian beam to apply to the output. Default:
0(no output beam — the \(a_{\ell m}\) are simply truncated atlmaxand synthesized). PassNoneto use the effective resolution beam: a Gaussian whose FWHM is proportional to the output pixel size, equivalent to:fwhm_out = hp.effective_resolution_fwhm(nside_out)
This is the scaling adopted by the Planck Collaboration to define the effective beam at each HEALPix resolution (160 arcmin at NSIDE 64, halving with each NSIDE doubling; see arXiv:1502.01588). Applying this beam suppresses Gibbs ringing around sharp features at the new pixel scale and is the recommended choice for science-grade maps of diffuse signals (CMB, Galactic foregrounds). Call
effective_resolution_fwhm()to look up the exact FWHM before downgrading.- beam_window_inarray-like or None, optional
Custom input beam transfer function to deconvolve, overriding
fwhm_in. Follows the format returned bygauss_beam:1D array of shape
(lmax+1,)for temperature-only.2D array of shape
(lmax+1, N)withN >= 2for polarized transforms, where column 0 is the temperature beam and column 1 is the polarization beam (E/B).
When provided,
fwhm_inis ignored. Default:None.- beam_window_outarray-like or None, optional
Custom output beam transfer function to apply, overriding
fwhm_out. Same format asbeam_window_in. When provided,fwhm_outis ignored. Default:None.- use_weightsbool, optional
If True, use HEALPix ring weights in
map2alm. Ring weights improve SHT accuracy at modest cost but are coarser than per-pixel weights. To enable ring weights you must also passuse_pixel_weights=False(the two weight schemes are mutually exclusive). Default:False— prefer the defaultuse_pixel_weights=Trueunless you specifically need ring-weighted behavior for compatibility with another tool.- datapathstr, optional
Directory where to find the HEALPix pixel-weight FITS files (
full_weights/healpix_full_weights_nside_NNNN.fits). IfNone(default), the files are downloaded automatically and cached via astropy.- use_pixel_weightsbool, optional
If True (default), use per-pixel
map2almweights for the most accurate SHT. Unlike plainmap2alm(which silently falls back if the weights file is missing),harmonic_ud_graderaises an exception when pixel weights are unavailable. Passuse_pixel_weights=Falseto disable this behavior explicitly and fall back to unweighted quadrature (or ring weights, ifuse_weights=True).- dtypedtype, optional
If provided, cast output map to this dtype.
- input_type{‘map’, ‘alm’}, optional
Whether
map_inis a pixel-space map ('map', default) or \(a_{\ell m}\) coefficients ('alm'). When'alm', themap2almstep is skipped: the transfer function is applied directly to the input \(a_{\ell m}\) and synthesized atnside_out. This is useful when you already have \(a_{\ell m}\) from a previous SHT and want to avoid a redundant forward transform. Wheninput_type='alm',nside_inmust be provided explicitly, and theiter,use_weights,use_pixel_weights, anddatapathparameters are ignored.
- Returns:
- map_outndarray
Output map(s) at
nside_out.
See also
ud_gradeChange resolution by pixel averaging (faster, but can introduce aliasing for band-limited signals).
gauss_beamCompute Gaussian beam transfer functions for use with
beam_window_in/beam_window_out.
Notes
When to use this function:
harmonic_ud_gradeis the recommended way to change a HEALPix map’s resolution when the signal is diffuse and well-described by a power spectrum (e.g. CMB, Galactic foregrounds). It avoids the aliasing and ringing artifacts thatud_grade(pixel averaging) can introduce for band-limited signals. For maps dominated by point sources or sharp features (e.g. source catalogs, masks),ud_grademay be more appropriate — see the comparison tutorial linked below.Quick-start examples:
Downgrade a CMB temperature map from NSIDE 256 to 64 with the effective resolution beam (recommended for diffuse signals). The output map will be smoothed to 160.0 arcmin (
effective_resolution_fwhm(64)):m_out = hp.harmonic_ud_grade(m_in, nside_out=64, fwhm_out=None)
Check the beam width before calling:
fwhm = hp.effective_resolution_fwhm(64, arcmin=True) # 160.0 arcmin
Downgrade with plain bandlimit truncation (no output smoothing):
m_out = hp.harmonic_ud_grade(m_in, nside_out=64)
Downgrade with explicit input beam and a custom output beam:
m_out = hp.harmonic_ud_grade( m_in, nside_out=64, fwhm_in=np.radians(30), fwhm_out=np.radians(60), )
Reconvolve a map from a 10-arcmin beam to a 30-arcmin beam at the same NSIDE (no resolution change):
m_out = hp.harmonic_ud_grade( m_in, nside_out=nside, fwhm_in=np.radians(10/60), fwhm_out=np.radians(30/60), )
Pass pre-computed \(a_{\ell m}\) to skip the forward SHT:
m_out = hp.harmonic_ud_grade( alm, nside_out=64, nside_in=256, input_type='alm', fwhm_out=None, )
Best practices:
Always specify ``fwhm_in`` when your input map has been smoothed with a known beam. If you don’t, the output will carry the input beam’s attenuation into the new resolution without correction, which can bias power-spectrum estimates.
Use the effective resolution beam (
fwhm_out=None) for science-grade maps of diffuse signals (CMB, Galactic foregrounds). This applies the Planck-scaled beam that suppresses Gibbs ringing at the new pixel scale (Planck Collaboration 2015 X). Plain bandlimit truncation (fwhm_out=0, the default) is appropriate when you need the sharpest possible output or when the downstream analysis handles ringing separately.Use ``pixwin=False`` to skip pixel-window correction entirely. This is appropriate for noise maps (where deconvolving the pixel window would incorrectly amplify high-\(\ell\) noise) or when you want only beam corrections without pixel-window effects.
Use ``input_type=’alm’`` when you already have \(a_{\ell m}\) from a previous analysis step (e.g. component separation, map-making). This avoids a redundant
map2almround-trip and is both faster and more accurate. Wheninput_type='alm',nside_inmust also be provided so that the correct input pixel window can be applied.Pixel weights and iteration: The default
use_pixel_weights=Truewith automatic iteration provides the most accurate \(a_{\ell m}\) for band-limited signals. To match a ring-weighted reference implementation instead, passuse_pixel_weights=False, use_weights=True.
Reconvolution at the same NSIDE: When
nside_out == nside_in, the function performs beam reconvolution — deconvolving the input beam and applying the output beam while keeping the same pixelization. This is similar to the HEALPix process_alm /alteralmfacility.Tutorials:
harmonic_ud_grade vs ud_grade comparison: demonstrates aliasing suppression, noise handling, and Gibbs-ringing tradeoffs with visual examples.
harmonic_ud_grade vs skytools change_resolution: side-by-side API comparison with
skytools.change_resolution, including custom beam transfer functions, pixel-window handling, andinput_type='alm'.