Index

Pixell.EnmapType

Map type, contains an AbstractArray and a WCS object, but behaves like the AbstractArray it contains for array operations. It only implements the subset of Base.Array operations which are common on maps. You should work with the data directly using enmap_instance.data if you need additional Array functions.

source
Pixell.dplanckFunction

The derivative of the planck spectrum with respect to temperature, evaluated at frequencies f and temperature T, in units of Jy/sr/K.

A blackbody has intensity $I = 2hf^3/c^2/(\exp{hf/kT}-1) = V/(\exp{x}-1)$ with $V = 2hf^3/c^2$, $x = hf/kT$.

\[dI/dx = -V/(\exp{x}-1)^2 * \exp{x}\]

\[\begin{aligned} dI/dT &= dI/dx * dx/dT \\ &= 2hf^3/c^2/(\exp{x}-1)^2*\exp{x} * hf/k / T^2 \\ &= 2*h^2*f^4/c^2/k/T^2 * \exp{x}/(\exp{x}-1)^2 \\ &= 2*x^4 * k^3*T^2/(h^2*c^2) * \exp{x}/(\exp{x}-1)^2 \\ &= \dots /(4* \sinh(x/2)^2) \end{aligned}\]

source
Pixell.fullsky_geometryMethod
fullsky_geometry([W=CarClenshawCurtis], res; shape = nothing, dims = ())

Generates a full-sky geometry.

Arguments:

  • proj=CarClenshawCurtis(): [optional] projection
  • res: resolution in radians. Passing a Number produces a square pixel.

Passing a tuple with (ΔRA, ΔDEC) produces a rectangular pixel.

Keywords

  • shape::NTuple=nothing: shape of the map. If not specified, will be computed.
  • dims::NTuple=(): additional dimensions to append to the shape, such as (3,) for IQU

to generate a map with (nx, ny, 3).

Returns:

  • shape::Tuple, wcs::W: a tuple containing the shape of the map and the WCS

Examples

julia> shape, wcs = fullsky_geometry(deg2rad(30/60))  # 30 arcmin pixel
((720, 361), WCSTransform(naxis=2,cdelt=[-0.5, 0.5],crval=[0.25, 0.0],crpix=[360.5, 181.0]))
source
Pixell.laxes_cylMethod
laxes_cyl(shape, wcs)

Get multipole axes for a cylindrical pixelization, with a compromise choice for handling the spherical geometry that changes the angle between meridians in different declinations.

source
Pixell.pix2skyFunction
pix2sky(m::Enmap, pixcoords)

Convert 1-indexed pixels to sky coordinates. While typically WCS can specify custom units, Pixell uses radians throughout.

Arguments:

  • m::Enmap: the map that provides a coordinate system
  • pixcoords: pixcoords should be a 2-d array where "pixcoords[:, i]" is the i-th set of coordinates, or a 1-d array representing a single set of coordinates.

Returns:

  • Array: same shape as pixcoords

Examples

julia> shape, wcs = fullsky_geometry(deg2rad(1))
julia> pix2sky(shape, wcs, [1.0, 1.0])
2-element Vector{Float64}:
 -3.141592653589793
 -1.5707963267948966
source
Pixell.pix2sky!Method
pix2sky!(m::Enmap, pixcoords, skycoords)

Convert 1-indexed pixels to sky coordinates, in-place. While typically WCS can specify custom units, Pixell uses radians throughout.

Arguments:

  • m::Enmap: the map that provides a coordinate system
  • pixcoords: pixel coordinates should be a 2-d array where "pixcoords[:, i]" is the i-th set of coordinates, or a 1-d array representing a single set of coordinates.
  • skycoords: output array for sky coordinates, must be same same as pixcoords

Returns:

  • Array: skycoords

Examples

julia> shape, wcs = fullsky_geometry(deg2rad(1))
       pixcoords =  100 .* rand(2,4096 * 2)
       skycoords =  similar(pixcoords)

julia> pix2sky!(shape, wcs, pixcoords, skycoords)
source
Pixell.pix2skyMethod
pix2sky(m::Enmap{T,N,AA,<:AbstractCARWCS}, ra_pixel, dec_pixel)

Compute the sky position of a single position on the sky.

Only implemented for CAR projections, so the input map is of type Enmap{T,N,AA,<:AbstractCARWCS}. This takes pixel indices for RA and DEC, and returns a tuple containing the corresponding RA and DEC.

Examples

julia> shape, wcs = fullsky_geometry(deg2rad(1))
julia> pix2sky(shape, wcs, 30.0, 80.0)
(151.0, -11.0)
source
Pixell.pixareamapMethod

Generate a similar Enmap whose pixel values are the areas of the pixels in steradians.

source
Pixell.rewindMethod
rewind(angles, period=2π, ref_angle=0)

Given angles or other cyclic coordinates with the specified period, such that the angle + period has the same meaning as the original angle, this function adds or subtracts multiples of the period such that the result has the same meaning, but now all angles lie in an interval of length the specified period, centered on the reference angle ref_angle.

source
Pixell.sky2pixFunction
sky2pix(m::Enmap, skycoords)

Convert sky coordinates to 1-indexed pixels. While typically WCS can specify custom units, Pixell uses radians throughout.

Arguments:

  • m::Enmap: the map to obtain the coordinates from
  • skycoords: skycoords should be a 2-d array where "skycoords[:, i]" is the i-th set of coordinates, or a 1-d array representing a single set of coordinates.

Returns:

  • Array: same shape as skycoords

Examples

julia> shape, wcs = fullsky_geometry(deg2rad(1))
julia> sky2pix(shape, wcs, [30.0, 50.0])
2-element Vector{Float64}:
 151.0
 141.0
source
Pixell.sky2pix!Function
sky2pix!(m::Enmap, skycoords, pixcoords)

Convert sky coordinates to 1-indexed pixels, in-place. While typically WCS can specify custom units, Pixell uses radians throughout.

Arguments:

  • m::Enmap: the map that provides a coordinate system
  • skycoords: sky coordinates should be a 2-d array where "skycoords[:, i]" is the i-th set of coordinates, or a 1-d array representing a single set of coordinates.
  • pixcoords: output array for pixel coordinates, must be same same as pixcoords

Returns:

  • Array: pixcoords

Examples

julia> shape, wcs = fullsky_geometry(deg2rad(1))
       skycoords =  similar(pixcoords)
       pixcoords =  100 .* rand(2,4096 * 2)
julia> sky2pix!(shape, wcs, skycoords, pixcoords)
source
Pixell.sky2pixMethod
sky2pix(m::Enmap{T,N,AA,<:AbstractCARWCS}, ra, dec)

Compute 1-indexed pixels into sky coordinates.

Only implemented for CAR (Clenshaw-Curtis variant) projections. Takes RA and DEC and returns a tuple containing the corresponding pixel indices. If vectors of RA and DEC are given, then vectors of pixel indices will be returned.

Examples

julia> shape, wcs = fullsky_geometry(deg2rad(1))
julia> sky2pix(shape, wcs, deg2rad(30.0), deg2rad(80.0))
(151.0, 171.0)
source