Power Spectrum
Power spectrum components $C_{\ell}$ are encoded as Vector{T}. Healpix.jl implements functions to perform sht operations on power spectra, e.g. to obtain a map or a set of 'Alm', as well as writing/reading a power spectrum from a FITS file.
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).
Loading and saving power spectrum components
Healpix.jl implements functions to read/write the components $C_{\ell}$ from/to a FITS files.
Healpix.readClFromFITS
— FunctionreadClFromFITS{T <: Real}(f::CFITSIO.FITSFile, t::Type{T}; col_num = 2) -> Vector{T}
readClFromFITS{T <: Real}(fileName::String, t::Type{T}; col_num = 2) -> Vector{T}
Read a set of C_ℓ coefficients from a FITS file.
Healpix.writeClToFITS
— FunctionwriteClToFITS(f::CFITSIO.FITSFile, Cl::Vector{T}) where {T <: Real}
writeClToFITS(fileName, Cl::Vector{T}; overwrite = true) where {T <: Real}
Write a set of C_ℓ coefficients to a FITS file.
Converting different power spectrum representation
It's often useful to represent, especially for plotting it, the power spectrum in the form of $D_{\ell}$. Healpix.jl implements a couple of functions to convert a power spectrum from/to such a representation.
Healpix.cl2dl
— Functioncl2dl(cl::AbstractVector{T}, lmin::Integer) where {T <: Real}
Convert a set of $C_{\ell}$ to $D_{\ell}$ power spectrum, where $D_{\ell} = \ell (\ell + 1) C_{\ell} / 2\pi$. The first components are set to zero if not present.
Arguments:
cl::AbstractVector{T}
: Array of C_ℓ componentslmin::Integer
: minimum l in the representation of the C_ℓ power spectrum
Returns:
Vector{T}
: Array of D_ℓ power spectrum components
Healpix.dl2cl
— Functiondl2cl(dl::AbstractVector{T}, lmin::Integer) where {T <: Real}
Convert a set of $D_{\ell}$ to $C_{\ell}$ power spectrum, where $C_{\ell} = 2\pi D_{\ell} / \ell (\ell + 1)$. The first components are set to zero if not present. The monopole component is set to zero in any case to avoid Inf values.
Arguments:
dl::AbstractVector{T}
: Array of D_ℓ componentslmin::Integer
: minimum l in the representation of the Dℓ power spectrum
Returns:
Vector{T}
: Array of C_ℓ power spectrum components
Synthesizing harmonic coefficients from power spectrum
Generate a Alm
instance with a random set of $a_{\ell m}$ coefficients. Each harmonic coefficient $a_{\ell m}$ is a realization of a gaussian distribution with zero mean and $C_{\ell}$ variance.
Healpix.synalm!
— Functionsynalm!(cl::Vector{T}, alm::Alm{ComplexF64, Vector{ComplexF64}}, rng::AbstractRNG) where {T <: Real}
synalm!(cl::Vector{T}, alm::Alm{ComplexF64, Vector{ComplexF64}}) where {T <: Real}
Generate a set of $a_{\ell m}$ from a given power spectra $C_{\ell}$. The output is written into the Alm
object passed in input.
Arguments:
cl::AbstractVector{T}
: The array representing the power spectrum components $C_{\ell}$,
starting from $\ell = 0$.
alm::Alm{Complex{T}}
: The array representing the spherical harmonics coefficients $a_{\ell m}$
we want to write the result into.
rng::AbstractRNG
: (optional) the RNG to be used for generating the $a_{\ell m}$. It allows
to set the seed beforehand guaranteeing the reproducibility of the process.
Healpix.synalm
— Functionsynalm(cl::Vector{T}, lmax::Integer, mmax::Integer, rng::AbstractRNG) where {T <: Real}
synalm(cl::Vector{T}, lmax::Integer, mmax::Integer) where {T <: Real}
synalm(cl::Vector{T}, lmax::Integer, rng::AbstractRNG) where {T <: Real}
synalm(cl::Vector{T}, lmax::Integer) where {T <: Real}
synalm(cl::Vector{T}, rng::AbstractRNG) where {T <: Real}
synalm(cl::Vector{T}) where {T <: Real}
Generate a set of $a_{\ell m}$ from a given power spectra $C_{\ell}$. The output is written into a new Alm
object of given lmax.
Arguments:
cl::AbstractVector{T}
: The array representing the power spectrum components $C_{\ell}$,
starting from $\ell = 0$.
lmax::Integer
: the maximum $ℓ$ coefficient, will default tolength(cl)-1
if not specified.mmax::Integer
: the maximum $m$ coefficient, will default tolmax
if not specified.rng::AbstractRNG
: (optional) the RNG to be used for generating the $a_{\ell m}$. It allows
to set the seed beforehand guaranteeing the reproducibility of the process.
Generating a map from power spectrum
Synthesize a set of Alm
through 'synalm' and generates a map from it through 'alm2map'.
Healpix.synfast!
— Functionsynfast!(cl::Vector{T}, map::HealpixMap{T, RingOrder}, lmax::Integer, rng::AbstractRNG) where {T <: Real}
synfast!(cl::Vector{T}, map::HealpixMap{T, RingOrder}, lmax::Integer) where {T <: Real}
synfast!(cl::Vector{T}, map::HealpixMap{T, RingOrder}, rng::AbstractRNG) where {T <: Real}
synfast!(cl::Vector{T}, map::HealpixMap{T, RingOrder}) where {T <: Real}
Generate a map from a given power spectra $C_{\ell}$. The result is saved into the HealpixMap
passed in input.
Arguments:
cl::AbstractVector{T}
: The array representing the power spectrum components $C_{\ell}$.map::HealpixMap{T, RingOrder}
: the map that will contain the result.lmax::Integer
: the maximum $ℓ$ coefficient, will default tolength(cl)-1
if not specified.rng::AbstractRNG
: (optional) the RNG to be used for generating the $a_{\ell m}$. It allows
to set the seed beforehand guaranteeing the reproducibility of the process.
Healpix.synfast
— Functionsynfast(cl::Vector{T}, nside::Integer, lmax::Integer, rng::AbstractRNG) where {T <: Real}
synfast(cl::Vector{T}, nside::Integer, lmax::Integer) where {T <: Real}
synfast(cl::Vector{T}, nside::Integer, rng::AbstractRNG) where {T <: Real}
synfast(cl::Vector{T}, nside::Integer) where {T <: Real}
Generate a HealpixMap
with given Nside, from a given power spectra $C_{\ell}$.
Arguments:
cl::AbstractVector{T}
: The array representing the power spectrum components $C_{\ell}$.nside::Integer
: nside of the map that will contain the result.lmax::Integer
: the maximum $ℓ$ coefficient, will default tolength(cl)
-1 if not specified.rng::AbstractRNG
: (optional) the RNG to be used for generating the $a_{\ell m}$. It allows
to set the seed beforehand guaranteeing the reproducibility of the process.
Computing the power spectrum from a map
Compute the (cross-) power spectrum of one (or two) 'HealpixMap'.
Healpix.anafast
— Functionanafast(map::HealpixMap{Float64, RingOrder, AA}; lmax=nothing, mmax=nothing, niter::Integer = 3) where {T <: Real,AA <: AbstractArray{T,1}} -> Vector{Float64}
anafast(map₁::HealpixMap{Float64, RingOrder, AA}, map₂::HealpixMap{Float64, RingOrder, AA}; lmax=nothing, mmax=nothing, niter::Integer = 3) where {T <: Real,AA <: AbstractArray{T,1}} -> Vector{Float64}
Computes the power spectrum of a Healpix map, or the cross-spectrum between two maps if map2
is given. No removal of monopole or dipole is performed. The input maps must be in ring-ordering.
Arguments:
map₁::HealpixMap{Float64, RingOrder, AA}
: the spherical harmonic coefficients of the first fieldmap₂::HealpixMap{Float64, RingOrder, AA}
: the spherical harmonic coefficients of the second field
Returns:
Array{T}
containing $C_{\ell}$, with the first element referring to ℓ=0.