Spherical harmonics

Starting from version 2.4, Healpix.jl implements generalized Fourier transformations through the libsharp library, to convert a map from its pixel-space representation to its decomposition in spherical harmonics. This has multiple applications, the most relevant being the analysis of Cosmic Microwave Background maps and the efficient computation of convolution operators.

Everything revolves around the Alm type, which encodes a set of spherical harmonics and is thus conceptually equivalent to the concept of a HealpixMap, only living in the harmonic space:

Healpix.AlmType

An array of harmonic coefficients (a_ℓm).

The type T is used for the value of each harmonic coefficient, and it must be a Number (one should however only use complex types for this). The type AA is used to store the array of coefficients; a typical choice is Vector.

A Alm type contains the following fields:

  • alm: the array of harmonic coefficients
  • lmax: the maximum value for $ℓ$
  • mmax: the maximum value for $m$
  • tval: maximum number of $m$ coefficients for the maximum $ℓ$

The $a_{\ell m}$ are stored by $m$: if $\ell_{max}$ is 16, the first 16 elements are $m=0$, $\ell=0-16$, then the following 15 elements are $m=1$, $\ell=1-16$, then $m=2$, $\ell=2-16$ and so on until the last element, the 153th, is $m=16$, $\ell=16$.

source

In the general case, the number of coefficients in a spherical harmonic expansion is infinite. For obvious reasons, Healpix.jl only allows to store band-limited expansions. The function numberOfAlms returns the number of floating-point numbers used to store the expansion, as a function of the maximum value for $\ell$ and $m$.

Healpix.numberOfAlmsFunction
numberOfAlms(lmax::Integer, mmax::Integer) -> Integer
numberOfAlms(lmax::Integer) -> Integer

Return the size of the array of complex numbers needed to store the a_ℓm coefficients in the range of ℓ and m specified by lmax and mmax. If mmax is not specified, it is assumed to be equal to lmax. If lmax and mmax are inconsistent or negative, a DomainError exception is thrown.

source

Converting between pixel space and harmonic space

Healpix.jl implements the four functions alm2map, map2alm, alm2map!, and map2alm! to convert a map from a pixel-space representation to the harmonic space and vice-versa. The functions ending with ! are mutating functions, which means that they assume that the result must be saved in a preallocated variable; they are space- and time-efficient and should be used when you want your code to be performant, or when you plan to apply the same operation several times (e.g., in a Monte Carlo simulation). The synthesis operation (alm2map) is generally referred to with the matrix operator $\mathrm{Y}$, while his inverse (map2alm) with $\mathrm{Y}^{-1}$.

Healpix.jl also implements the two adjoint functions adjoint_alm2map! and adjoint_map2alm!, represented by $\mathrm{Y}^{\mathrm{T}}$ and $(\mathrm{Y}^{-1})^\mathrm{T}$ respectively. While the synthesis operator on a general scalar field $f(\theta, \phi)$ can be defined through an exact summation as $f(\theta, \phi) = \mathrm{Y} \, a_{\ell m} \quad \text{where} \quad f(\theta, \phi) = \sum_{\ell=0}^{\infty} \sum_{m=-\ell}^{\ell} a_{\ell m} Y_{\ell m} (\theta, \phi)$. The analysis operator is defined through an integral operator as $a_{\ell m} = \mathrm{Y}^{-1} f(\theta, \phi) \quad \text{where} \quad a_{\ell m} = \int_0^{2\pi} \int_0^\pi Y^*_{\ell m}(\theta, \phi)\, f(\theta, \phi) \sin\theta \, d\theta \,d\phi$. Though, in the real case wherein maps are pixelized, the latter ends up being approximated through a summation over the pixels. Here is where the adjoint of the synthesis operator, $\mathrm{Y}^{\mathrm{T}}$, comes into play. It is defined through: $ \mathrm{Y}^{\mathrm{T}} f(\theta, \phi) \equiv \sum{i = 1}^{N{\mathrm{pix}}} Y^*{\ell m,\, i} \, fi,$ which is an exact operation. Note that the latter does not give directly the $a_{\ell m}$ coefficients, since $\mathrm{Y}^{-1} \simeq \mathrm{W}\, \mathrm{Y}^{\mathrm{T}}$, where $\mathrm{W}$ is a diagonal matrix whose non-zero elements are approximately constant and equal to $4 \pi / N_{\mathrm{pix}}$, depending on the map pixelization. The latter realtion is also useful to obtain the adjoint of the analysis operator: $(\mathrm{Y}^{-1})^\mathrm{T} = \mathrm{W}^{\mathrm{T}}\,\mathrm{Y} = \mathrm{W}\,\mathrm{Y}$.

Here is an example:

using Random

# Ensure reproducibility by using a fixed seed
Random.seed!(1234)

nside = 8
m = HealpixMap{Float32,RingOrder}(nside)

# Initialize the pixels to random values in the 0…1 range
for i in 1:length(m)
    m[i] = rand(Float32)
end

alm = map2alm(m)

# Go back to pixel space
newm = alm2map(alm, nside)
768-element HealpixMap{Float64, RingOrder, Vector{Float64}}:
  0.16526318742679227
  0.36187810449645147
 -0.0013815388023585007
  0.4008476600279501
  0.2530007233310308
  0.7785956979193497
  0.43071187199759764
  0.937945224190703
  0.7554174116102923
  0.603273537633899
  ⋮
  0.4434725912890835
  0.4226379392252659
 -0.2052063921524328
  0.41659016756519374
  0.6860000214208772
  0.7211351759306384
  0.6707478409770004
  0.47791307158759316
  0.7109101053740822

The variable newm is a map that is close enough to m, yet it is not exactly the same because of the approximations done by both map2alm and alm2map.

Healpix.map2alm!Function
map2alm!(map::HealpixMap{Float64, RingOrder, Array{Float64, 1}}, alm::Alm{ComplexF64, Array{ComplexF64, 1}}; niter::Integer=3)
map2alm!(map::PolarizedHealpixMap{Float64, RingOrder, Array{Float64, 1}}, alm::Array{Alm{ComplexF64, Array{ComplexF64, 1}},1};
    niter::Integer=3)

This function performs a spherical harmonic transform on the map and places the results in the passed alm object. This function requires types derived from Float64, since it is done in-place.

Arguments

  • map: the map that must be decomposed in spherical harmonics. It can either be a HealpixMap{Float64, RingOrder} type (scalar map) or a PolarizedHealpixMap{Float64, RingOrder} type (polarized map).

  • alm::Alm{ComplexF64, Array{ComplexF64, 1}}: the spherical harmonic coefficients to be written to.

Keywords

  • niter::Integer: number of iterations of SHTs to perform, to enhance accuracy
source
Healpix.map2almFunction
map2alm(map::HealpixMap{Float64, RingOrder, AA};
    lmax=nothing, mmax=nothing, niter::Integer=3)
map2alm(m::HealpixMap{T, RingOrder, AA}; lmax=nothing, mmax=nothing,
    niter::Integer=3) where {T <: Real, AA <: AbstractArray{T, 1} }

Compute the spherical harmonic coefficients of a map. To enhance precision, more iterations of the transforms can be performed by passing a nonzero niter. The underlying SHT library libsharp performs all calculations using Cdouble types, so all inputs are converted to types based on Float64.

Arguments

  • map: the map to decompose in spherical harmonics. It can either be a HealpixMap{T, RingOrder, AA} type (scalar map) or a PolarizedHealpixMap{T, RingOrder, AA} type (polarized map).

Keywords

  • lmax::Integer: the maximum ℓ coefficient, will default to 3*nside-1 if not specified.

  • mmax::Integer: the maximum m coefficient

  • niter::Integer: number of SHT iterations, to enhance precision. Defaults to 3

Returns

  • Alm{ComplexF64, Array{ComplexF64, 1}}: the spherical harmonic coefficients corresponding to the map
source
Healpix.alm2map!Function
alm2map!(alm::Alm{ComplexF64, Array{ComplexF64, 1}}, map::HealpixMap{Float64, RingOrder, Array{Float64, 1}})
alm2map!(alm::Array{Alm{ComplexF64, Array{ComplexF64, 1}},1}, map::PolarizedHealpixMap{Float64, RingOrder, Array{Float64, 1}})

This function performs a spherical harmonic transform on the map and places the results in the passed alm object. This function requires types derived from Float64, since it is done in-place.

Arguments

  • alm::Alm{ComplexF64, Array{ComplexF64, 1}}: the spherical harmonic coefficients to perform the spherical harmonic transform on.

  • map: the map that will contain the result. It can either be a HealpixMap{Float64, RingOrder, Array{Float64, 1}} type (scalar map) or a PolarizedHealpixMap{Float64, RingOrder, Array{Float64, 1}} (polarized map).

source
Healpix.alm2mapFunction
alm2map(alm::Alm{ComplexF64, Array{Float64, 1}}, nside::Integer)
alm2map(alm::Alm{T, Array{T, 1}}, nside::Integer) where T
alm2map(alm::Array{Alm{ComplexF64, Array{ComplexF64, 1}},1}, nside::Integer)
alm2map(alms::Array{Alm{T, Array{T, 1}},1}, nside::Integer) where T

Compute a map from spherical harmonic coefficients. The underlying SHT library libsharp performs all calculations in Cdouble, so all inputs are converted to types based on Float64.

Arguments

  • alm: the spherical harmonic coefficients to transform. If of type Alm{T, Array{T, 1}}, we assume a spin-0 spherical harmonic transform. If an array of Alm is passed, we assume that the components correspond to T, E, and B coefficients.

Keywords

  • nside::Integer: Healpix resolution parameter

Returns

  • HealpixMap{Float64, RingOrder, Array{Float64, 1}} or PolarizedHealpixMap{Float64, RingOrder, Array{Float64, 1}} depending on if the input alm is of type Alm{T, Array{T, 1}} or Array{Alm{T, Array{T, 1}}} respectively.
source
Healpix.adjoint_alm2map!Function
adjoint_alm2map!(map::HealpixMap{Float64, RingOrder, Array{Float64, 1}}, alm::Alm{ComplexF64, Array{ComplexF64, 1}})
adjoint_alm2map!(map::PolarizedHealpixMap{Float64, RingOrder, Array{Float64, 1}}, alm::Array{Alm{ComplexF64, Array{ComplexF64, 1}},1})

This function performs a spherical harmonic transform Yᵀ on the map and places the results in the passed alm object. This function requires types derived from Float64, since it is done in-place.

Arguments

  • map: the map that must be decomposed in spherical harmonics. It can either be a HealpixMap{Float64, RingOrder} type (scalar map) or a PolarizedHealpixMap{Float64, RingOrder} type (polarized map).

  • alm::Alm{ComplexF64, Array{ComplexF64, 1}}: the spherical harmonic coefficients to be written to.

source
Healpix.adjoint_map2alm!Function
adjoint_map2alm!(alm::Alm{ComplexF64, Array{ComplexF64, 1}}, map::HealpixMap{Float64, RingOrder, Array{Float64, 1}})
adjoint_map2alm!(alm::Array{Alm{ComplexF64, Array{ComplexF64, 1}},1}, map::PolarizedHealpixMap{Float64, RingOrder, Array{Float64, 1}})

This function performs the spherical harmonic transform (Y^-1)^t = W^t Y = W Y on the map and places the results in the passed alm object. This function requires types derived from Float64, since it is done in-place.

Arguments

  • alm::Alm{ComplexF64, Array{ComplexF64, 1}}: the spherical harmonic coefficients to perform the spherical harmonic transform on.

  • map: the map that will contain the result. It can either be a HealpixMap{Float64, RingOrder, Array{Float64, 1}} type (scalar map) or a PolarizedHealpixMap{Float64, RingOrder, Array{Float64, 1}} (polarized map).

source

From harmonic coefficients to the power spectrum

You can use the function alm2cl to convert a set of $a_{\ell m}$ coefficients into the components $C_\ell$ of the power spectrum. The pixelization also induces a transfer function, which can be obtained from pixwin. A simple Gaussian beam window function in the asymptotic small-beam limit can be computed with gaussbeam.

Healpix.alm2clFunction
alm2cl(alm::Alm{Complex{T}}) where {T <: Number} -> Vector{T}
alm2cl(alm₁::Alm{Complex{T}}, alm₂::Alm{Complex{T}}) where {T <: Number} -> Vector{T}

Compute $C_{\ell}$ from the spherical harmonic coefficients of one or two fields.

Arguments

  • alm₁::Alm{Complex{T}}: the spherical harmonic coefficients of the first field
  • alm₂::Alm{Complex{T}}: the spherical harmonic coefficients of the second field

Returns

  • Array{T} containing $C_{\ell}$, with the first element referring to ℓ=0.
source
Healpix.pixwinFunction
pixwin(nside; pol=false)

Arguments:

  • nside: HEALPix resolution parameter

Keywords

  • pol::Bool=false: if true, also return polarization pixel window

Returns:

  • Vector{Float64} pixel window function. If pol=true, returns a Tuple of the temperature and polarization pixel windows.

Examples

julia> pixwin(4)
17-element Vector{Float64}:
 1.0000000000020606
 0.9942340766588788
 ⋮
 0.4222841034207188
source
Healpix.gaussbeamFunction
gaussbeam(fwhm::T, lmax::Int; pol=false) where T

Compute the Gaussian beam window function $B_{\ell}$ given the FWHM of the beam in radians, where $C_{\ell, \mathrm{measured}} = B_{\ell}^2 C_{\ell}$. This beam is valid in the limit of $\sigma^2 \ll 0$, which is the case for all high-resolution CMB experiments.

Arguments

  • fwhm::T: FWHM of the Gaussian beam in radians
  • lmax::Int: maximum multipole ℓ
  • pol=false: if false, returns the spin-0 beam for i.e. intensity. if true, returns the spin-2 beam

Returns

  • Array{T,1} containing $B_{\ell}$, with the first element referring to ℓ=0.
source

Algebraic operations in harmonic space

Healpix.jl provides overloads of the Base functions +, -, *, /, as well as LinearAlgebra.dot, allowing to carry out these fundamental operations element-wise in harmonic space directly.

For example, an element-wise sum between two Alm objects can be performed as follows:

using Healpix

#just two constant Alm objects
myalm1 = Healpix.Alm(5,5, ones(ComplexF64, Healpix.numberOfAlms(5)))
myalm1 = Healpix.Alm(5,5, ones(ComplexF64, Healpix.numberOfAlms(5)))

alm_sum = myalm1 + myalm2 #each element will be = 2 + 0im

Multiplying or dividing a set of Alm by a generic function of $\ell$ or a constant

The operators * and / can be used to multiply or divide an Alm by an $\ell$-dependent generic function $f_\ell$ (or just a constant, of type Number).

In this case a new instance of Alm type will be returned. To perform a more efficient in-place operation refer to almxfl!, as shown in this brief example:

using Healpix

#just two constant Alm objects
myalm = Healpix.Alm(5,5, ones(ComplexF64, Healpix.numberOfAlms(5)))
myf_l = ones(Healpix.numberOfAlms(5)) .* 2

#will return a new object:
myalm*myf_l
myalm/myf_l

#will overwrite myalm:
almxfl!(myalm, myf_l)
almxfl!(myalm, 1.0 ./ myf_l) #division

In either case the call to such operator consists in a shortcut to the function almxfl.

Healpix.almxflFunction
almxfl(alm::Alm{Complex{T}}, fl::AbstractVector{T}) where {T <: Number} -> Alm{T}

Multiply an aℓm by a vector bℓ representing an ℓ-dependent function, without changing the a_ℓm passed in input.

ARGUMENTS

  • alms::Alm{Complex{T}}: The Alm object containing the spherical harmonics coefficients
  • fl::AbstractVector{T}: The array containing the factors fℓ to be multiplied by aℓm

#RETURNS

  • Alm{Complex{T}}: The result of aℓm * fℓ.
source
Healpix.almxfl!Function
almxfl!(alm::Alm{Complex{T}}, fl::AA) where {T <: Number,AA <: AbstractArray{T,1}}

Multiply IN-PLACE an aℓm by a vector bℓ representing an ℓ-dependent function.

ARGUMENTS

  • alms::Alm{Complex{T}}: The Alm object containing the spherical harmonics coefficients
  • fl::AbstractVector{T}: The array containing the factors fℓ to be multiplied by aℓm
source
Base.:+Function

+(alm₁::Alm{Complex{T}}, alm₂::Alm{Complex{T}}) where {T <: Number}

Perform the element-wise sum in a_ℓm space.
A new `Alm` object is returned.
source
Base.:-Function

-(alm₁::Alm{Complex{T}}, alm₂::Alm{Complex{T}}) where {T <: Number}

Perform the element-wise subtraction in a_ℓm space.
A new `Alm` object is returned.
source
Base.:*Function

*(alm₁::Alm{Complex{T}}, alm₂::Alm{Complex{T}}) where {T <: Number}

Perform the element-wise product in a_ℓm space.
A new `Alm` object is returned.
source

*(alm₁::Alm{Complex{T}}, fl::AA) where {T <: Number, AA <: AbstractArray{T,1}}

Perform the product of an `Alm` object by a function of ℓ in a_ℓm space.
Note: this consists in a shortcut of [`almxfl`](@ref),
therefore a new `Alm` object is returned.
source

*(alm₁::Alm{Complex{T}}, fl::AbstractVector{T}) where {T <: Number}

Perform the element-wise product of an `Alm` object by a constant in a_ℓm space.
source
Base.:/Function

/(alm₁::Alm{Complex{T}}, alm₂::Alm{Complex{T}}) where {T <: Number}

Perform an element-wise division in a_ℓm space between two `Alm`s.
A new `Alm` object is returned.
source

/(alm₁::Alm{Complex{T}}, fl::AA) where {T <: Number,AA <: AbstractArray{T,1}}

Perform an element-wise division by a function of ℓ in a_ℓm space.
A new `Alm` object is returned.
source

/(alm₁::Alm{Complex{T}}, c::Number) where {T <: Number}

Perform an element-wise division by a constant in a_ℓm space.
A new `Alm` object is returned.
source

Dot product

Healpix.jl implements an overload of the operator LinearAlgebra.dot (along with its shortcut ) to perform a dot product directly in harmonic space.

using Healpix
myalm = Healpix.Alm(5,5, ones(ComplexF64, Healpix.numberOfAlms(5)))

dot_res = myalm ⋅ myalm #equivalent to dot(myalm, myalm)
Missing docstring.

Missing docstring for LinearAlgebra.dot. Check Documenter's build log for details.

Loading and saving harmonic coefficients

Healpix.readAlmFromFITSFunction
readAlmFromFITS{T <: Complex}(f::CFITSIO.FITSFile, t::Type{T}) -> Alm{T}
readAlmFromFITS{T <: Complex}(fileName::String, t::Type{T}) -> Alm{T}

Read a set of a_ℓm coefficients from a FITS file. If the code fails, CFITSIO will raise an exception. (Refer to the CFITSIO library for more information.)

source
Healpix.writeAlmToFITSFunction
writeAlmToFITS(f::CFITSIO.FITSFile, alm::Alm{Complex{T}}) where {T <: Number}
writeAlmToFITS(fileName, alm::Alm{Complex{T}}; overwrite = true) where {T <: Number}

Write a set of a_ℓm coefficients into a FITS file. If the code fails, CFITSIO will raise an exception. (Refer to the CFITSIO library for more information.) In the fits file the alms are written with explicit index scheme, $\mathrm{index} = \ell^2 + \ell + m + 1$, possibly out of order (check almExplicitIndex).

source

Alm Indexing

You can use almExplicitIndex to compute the so-called explicit indexing. It is exploited for instance in readAlmFromFITS and writeAlmToFITS.

Healpix.almExplicitIndexFunction
almExplicitIndex(lmax) -> Vector{Int}
almExplicitIndex(lmax, mmax) -> Vector{Int}
almExplicitIndex(alm::Alm{T}) where {T} -> Vector{Int}

Compute the explicit index scheme, i.e. $\mathrm{index} = \ell^2 + \ell + m + 1$ up to a certain $ℓ$ and $m$ if specified, or taken from the Alm passed. If not passed, mmax is defaulted to lmax. If lmax and mmax are inconsistent or negative, a DomainError exception is thrown.

source

The following functions can be used, in an analogous way as eachindex, in the case of arrays, to obtain sets of indexes or values of $\ell$ and $m$.

On the same line as eachindex, these can be very useful when implementing for-cycles over Alm objects. Here is an example of how to exploit each_ell_m to print explicitly the major-m ordering of a set of complex-stored Alm:

using Random

# Initialize a random set of alm
alm = Alm(4,4, randn(ComplexF64, numberOfAlms(4,4)))

#print a_lm values knowing each corresponding l and m values
function print_alm(alm)
    i=1
    for (l,m) in each_ell_m(alm)
        a_lm = alm.alm[i]
        print("ℓ = $l, |m| = $m: a_ℓm = $a_lm \n")
        i+=1
    end
end

print_alm(alm)
ℓ = 0, |m| = 0: a_ℓm = -0.05890341708007019 - 0.001451699516140343im
ℓ = 1, |m| = 0: a_ℓm = 0.15348960130527967 - 0.3680136602915076im
ℓ = 2, |m| = 0: a_ℓm = -0.1404181071120885 + 0.2104651891571884im
ℓ = 3, |m| = 0: a_ℓm = -0.610648678062473 - 0.847789396785335im
ℓ = 4, |m| = 0: a_ℓm = -0.7766174857568857 + 0.6342153110699746im
ℓ = 1, |m| = 1: a_ℓm = 0.9480101650016622 - 1.136381772622842im
ℓ = 2, |m| = 1: a_ℓm = -0.04247218587309066 + 0.0861167580589787im
ℓ = 3, |m| = 1: a_ℓm = -0.7852494354333247 - 0.2729763545696586im
ℓ = 4, |m| = 1: a_ℓm = 0.1364808621806462 + 1.2762183868287824im
ℓ = 2, |m| = 2: a_ℓm = 0.24067100771247776 - 0.5963430511551708im
ℓ = 3, |m| = 2: a_ℓm = 0.6379917026290011 - 1.2602205519112877im
ℓ = 4, |m| = 2: a_ℓm = 1.3189206031085603 - 0.3449339018723129im
ℓ = 3, |m| = 3: a_ℓm = 0.24697743879940517 - 0.2562873731354161im
ℓ = 4, |m| = 3: a_ℓm = 0.8078786506105151 + 0.7934122673017472im
ℓ = 4, |m| = 4: a_ℓm = -0.10265442319188665 + 0.7512937490460659im
Healpix.each_ellFunction

eachell(alm::Alm{Complex{T}}, m::Integer) where {T <: Number} -> Vector{Int} eachell(alm::Alm{Complex{T}}, ms::AbstractArray{I, 1}) where {T <: Number, I <: Integer} -> Vector{Int}

Returns an array of all the allowed ℓ values in `alm` for the given `m`.
source
Healpix.each_ell_idxFunction

eachellidx(alm::Alm{Complex{T}}, m::Integer) where {T <: Number} -> Vector{Int} eachellidx(alm::Alm{Complex{T}}, ms::AbstractArray{I, 1}) where {T <: Number, I <: Integer} -> Vector{Int}

Returns an array of the indexes of the harmonic coefficients in `alm` corresponding
to all the ℓ values for the given m value(s).
source
Healpix.each_mFunction

eachm(alm::Alm{Complex{T}}, l::Integer) where {T <: Number} -> Vector{Int} eachm(alm::Alm{Complex{T}}, ls::AbstractArray{I, 1}) where {T <: Number, I <: Integer} -> Vector{Int}

Returns an array containing all the allowed m values in `alm` for the given ℓ value(s).
source
Healpix.each_m_idxFunction

eachmidx(alm::Alm{Complex{T}}, l::Integer) where {T <: Number} -> Vector{Int} eachmidx(alm::Alm{Complex{T}}, ls::AbstractArray{I, 1}) where {T <: Number, I <: Integer} -> Vector{Int}

Returns an array of the indexes of the harmonic coefficients in `alm` corresponding
to all the allowed m values for the given ℓ value(s).
source
Healpix.each_ell_mFunction

eachellm(alm::Alm{Complex{T}}) where {T <: Number} -> Vector{Int}

Returns an array of tuples `(l, m)` of all the ℓ and m values of `alm` in
m-major order (the same order as how the harmonic coefficients are stored in `Alm` objects).
source

Full Pixel Weights

The default map2alm uses iteration to obtain an accurate transform. One can instead apply a pixel weight to compute an accurate transform in a single pass, like quadrature. The easiest way to the pixel weight files is to run

git clone --depth 1 https://github.com/healpy/healpy-data

These weights are in a compressed format that is read with readFullWeights and multiplied into a map with applyFullWeights!.

nside = 32
compressed_weights = Healpix.readFullWeights(
    "healpix_full_weights_nside_$(lpad(nside,4,'0')).fits")
m = Healpix.HealpixMap{Float64,Healpix.RingOrder}(ones(Healpix.nside2npix(nside)))
Healpix.applyFullWeights!(m, compressed_weights)
alm = Healpix.map2alm(m; niter=0)

The subsequent map2alm only needs niter=0.

Healpix.readFullWeightsFunction
readFullWeights(filename::String)

Read the set of pixel weights used to compute the generalized Fourier transform of a map.

These weights are usually precomputed; you can download the ones available in the Healpy repository using the following command:

git clone --depth 1 https://github.com/healpy/healpy-data

Arguments:

  • filename::String: filename of the full pixel weights

Returns:

  • Vector{Float64}: contains the compressed pixel weights
source
Healpix.applyFullWeights!Function
applyFullWeights!(m::HealpixMap{T, RingOrder}, [wgt::Vector{T}]) where T

Apply a pixel weighting to a map for more accurate SHTs. Note that this only helps for lmax<=1.5*Nside. If this is not the case, the pixel weights may do more harm than good.

Pixel weights are automatically downloaded if not specified.

Arguments:

  • m::HealpixMap{T, RingOrder}: map to modify
  • wgt::Vector{T} (optional): compressed pixel weights. If not specified, this routine will look for weights in artifacts.
source
applyFullWeights!(m::PolarizedHealpixMap{T, RingOrder}, [wgt::Vector{T}]) where T

Apply a pixel weighting to a polarized map for more accurate SHTs.

Arguments:

  • m::PolarizedHealpixMap{T, RingOrder}: map to modify
  • wgt::Vector{T} (optional): compressed pixel weights. If not specified, an artifact will be sought.
source