Index
Pixell.AbstractDT
Pixell.CarClenshawCurtis
Pixell.Enmap
Pixell.FFTLogPlan
Pixell.Gnomonic
Healpix.alm2map
Healpix.map2alm
Pixell.create_sht_band
Pixell.distance_transform
Pixell.distance_transform
Pixell.dplanck
Pixell.first_last_rings_in_fullsky
Pixell.fullringnum
Pixell.fullringsize
Pixell.fullsky_geometry
Pixell.getlmax
Pixell.iscyl
Pixell.laxes_cyl
Pixell.make_cc_geom_info
Pixell.metric
Pixell.pix2sky
Pixell.pix2sky
Pixell.pix2sky!
Pixell.pixareamap
Pixell.pixareamap!
Pixell.rewind
Pixell.sky2pix
Pixell.sky2pix
Pixell.sky2pix!
Pixell.skyarea
Pixell.AbstractDT
— TypeDistance Transform (DT), especially Spherical Distance Transforms (SDT)
Pixell.CarClenshawCurtis
— TypeFast custom WCS structure.
Pixell.Enmap
— TypeMap 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.
Pixell.FFTLogPlan
— TypeFFTLog from Hamilton 2000. Construct this type with plan_fftlog
.
Pixell.Gnomonic
— TypeFast custom WCS structure.
Healpix.alm2map
— Methodperform inverse SHT of an intensity map
Healpix.map2alm
— Methodperform forward SHT of an intensity map
Pixell.create_sht_band
— MethodComplete rings for spherical harmonic transforms. Also converts to AbstractArray{Float64,2}.
Pixell.distance_transform
— MethodBrute-force spherical distance transform. ~ O(Nₚᵢₓ × N₀). Doesn't handle wrapping.
Pixell.distance_transform
— MethodCompute an exact distance transform ~ O(Nₚᵢₓ), based on Mullikin 1992 and Danielsson 1980
Pixell.dplanck
— FunctionThe 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}\]
Pixell.first_last_rings_in_fullsky
— MethodGet index of ring given angle
Pixell.fullringnum
— MethodNumber of rings in a fullsky version of this WCS.
Pixell.fullringsize
— MethodNumber of pixels in full CAR ring of this WCS.
Pixell.fullsky_geometry
— Methodfullsky_geometry([W=CarClenshawCurtis], res; shape = nothing, dims = ())
Generates a full-sky geometry.
Arguments:
proj=CarClenshawCurtis()
: [optional] projectionres
: 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]))
Pixell.getlmax
— MethodGenerate an estimate of the Nyquist frequency for a map
Pixell.iscyl
— MethodCheck if a WCS is a cylindrical pixelization
Pixell.laxes_cyl
— Methodlaxes_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.
Pixell.make_cc_geom_info
— MethodClenshaw-Curtis quadrature weights, make a Libsharp GeomInfo
Pixell.metric
— MethodMetric on the sphere for RA (α) and DEC (δ)
Pixell.pix2sky
— Functionpix2sky(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 systempixcoords
: 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
Pixell.pix2sky!
— Methodpix2sky!(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 systempixcoords
: 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)
Pixell.pix2sky
— Methodpix2sky(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)
Pixell.pixareamap!
— MethodIn-place write to pixels the areas of those pixels in steradians.
Pixell.pixareamap
— MethodGenerate a similar Enmap whose pixel values are the areas of the pixels in steradians.
Pixell.rewind
— Methodrewind(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
.
Pixell.sky2pix
— Functionsky2pix(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 fromskycoords
: 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
Pixell.sky2pix!
— Functionsky2pix!(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 systemskycoords
: 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)
Pixell.sky2pix
— Methodsky2pix(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)
Pixell.skyarea
— MethodArea of a patch given shape and WCS, in steradians.