Library

Library

Pixel Functions

LibHealpix.nside2npixFunction.
nside2npix(nside)

Compute the number of pixels in a Healpix map with the given value of nside.

Arguments:

  • nside - the Healpix resolution parameter

Usage:

julia> nside2npix(4)
192

julia> nside2npix(256)
786432

See Also: npix2nside, nside2nring

source
LibHealpix.npix2nsideFunction.
npix2nside(npix)

Compute the value of the nside parameter for a Healpix map with the given number of pixels.

Arguments:

  • npix - the number of pixels in the map

Usage:

julia> npix2nside(192)
4

julia> npix2nside(786432)
256

See Also: nside2npix, nside2nring

source
nside2nring(nside)

Compute the number of equal latitude rings in the Healpix map with the given value of nside.

Arguments:

  • nside - the Healpix resolution parameter

Usage:

julia> nside2nring(4)
15

julia> nside2nring(256)
1023

See Also: nside2npix, npix2nside

source
LibHealpix.ang2vecFunction.
ang2vec(theta, phi)

Compute the Cartesian unit vector to the spherical coordinates $(θ, ϕ)$.

Arguments:

  • theta - the inclination angle $θ$

  • phi - the azimuthal angle $ϕ$

Usage:

julia> ang2vec(0, 0)
3-element LibHealpix.UnitVector:
 0.0
 0.0
 1.0

julia> ang2vec(π/2, π/2)
3-element LibHealpix.UnitVector:
 6.12323e-17
 1.0
 6.12323e-17

See Also: vec2ang

source
LibHealpix.vec2angFunction.
vec2ang(vec)

Compute the spherical coordinates $(θ, ϕ)$ from the given unit vector.

Arguments:

  • vec - the input Cartesian unit vector

Usage:

julia> vec2ang([1, 0, 0])
(1.5707963267948966, 0.0)

julia> vec2ang([0, 1, 0])
(1.5707963267948966, 1.5707963267948966)

julia> vec2ang([0, 0, 1])
(0.0, 0.0)

See Also: ang2vec

source
LibHealpix.nest2ringFunction.
nest2ring(nside, ipix)

Convert the given pixel index from the nested to the ring indexing scheme.

Arguments:

  • nside - the Healpix resolution parameter

  • ipix - the pixel index (nested scheme)

Usage:

julia> nest2ring(256, 1)
391809

julia> nest2ring(256, 2)
390785

See Also: ring2nest

source
LibHealpix.ring2nestFunction.
ring2nest(nside, ipix)

Convert the given pixel index from the ring to the nested indexing scheme.

Arguments:

  • nside - the Healpix resolution parameter

  • ipix - the pixel index (ring scheme)

Usage:

julia> ring2nest(256, 1)
65536

julia> ring2nest(256, 2)
131072

See Also: nest2ring

source
ang2pix_nest(nside, theta, phi)

Compute the pixel index (in the nested scheme) that contains the point on the sphere given by the spherical coordinates $(θ, ϕ)$.

Arguments:

  • nside - the Healpix resolution parameter

  • theta - the inclination angle $θ$

  • phi - the azimuthal angle $ϕ$

Usage:

julia> ang2pix_nest(256, 0, 0)
65536

julia> ang2pix_nest(256, π/2, π/2)
354987

See Also: ang2pix_ring, pix2ang_nest, pix2ang_ring

source
ang2pix_ring(nside, theta, phi)

Compute the pixel index (in the ring scheme) that contains the point on the sphere given by the spherical coordinates $(θ, ϕ)$.

Arguments:

  • nside - the Healpix resolution parameter

  • theta - the inclination angle $θ$

  • phi - the azimuthal angle $ϕ$

Usage:

julia> ang2pix_ring(256, 0, 0)
1

julia> ang2pix_ring(256, π/2, π/2)
392961

See Also: ang2pix_nest, pix2ang_nest, pix2ang_ring

source
pix2ang_nest(nside, ipix)

Compute the spherical coordinates $(θ, ϕ)$ corresponding to the given pixel center.

Arguments:

  • nside - the Healpix resolution parameter

  • ipix - the pixel index (nested scheme)

Usage:

julia> pix2ang_nest(256, 1)
(1.5681921571847817, 0.7853981633974483)

julia> pix2ang_nest(256, 2)
(1.5655879699137618, 0.7884661249732196)

See Also: pix2ang_ring, ang2pix_nest, ang2pix_ring

source
pix2ang_ring(nside, ipix)

Compute the spherical coordinates $(θ, ϕ)$ corresponding to the given pixel center.

Arguments:

  • nside - the Healpix resolution parameter

  • ipix - the pixel index (ring scheme)

Usage:

julia> pix2ang_ring(256, 1)
(0.0031894411211228764, 0.7853981633974483)

julia> pix2ang_ring(256, 2)
(0.0031894411211228764, 2.356194490192345)

See Also: pix2ang_nest, ang2pix_nest, ang2pix_ring

source
vec2pix_nest(nside, vec)

Compute the pixel index (in the nested scheme) that contains the point on the sphere given by the Cartesian unit vector.

Arguments:

  • nside - the Healpix resolution parameter

  • vec - the input Cartesian unit vector

Usage:

julia> vec2pix_nest(256, [1, 0, 0])
289451

julia> vec2pix_nest(256, [0, 1, 0])
354987

julia> vec2pix_nest(256, [0, 0, 1])
65536

See Also: vec2pix_ring, pix2vec_nest, pix2vec_ring

source
vec2pix_ring(nside, vec)

Compute the pixel index (in the ring scheme) that contains the point on the sphere given by the Cartesian unit vector.

Arguments:

  • nside - the Healpix resolution parameter

  • vec - the input Cartesian unit vector

Usage:

julia> vec2pix_ring(256, [1, 0, 0])
392705

julia> vec2pix_ring(256, [0, 1, 0])
392961

julia> vec2pix_ring(256, [0, 0, 1])
1

See Also: vec2pix_nest, pix2vec_nest, pix2vec_ring

source
pix2vec_nest(nside, ipix)

Compute the Cartesian unit vector corresponding to the given pixel center.

Arguments:

  • nside - the Healpix resolution parameter

  • ipix - the pixel index (nested scheme)

Usage:

julia> pix2vec_nest(256, 1)
3-element LibHealpix.UnitVector:
 0.707104
 0.707104
 0.00260417

julia> pix2vec_nest(256, 2)
3-element LibHealpix.UnitVector:
 0.704925
 0.709263
 0.00520833

See Also: pix2vec_ring, vec2pix_nest, vec2pix_ring

source
pix2vec_ring(nside, ipix)

Compute the Cartesian unit vector corresponding to the given pixel center.

Arguments:

  • nside - the Healpix resolution parameter

  • ipix - the pixel index (ring scheme)

Usage:

julia> pix2vec_ring(256, 1)
3-element LibHealpix.UnitVector:
 0.00225527
 0.00225527
 0.999995

julia> pix2vec_ring(256, 2)
3-element LibHealpix.UnitVector:
 -0.00225527
  0.00225527
  0.999995

See Also: pix2vec_nest, vec2pix_nest, vec2pix_ring

source

Healpix Maps

abstract type HealpixMap{T<:Number} <: AbstractVector{T}

This abstract type represents a Healpix equal-area pixelization of the sphere.

Subtypes:

  • RingHealpixMap - a HealpixMap where pixels are ordered along rings of constant latitude. This ordering should be used for performing spherical harmonic transforms.

  • NestHealpixMap - a HealpixMap where nearby pixels also tend to be nearby in memory.

source
struct RingHealpixMap{T<:Number} <: HealpixMap{T}

This type represents a Healpix equal-area pixelization of the sphere where pixels are ordered along rings of constant latitude.

Fields:

  • nside - the Healpix resolution parameter

  • pixels - the list of pixel values

Constructors:

RingHealpixMap(T, nside)

Construct a RingHealpixMap with the element type T and resolution parameter nside. All of the pixels will be set to zero initially.

RingHealpixMap(pixels)

Construct a RingHealpixMap with the given list of pixel values. The resolution parameter nside will be inferred from the number of pixels. However a LibHealpixException will be thrown if given an invalid number of pixels.

RingHealpixMap(nside, pixels)

Construct a RingHealpixMap with the given resolution parameter nside and initial list of pixel values. This constructor is cheaper than RingHealpixMap(pixels) if the correct value of nside is already known.

Usage:

julia> map = RingHealpixMap(Float64, 256)
       for idx = 1:length(map)
           map[idx] = randn()
       end
       map + map == 2map
true

See also: HealpixMap, NestHealpixMap, Alm

source
struct NestHealpixMap{T<:Number} <: HealpixMap{T}

This type represents a Healpix equal-area pixelization of the sphere where nearby pixels also tend to be nearby in memory.

Fields:

  • nside - the Healpix resolution parameter

  • pixels - the list of pixel values

Constructors:

NestHealpixMap(T, nside)

Construct a NestHealpixMap with the element type T and resolution parameter nside. All of the pixels will be set to zero initially.

NestHealpixMap(pixels)

Construct a NestHealpixMap with the given list of pixel values. The resolution parameter nside will be inferred from the number of pixels. However a LibHealpixException will be thrown if given an invalid number of pixels.

NestHealpixMap(nside, pixels)

Construct a NestHealpixMap with the given resolution parameter nside and initial list of pixel values. This constructor is cheaper than NestHealpixMap(pixels) if the correct value of nside is already known.

Usage:

julia> map = NestHealpixMap(Float64, 256)
       for idx = 1:length(map)
           map[idx] = randn()
       end
       map + map == 2map
true

See also: HealpixMap, RingHealpixMap, Alm

source
writehealpix(filename, map)

Write the HealpixMap to disk as a FITS image.

Arguments:

  • filename - the name of the output file (eg. "/path/to/healpix.fits")

  • map - the Healpix map to write

Keyword Arguments:

  • coordsys - the coordinate system of the map (one of "G" galactic, "E" ecliptic, or "C" celestial)

  • replace - if set to true, the output file will be automatically overwritten if it exists

See also: readhealpix

source
readhealpix(filename)

Read a HealpixMap (stored as a FITS image) from disk.

Arguments:

  • filename - the name of the input file (eg. "/path/to/healpix.fits")

See also: writehealpix

source
LibHealpix.ang2pixFunction.
ang2pix(map, theta, phi)

Compute the pixel index that contains the point on the sphere given by the spherical coordinates $(θ, ϕ)$.

Arguments:

  • map - the input Healpix map

  • theta - the inclination angle $θ$ (in radians)

  • phi - the azimuthal angle $ϕ$ (in radians)

Usage:

julia> ang2pix(RingHealpixMap(Float64, 256), π/2, π/2)
392961

julia> ang2pix(NestHealpixMap(Float64, 256), π/2, π/2)
354987

See Also: pix2ang, ang2pix_nest, ang2pix_ring

source
LibHealpix.pix2angFunction.
pix2ang(map, ipix)

Compute the spherical coordinates $(θ, ϕ)$ corresponding to the given pixel center.

Arguments:

  • map - the input Healpix map

  • ipix - the pixel index

Usage:

julia> pix2ang(RingHealpixMap(Float64, 256), 1)
(0.0031894411211228764, 0.7853981633974483)

julia> pix2ang(NestHealpixMap(Float64, 256), 1)
(1.5681921571847817, 0.7853981633974483)

See Also: ang2pix, pix2ang_nest, pix2ang_ring

source
LibHealpix.vec2pixFunction.
vec2pix(map, vec)

Compute the pixel index that contains the point on the sphere given by the Cartesian unit vector.

Arguments:

  • map - the input Healpix map

  • vec - the input Cartesian unit vector

Usage:

julia> vec2pix(RingHealpixMap(Float64, 256), [0, 0, 1])
1

julia> vec2pix(NestHealpixMap(Float64, 256), [0, 0, 1])
65536

See Also: pix2vec, vec2pix_nest, vec2pix_ring

source
LibHealpix.pix2vecFunction.
pix2vec(map, ipix)

Compute the Cartesian unit vector corresponding to the given pixel center.

Arguments:

  • map - the input Healpix map

  • ipix - the pixel index (nested scheme)

Usage:

julia> pix2vec(RingHealpixMap(Float64, 256), 1)
3-element LibHealpix.UnitVector:
 0.00225527
 0.00225527
 0.999995

julia> pix2vec(NestHealpixMap(Float64, 256), 1)
3-element LibHealpix.UnitVector:
 0.707104
 0.707104
 0.00260417

See Also: vec2pix, pix2vec_nest, pix2vec_ring

source
LibHealpix.UNSEENFunction.
LibHealpix.UNSEEN()

Get the sentinal value indicating a blind or masked pixel.

source
LibHealpix.interpolate(map, theta, phi)

Linearly interpolate the Healpix map at the given spherical coordinates $(θ, ϕ)$.

Arguments:

  • map - the input Healpix map

  • theta - the inclination angle $θ$ (in radians)

  • phi - the azimuthal angle $ϕ$ (in radians)

Usage:

julia> healpixmap = RingHealpixMap(Float64, 256)
       for idx = 1:length(healpixmap)
           healpixmap[idx] = idx
       end
       LibHealpix.interpolate(healpixmap, 0, 0)
2.5

See Also: ang2pix

source
LibHealpix.query_discFunction.
query_disc(nside, ordering, theta, phi, radius; inclusive=true)
query_disc(map, theta, phi, radius; inclusive=true)

Return a list of all pixels contained within a circular disc of the given radius.

Arguments:

  • nside - the Healpix resolution parameter

  • ordering - the ordering of the Healpix map (either LibHealpix.ring or LibHealpix.nest

  • theta - the inclination angle $θ$ (in radians)

  • phi - the azimuthal angle $ϕ$ (in radians)

  • radius - the radius of the disc (in radians)

  • map - the input Healpix map (nside and ordering will be inferred from the map)

Keyword Arguments:

  • inclusive - if set to `true pixels partially contained within the disc will be included, otherwise they are excluded

Usage:

julia> query_disc(512, LibHealpix.ring, 0, 0, deg2rad(10/60), inclusive=false)
4-element Array{Int32,1}:
 1
 2
 3
 4

julia> query_disc(512, LibHealpix.ring, 0, 0, deg2rad(10/60), inclusive=true) |> length
24
source

Spherical Harmonic Coefficients

LibHealpix.AlmType.
struct Alm{T<:Number} <: AbstractVector{T}

This type holds a vector of spherical harmonic coefficients.

Fields:

  • lmax - the maximum value for the $l$ quantum number

  • mmax - the maximum value for the $m$ quantum number (note that $m ≤ l$)

  • coefficients - the list of spherical harmonic coefficients

Constructors:

Alm(T, lmax, mmax)

Construct an Alm object that will store all spherical harmonic coefficients with element type T, $l ≤ lₘₐₓ$, and $m ≤ mₘₐₓ$. All of the coefficients will be initialized to zero.

Alm(lmax, mmax, coefficients)

Construct an Alm object with the given list of initial coefficients corresponding to $l ≤ lₘₐₓ$, and $m ≤ mₘₐₓ$. A LibHealpixException will be thrown if too many or too few coefficients are provided.

Usage:

julia> alm = Alm(Complex128, 10, 10)
       for (l, m) in lm(alm)
           @lm alm[l, m] = l * m
       end
       @lm(alm[10, 5]) == 50
true
Note

The lm function is used to iterate over the spherical harmonic quantum numbers $l$ and $m$.

Note

The @lm macro is used to index into an Alm object when given the spherical harmonic quantum numbers $l$ and $m$.

See also: RingHealpixMap, NestHealpixMap, lm, @lm

source
LibHealpix.lmFunction.
lm(lmax, mmax)
lm(alm)

Construct an interator for iterating over all possible values of the spherical harmonic quantum numbers $l ≤ lₘₐₓ$ and $m ≤ mₘₐₓ$.

Arguments:

  • lmax - the maximum value of $l$

  • mmax - the maximum value of $m$

  • alm - if an Alm object is provided, lmax and mmax will be inferred from the corresponding fields

Usage:

julia> for (l, m) in lm(2, 1)
           @show l, m
       end
(l, m) = (0, 0)
(l, m) = (1, 0)
(l, m) = (2, 0)
(l, m) = (1, 1)
(l, m) = (2, 1)

See Also: Alm, @lm

source
LibHealpix.@lmMacro.
@lm

This macro is used to index an Alm object when given the values for quantum numbers $l$ and $m$.

Usage

julia> alm = Alm(Int, 2, 1)
       for (l, m) in lm(alm)
           @lm alm[l, m] = l + m
       end

julia> @lm alm[1, 1]
2

julia> @lm alm[1, :] # all coefficients with l == 1
2-element Array{Int64,1}:
 1
 2

julia> @lm alm[:, 1] # all coefficients with m == 1
2-element Array{Int64,1}:
 2
 3

Background

Alm implements the AbstractVector interface which allows the type to be used in place of a standard Vector in many cases. This generally makes sense because Alm is simply a wrapper around a standard Vector.

However, one consequence of being an AbstractVector is that the two-element getindex function already has a meaning and therefore alm[l, m] cannot be used to mean "give me the coefficient corresponding to the quantum numbers $l$ and $m$". Instead @lm alm[l, m] calls a separate function that does give you the coefficient for $l$ and $m$.

See Also: Alm, lm

source

Spherical Harmonic Transforms

LibHealpix.map2almFunction.
map2alm(map, lmax, mmax)

Compute the spherical harmonic coefficients of the given Healpix map by means of a spherical harmonic transform.

Arguments:

  • map - the input Healpix map (must be ring ordered)

  • lmax - the maximum value for the $l$ quantum number

  • mmax - the maximum value for the $m$ quantum number

Keyword Arguments:

  • iterations - the number of iterations to perform

Note

Set iterations to something greater than 0 if more precision is required.

Usage:

julia> map = RingHealpixMap(Float64, 4)
       map[:] = 1
       alm = map2alm(map, 1, 1)
       @lm alm[0, 0]
3.5449077018110318 + 0.0im

See Also: alm2map

source
LibHealpix.alm2mapFunction.
alm2map(alm, nside)

Compute the Healpix map corresponding to the given spherical harmonic coefficients by means of an inverse spherical harmonic transform.

Arguments:

  • alm - the input list of spherical harmonic coefficients

  • nside - the resolution of the output Healpix map

Usage:

julia> alm = Alm(Complex128, 1, 1)
       @lm alm[0, 0] = 1
       map = alm2map(alm, 1)
12-element LibHealpix.RingHealpixMap{Float64}:
 0.282095
 0.282095
 0.282095
 0.282095
 0.282095
 0.282095
 0.282095
 0.282095
 0.282095
 0.282095
 0.282095
 0.282095

See Also: map2alm

source

Visualization

LibHealpix.mollweideFunction.
mollweide(map, size=(512, 1024))

Create a Mollweide projected image of the given Healpix map.

Note

The image values will be set to zero outside of the projection area.

Arguments:

  • map - the input Healpix map

  • size - the size of the output image

Usage:

julia> map = RingHealpixMap(Int, 1)
       map[:] = 1
       mollweide(map, (10, 20))
10×20 Array{Int64,2}:
 0  0  0  0  0  0  1  1  1  1  1  1  1  1  0  0  0  0  0  0
 0  0  0  1  1  1  1  1  1  1  1  1  1  1  1  1  1  0  0  0
 0  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  0
 1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1
 1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1
 1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1
 1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1
 0  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  1  0
 0  0  0  1  1  1  1  1  1  1  1  1  1  1  1  1  1  0  0  0
 0  0  0  0  0  0  1  1  1  1  1  1  1  1  0  0  0  0  0  0
source