Welcome to healpy documentation!

Tutorial

Healpy tutorial

Creating and manipulating maps

Maps are simply numpy arrays, where each array element refers to a location in the sky as defined by the Healpix pixelization schemes (see the healpix website).

Note: Running the code below in a regular Python session will not display the maps; it’s recommended to use IPython:

% ipython

...then select the appropriate backend display:

>>> %matplotlib inline # for IPython notebook
>>> %matplotlib qt     # using Qt (e.g. Windows)
>>> %matplotlib osx    # on Macs
>>> %matplotlib gtk    # GTK

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")
_images/moll_nside32_ring.png

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")
_images/moll_nside32_nest.png

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 also supports 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.

_images/wmap_histeq_ecl.png

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 in the UNSEEN value 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)
_images/wmap_powspec.png

Gaussian beam map smoothing is provided by smoothing():

>>> wmap_map_I_smoothed = hp.smoothing(wmap_map_I, fwhm=np.radians(1.))
>>> 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 for both is include with Healpy and is built automatically, so you do not need to install them yourself.

Quick installation with Pip

The quickest way to install Healpy is with pip (>= 1.4.2), which automatically fetches the latest version of Healpy and any missing dependencies:

pip install --user healpy

If you have installed with pip, you can keep your installation up to date by upgrading from time to time:

pip install --user --upgrade healpy

Even quicker installation on Mac OS with MacPorts

If you are using a Mac and have the MacPorts package manager, it’s even easer to install Healpy with:

sudo port install py27-healpy

Binary apt-get style packages are also available in the development versions of Debian (sid) and Ubuntu (utopic).

Almost-as-quick installation from official source release

Healpy is also available in the Python Package Index (PyPI). You can download it with:

curl -O https://pypi.python.org/packages/source/h/healpy/healpy-1.7.4.tar.gz

and build it with:

tar -xzf healpy-1.7.4.tar.gz
pushd healpy-1.7.4
python setup.py install --user
popd

If everything goes fine, you can test it:

python
>>> import matplotlib.pyplot as plt
>>> import numpy as np
>>> import healpy as hp
>>> hp.mollview(np.arange(12))
>>> plt.show()

or run the test suite with nose:

cd healpy-1.7.4 && python setup.py test

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.

Development install

Developers building from a snapshot of the github repository need:

  • autoconf and libtool (in Debian or Ubuntu: sudo apt-get install autoconf automake libtool pkg-config)
  • cython > 0.16
  • run git submodule init and git submodule update to get the bundled 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

then 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(map_in, *args, **kwds) Smooth a map with a Gaussian symmetric beam.
smoothalm(alms[, fwhm, sigma, pol, mmax, ...]) 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, verbose]) 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.
orthview([map, fig, rot, coord, unit, ...]) Plot an healpix map (given as an array) in Orthographic 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.

Geometrical helpers

vec2dir(vec[, vy, vz, lonlat]) Transform a vector to angle given by theta,phi.
dir2vec(theta[, phi, lonlat]) Transform a direction theta,phi to a unit vector.
angdist(dir1, dir2[, lonlat]) Returns the angular distance between dir1 and dir2.

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)

Methods

MollweideAxes(*args, **kwds) Define a mollweide Axes to handle mollweide projection.
HpxMollweideAxes(*args, **kwds)

Methods

CartesianAxes(*args, **kwds) Define a cylindrical Axes to handle cylindrical projection.
HpxCartesianAxes(*args, **kwds)

Methods

zoomtool – Interactive visualisation

Interactive map visualization

mollzoom

Indices and tables