Healpy Tools¶
The healpytools provides an extension for the healpy framework (https://healpy.readthedocs.io). There are various functionalities on top of healpy, e.g. sample random directions in pixel or create distributions on the sphere (dipole, fisher).
Extensions to healpy that covers: -> pixel, vector, angular transformations -> drawing vectors uniformly in pixel -> various probability density functions (exposure, fisher, dipole)
-
healpytools.
ang2pix
(nside, phi, theta, nest=False)[source]¶ Convert spherical angle (astrotools definition) to HEALpixel ipix Substitutes hp.ang2pix
- Parameters
nside – nside of the healpy pixelization
phi – longitude in astrotools definition
theta – latitude in astrotools definition
nest – set True in case you work with healpy’s nested scheme
- Returns
pixel number(s)
-
healpytools.
ang2vec
(phi, theta)[source]¶ Substitutes healpy.ang2vec() to use our angle convention
- Parameters
phi – longitude, range (pi, -pi), 0 points in x-direction, pi/2 in y-direction
theta – latitude, range (pi/2, -pi/2), pi/2 points in z-direction
- Returns
vector of shape (3, n)
-
healpytools.
angle
(nside, ipix, jpix, nest=False, each2each=False)[source]¶ Give the angular distance between two pixel arrays.
- Parameters
nside – nside of the healpy pixelization
ipix – healpy pixel i (either int or array like int)
jpix – healpy pixel j (either int or array like int)
nest – use the nesting scheme of healpy
each2each – if true, calculates every combination of the two lists v1, v2
- Returns
angular distance in radians
-
healpytools.
check_problematic_pixel
(nside, ipix, a0, zmax, deviation=0.5, coord_system='gal')[source]¶ Checks input pixel for exposure deviation within the corner points from more than certain threshold (default: 0.5).
- Parameters
nside – nside of the healpy pixelization
ipix – pixel number(s)
a0 – latitude of detector (-90, 90) in degrees (default: Auger)
zmax – maximum acceptance zenith angle (0, 90) degrees
deviation – maximum deviation between exposure values in pixel corners
coord_system – choose between different coordinate systems - gal, eq, sgal, ecl
-
healpytools.
dipole_pdf
(nside, a, v, y=None, z=None, pdf=True)[source]¶ Probability density function of a dipole. Returns 1 + a * cos(theta) for all pixels in hp.nside2npix(nside).
- Parameters
nside – nside of the healpy map
a – amplitude of the dipole, 0 <= a <= 1, automatically clipped
v – either (x, y, z) vector of the pixel center(s) or only x-coordinate
y – y-coordinate(s) of the center
z – z-coordinate(s) of the center
pdf – if false, return unnormalized map (modulation around 1)
- Returns
weights
-
healpytools.
exposure_pdf
(nside=64, a0=- 35.25, zmax=60, coord_system='gal', pdf=True)[source]¶ Exposure (type: density function) of an experiment located at equatorial declination a0 and measuring events with zenith angles up to zmax degrees.
- Parameters
nside – nside of the output healpy map
a0 – equatorial declination [deg] of the experiment (default: AUGER, a0=-35.25 deg)
zmax – maximum zenith angle [deg] for the events
coord_system – choose between different coordinate systems - gal, eq, sgal, ecl
pdf – if false, return relative exposure
- Returns
weights of the exposure map
-
healpytools.
fisher_pdf
(nside, v, y=None, z=None, k=None, threshold=4, sparse=False, pdf=True)[source]¶ Probability density function of a fisher distribution of healpy pixels with mean direction (x, y, z) and concentration parameter kappa; normalized to 1.
- Parameters
nside – nside of the healpy map
v – either (x, y, z) vector of the center or only x-coordinate
y – y-coordinate of the center
z – z-coordinate of the center
k – kappa for the fisher distribution, 1 / sigma**2
threshold – Threshold in sigma up to where the distribution should be calculated
sparse – returns the map in the form (pixels, weights); this may be meaningfull for small distributions
pdf – if false, return values with maximum 1
- Returns
pixels, weights at pixels if sparse else a full map with length nside2npix(nside)
-
healpytools.
norder2npix
(norder)[source]¶ Give the number of pixel for the given HEALpix order.
- Parameters
norder – norder of the healpy pixelization
- Returns
npix, number of pixels of the healpy pixelization
-
healpytools.
norder2nside
(norder)[source]¶ Give the HEALpix nside parameter for the given HEALpix order.
- Parameters
norder – norder of the healpy pixelization
- Returns
nside, nside of the healpy pixelization
-
healpytools.
npix2norder
(npix)[source]¶ Give the HEALpix order for the given number of pixel.
- Parameters
npix – number of pixels of the healpy pixelization
- Returns
norder, norder of the healpy pixelization
-
healpytools.
nside2norder
(nside)[source]¶ Give the HEALpix order for the given HEALpix nside parameter.
- Parameters
nside – nside of the healpy pixelization
- Returns
norder, norder of the healpy pixelization
-
healpytools.
pix2ang
(nside, ipix, nest=False)[source]¶ Convert HEALpixel ipix to spherical angles (astrotools definition) Substitutes hp.pix2ang
- Parameters
nside – nside of the healpy pixelization
ipix – pixel number(s)
nest – set True in case you work with healpy’s nested scheme
- Returns
angles (phi, theta) in astrotools definition
-
healpytools.
pix2map
(nside, ipix)[source]¶ Converts healpy pixel to healpix map
- Parameters
nside – nside of the healpy pixelization
ipix – pixel number(s) in nside scheme
- Returns
healpix map of size npix
-
healpytools.
pix2vec
(nside, ipix, nest=False)[source]¶ Convert HEALpixel ipix to cartesian vector Substitutes hp.pix2vec
- Parameters
nside – nside of the healpy pixelization
ipix – pixel number(s)
nest – set True in case you work with healpy’s nested scheme
- Returns
vector of the pixel center(s)
-
healpytools.
query_disc
(nside, v, radius, nest=False, **kwargs)[source]¶ Substitutes hp.query_disc but supports also pixel number input
- Parameters
nside – nside of the healpy pixelization
v – either (x, y, z) vector of the pixel center or healpy pixel
radius – radius of disk in radians
nest – set True in case you work with healpy’s nested scheme
- Returns
pixels thay lie within radius from center v
-
healpytools.
rand_exposure_vec_in_pix
(nside, ipix, a0=- 35.25, zmax=60, coord_system='gal', deviation=0.5, nest=False)[source]¶ Draw vectors from a distribution within a HEALpixel that follow the exposure distribution within the pixel. It is much slower than rand_vec_in_pix() and should therefore only be used for problematic pixels (close to zero exposure).
- Parameters
nside – nside of the healpy pixelization
ipix – pixel number(s)
a0 – latitude of detector (-90, 90) in degrees (default: Auger)
zmax – maximum acceptance zenith angle (0, 90) degrees
coord_system – choose between different coordinate systems - gal, eq, sgal, ecl
deviation – maximum relative deviation between exposure values in pixel corners
nest – set True in case you work with healpy’s nested scheme
- Returns
vectors containing events from the pixel(s) specified in ipix
-
healpytools.
rand_pix_from_map
(healpy_map, n=1)[source]¶ Draw n random pixels from a HEALpix map.
- Parameters
healpy_map – healpix map (not necessarily normalized)
n – number of pixels that are drawn from the map
- Returns
an array of pixels with size n, that are drawn from the map
-
healpytools.
rand_vec_from_map
(healpy_map, n=1, nest=False)[source]¶ Draw n random vectors from a HEALpix map.
- Parameters
healpy_map – healpix map (not necessarily normalized)
n – number of pixels that are drawn from the map
nest – set True in case you work with healpy’s nested scheme
- Returns
an array of vectors with size n, that are drawn from the map
-
healpytools.
rand_vec_in_pix
(nside, ipix, nest=False)[source]¶ Draw vectors from a uniform distribution within a HEALpixel.
- Parameters
nside – nside of the healpy pixelization
ipix – pixel number(s)
nest – set True in case you work with healpy’s nested scheme
- Returns
vectors containing events from the pixel(s) specified in ipix
-
healpytools.
rotate_map
(healpy_map, rotation_axis, rotation_angle)[source]¶ Perform rotation of healpy map for given rotation axis and angle(s).
- Parameters
healpy_map – healpix map to be rotated
rotation_axis – rotation axis, either np.array([x, y, z]) or ndarray with shape (3, n)
rotation_angle – rotation angle in radians, either float or array size n
- Returns
rotated healpy map, same shape as input healpy map or shape (n, npix)
-
healpytools.
skymap_mean_quantile
(skymap, quantile=0.68)[source]¶ Calculates mean direction and quantile from an arbitrary healpy skymap
- Parameters
skymap – healpy skymap of valid healpy size
quantile – quantile for which the scatter angle is calculated
- Returns
either count, frequency, mean or rms maps
-
healpytools.
statistic
(nside, v, y=None, z=None, statistics='count', vals=None)[source]¶ Create HEALpix map of count, frequency or mean or rms value.
- Parameters
nside – nside of the healpy pixelization
v – either (x, y, z) vector of the pixel center(s) or only x-coordinate
y – y-coordinate(s) of the center
z – z-coordinate(s) of the center
statistics – keywords ‘count’, ‘frequency’, ‘mean’ or ‘rms’ possible
vals – values (array like) for which the mean or rms is calculated
- Returns
either count, frequency, mean or rms maps
-
healpytools.
vec2ang
(v)[source]¶ Substitutes healpy.vec2ang() to use our angle convention.
- Parameters
v – vector(s) of shape (3, n)
- Returns
tuple consisting of
phi [range (pi, -pi), 0 points in x-direction, pi/2 in y-direction],
theta [range (pi/2, -pi/2), pi/2 points in z-direction]
-
healpytools.
vec2pix
(nside, v, y=None, z=None, nest=False)[source]¶ Convert HEALpixel ipix to spherical angles (astrotools definition) Substitutes hp.vec2pix
- Parameters
nside – nside of the healpy pixelization
v – either (x, y, z) vector of the pixel center(s) or only x-coordinate
y – y-coordinate(s) of the center
z – z-coordinate(s) of the center
nest – set True in case you work with healpy’s nested scheme
- Returns
vector of the pixel center(s)