Welcome to healpy documentation!¶
Tutorial¶
Healpy tutorial¶
Creating and manipulating maps¶
Maps are simply numpy arrays, each array element refers to a location in the sky as defined by the Healpix pixelization schemes, see the healpix website.
The resolution of the map is defined by the NSIDE parameter, the nside2npix() function gives the number of pixel NPIX of the map:
>>> import numpy as np
>>> import healpy as hp
>>> NSIDE = 32
>>> m = np.arange(hp.nside2npix(NSIDE))
>>> hp.mollview(m, title="Mollview image RING")

Healpix supports two different ordering schemes, RING or NESTED, by default healpy maps are in *RING* ordering.
In order to work with NESTED ordering, all map related functions support the nest keyword, for example:
>>> hp.mollview(m, nest=True, title="Mollview image NESTED")

Reading and writing maps to file¶
Maps are read with the read_map() function:
>>> wmap_map_I = hp.read_map('../healpy/test/data/wmap_band_imap_r9_7yr_W_v4.fits')
by default input maps are converted to *RING* ordering, if they are in NESTED ordering. You can otherwise specify nest=True to retrieve a map is NESTED ordering, or nest=None to keep the ordering unchanged.
By default read_map() loads the first column, for reading other columns you can specify the field keyword.
write_map() writes a map to disk in FITS format, if the input map is a list of 3 maps, they are written to a single file as I,Q,U polarization components:
>>> hp.write_map("my_map.fits", wmap_map_I)
Visualization¶
Mollweide projection with mollview() is the most common visualization tool for HEALPIX maps, it supports also coordinate transformation:
>>> hp.mollview(wmap_map_I, coord=['G','E'], title='Histogram equalized Ecliptic', unit='mK', norm='hist', min=-1,max=1, xsize=2000)
>>> hp.graticule()
coord does galactic to ecliptic coordinate transformation, norm=’hist’ sets a histogram equalized color scale and xsize increases the size of the image. graticule() adds meridians and parallels.

gnomview() instead provides gnomonic projection around a position specified by rot:
>>> hp.gnomview(wmap_map_I, rot=[0,0.3], title='GnomView', unit='mK', format='%.2g')
shows a projection of the galactic center, xsize and ysize change the dimension of the sky patch.
mollzoom() is a powerful tool for interactive inspection of a map, it provides a mollweide projection where you can click to set the center of the adjacent gnomview panel.
Masked map, partial maps¶
By convention HEALPIX uses -1.6375e+30 to mark invalid or unseen pixels, this is stored in healpy as the constant UNSEEN().
All healpy functions automatically deal with maps with UNSEEN pixels, for example mollview() marks in grey that sections of a map.
There is an alternative way of dealing with UNSEEN pixel based on the numpy MaskedArray class, ma() loads a map as a masked array:
>>> mask = hp.read_map('../healpy/test/data/wmap_temperature_analysis_mask_r9_7yr_v4.fits').astype(np.bool)
>>> wmap_map_I_masked = hp.ma(wmap_map_I)
>>> wmap_map_I_masked.mask = np.logical_not(mask)
by convention the mask is 0 where the data are masked, while numpy defines data masked when the mask is True, so it is necessary to flip the mask.
>>> hp.mollview(wmap_map_I_masked.filled())
filling a masked array fills the UNSEEN value in and return a standard array that can be used by mollview. compressed() instead removes all the masked pixels and returns a standard array that can be used for examples by the matplotlib hist() function:
>>> import matplotlib.pyplot as plt
>>> plt.hist(wmap_map_I_masked.compressed(), bins = 1000)
Spherical harmonic transforms¶
healpy provides bindings to the C++ HEALPIX library for performing spherical harmonic transforms. anafast() computes the angular power spectrum of a map:
>>> LMAX = 1024
>>> cl = hp.anafast(wmap_map_I_masked.filled(), lmax=LMAX)
the relative ell array is just:
>>> ell = np.arange(len(cl))
therefore we can plot a normalized CMB spectrum and write it to disk:
>>> plt.figure()
>>> plt.plot(ell, ell * (ell+1) * cl)
>>> plt.xlabel('ell'); plt.ylabel('ell(ell+1)cl'); plt.grid()
>>> hp.write_cl('cl.fits', cl)

Gaussian beam map smoothing is provided by smoothing():
>>> wmap_map_I_smoothed = hp.smoothing(wmap_map_I, fwhm=60, arcmin=True)
>>> hp.mollview(wmap_map_I_smoothed, min=-1, max=1, title='Map smoothed 1 deg')
Installation¶
Installation procedure for Healpy¶
Requirements¶
Healpy depends on the Healpix C++ and cfitsio C libraries. Source code is include with Healpy and you do not have to install them separately.
Building against external Healpix and cfitsio¶
Healpy uses pkg-config to detect the presence of the Healpix and cfitsio libraries. pkg-config is available on most systems. If you do not have pkg-config installed, then Healpy will download and use (but not install) a Python clone called pykg-config.
If you want to provide your own external builds of Healpix and cfitsio, then download the following packages:
If you are going to install the packages in a nonstandard location (say, –prefix=/path/to/local), then you should set the environment variable PKG_CONFIG_PATH=/path/to/local/lib/pkgconfig when building. No other environment variable settings are necessary, and you do not need to set PKG_CONFIG_PATH to use Healpy after you have built it.
Then, unpack each of the above packages and build them with the usual ‘configure; make; make install’ recipe.
Installation¶
healpy is available on pipy, you can install it with:
pip install healpy
otherwise, you can download a source tarball from:
https://pypi.python.org/pypi/healpy
DO NOT DOWNLOAD from github, github does not include the dependencies.
and build it with:
cd healpy
python setup.py build
sudo python setup.py install
If everything goes fine, you can test it:
cd build/lib*
python
>>> import matplotlib.pyplot as plt
>>> import numpy as np
>>> import healpy as H
>>> H.mollview(np.arange(12))
>>> plt.show()
or run the test suite with nose:
nosetests -v
If the plot looks good, you can install:
sudo python setup.py install # install in default location, need root rights
or:
python setup.py install --install-lib=~/Softs/Python # will install healpy in directory ~/Softs/Python, which then must be in your PYTHONPATH
or:
python setup.py install --user # will install it in your User python directory (python >= 2.6)
Known issues¶
- Incompatibility with cfitisio from HEASOFT: due to a conflict of header file names it is currently not possible to use the cfitsio library provided with the HEASOFT package for compilation of Healpix C++. HEASOFT’s include directory contains a file called “rotmatrix.h” which clashes with Healpix’s own rotmatrix.h.
- Compilation problems in the C++ package: some gcc versions (we have reports for 4.4.5 and 4.4.6) crash with an internal compiler error during compilation of libsharp. Unfortunately we have not found a workaround for this compiler problem. To our knowledge, it has been fixed in gcc 4.4.7 and in the 4.5.x and newer versions.
Development install¶
- Developers building from a snapshot of the github repository need:
- autoconf (in Ubuntu: sudo apt-get install autoconf automake libtool pkg-config)
- cython > 0.14
- run git submodule init and git submodule update to get the healpix sources
the best way to install healpy if you plan to develop is to build the C++ extensions in place with:
python setup.py build_ext --inplace
the add the healpy/healpy folder to your PYTHONPATH
Clean¶
When you run “python setup.py”, temporary build products are placed in the “build” directory. If you want to clean out and remove the “build” directory, then run:
python setup.py clean --all
Reference¶
sphtfunc – Spherical harmonic transforms¶
From map to spherical harmonics¶
anafast(map1[, map2, nspec, lmax, mmax, ...]) | Computes the power spectrum of an Healpix map, or the cross-spectrum between two maps if map2 is given. |
map2alm(maps[, lmax, mmax, iter, pol, ...]) | Computes the alm of an Healpix map. |
From spherical harmonics to map¶
synfast(cls, nside[, lmax, mmax, alm, pol, ...]) | Create a map(s) from cl(s). |
alm2map(alms, nside[, lmax, mmax, pixwin, ...]) | Computes an Healpix map given the alm. |
alm2map_der1(alm, nside[, lmax, mmax]) | Computes an Healpix map and its first derivatives given the alm. |
Spherical harmonic transform tools¶
smoothing(*args, **kwds) | Smooth a map with a Gaussian symmetric beam. |
smoothalm(alms[, fwhm, sigma, invert, pol, ...]) | Smooth alm with a Gaussian symmetric beam function. |
alm2cl(alms1[, alms2, lmax, mmax, lmax_out, ...]) | Computes (cross-)spectra from alm(s). |
synalm(cls[, lmax, mmax, new]) | Generate a set of alm given cl. |
almxfl(alm, fl[, mmax, inplace]) | Multiply alm by a function of l. |
pixwin(nside[, pol]) | Return the pixel window function for the given nside. |
Alm() | This class provides some static methods for alm index computation. |
Other tools¶
gauss_beam(fwhm[, lmax, pol]) | Gaussian beam window function |
visufunc – Visualisation¶
Map projections¶
mollview([map, fig, rot, coord, unit, ...]) | Plot an healpix map (given as an array) in Mollweide projection. |
gnomview([map, fig, rot, coord, unit, ...]) | Plot an healpix map (given as an array) in Gnomonic projection. |
cartview([map, fig, rot, zat, coord, unit, ...]) | Plot an healpix map (given as an array) in Cartesian projection. |
Graticules¶
graticule([dpar, dmer, coord, local]) | Draw a graticule on the current Axes. |
delgraticules() | Delete all graticules previously created on the Axes. |
Tracing lines or points¶
projplot(*args, **kwds) | projplot is a wrapper around matplotlib.Axes.plot() to take into account the |
projscatter(*args, **kwds) | Projscatter is a wrapper around matplotlib.Axes.scatter() to take into account the spherical projection. |
projtext(*args, **kwds) | Projtext is a wrapper around matplotlib.Axes.text() to take into account the spherical projection. |
rotator – Rotation and geometry functions¶
Rotation¶
Rotator([rot, coord, inv, deg, eulertype]) | Rotation operator, including astronomical coordinate systems. |
rotateVector(rotmat, vec[, vy, vz, do_rot]) | Rotate a vector (or a list of vectors) using the rotation matrix given as first argument. |
rotateDirection(rotmat, theta[, phi, ...]) | Rotate the vector described by angles theta,phi using the rotation matrix given as first argument. |
projector – Spherical projections¶
Basic classes¶
SphericalProj([rot, coord, flipconv]) | This class defines functions for spherical projection. |
GnomonicProj([rot, coord, xsize, ysize, reso]) | This class provides class methods for Gnomonic projection. |
MollweideProj([rot, coord, xsize]) | This class provides class methods for Mollweide projection. |
CartesianProj([rot, coord, xsize, ysize, ...]) | This class provides class methods for Cartesian projection. |
projaxes – Axes for projection¶
Basic classes¶
SphericalProjAxes(ProjClass, *args, **kwds) | Define a special Axes to take care of spherical projection. |
GnomonicAxes(*args, **kwds) | Define a gnomonic Axes to handle gnomonic projection. |
HpxGnomonicAxes(*args, **kwds) | |
MollweideAxes(*args, **kwds) | Define a mollweide Axes to handle mollweide projection. |
HpxMollweideAxes(*args, **kwds) | |
CartesianAxes(*args, **kwds) | Define a cylindrical Axes to handle cylindrical projection. |
HpxCartesianAxes(*args, **kwds) |