healpy.pixelfunc.ang2pix

Contents

healpy.pixelfunc.ang2pix#

healpy.pixelfunc.ang2pix(nside, theta, phi, nest=False, lonlat=False, latauto=False, latbounce=True)#

ang2pix : nside,theta[rad],phi[rad],nest=False,lonlat=False,latauto=False,latbounce=True -> ipix (default:RING)

Parameters:
nsideint, scalar or array-like

The healpix nside parameter, must be a power of 2, less than 2**30

theta, phifloat, scalars or array-like

Angular coordinates of a point on the sphere

nestbool, optional

if True, assume NESTED pixel ordering, otherwise, RING pixel ordering

lonlatbool

If True, input angles are assumed to be longitude and latitude in degree, otherwise, they are co-latitude and longitude in radians.

latauto: bool, optional

If True, automatically adjust latitudes to be within [-90, 90] range.

latbounce: bool, optional

If True, longitude is not affected, as if scanning and bouncing back at the pole. If False, longitude is adjusted by 180 degrees, as if scanning through the pole. Needs latauto=True.

Returns:
pixint or array of int

The healpix pixel numbers. Scalar if all input are scalar, array otherwise. Usual numpy broadcasting rules apply.

See also

pix2ang, pix2vec, vec2pix

Examples

Note that some of the test inputs below that are on pixel boundaries such as theta=pi/2, phi=pi/2, have a tiny value of 1e-15 added to them to make them reproducible on i386 machines using x87 floating point instruction set (see healpy/healpy#528).

>>> import healpy as hp
>>> hp.ang2pix(16, np.pi/2, 0)
1440
>>> print(hp.ang2pix(16, [np.pi/2, np.pi/4, np.pi/2, 0, np.pi], [0., np.pi/4, np.pi/2 + 1e-15, 0, 0]))
[1440  427 1520    0 3068]
>>> print(hp.ang2pix(16, np.pi/2, [0, np.pi/2 + 1e-15]))
[1440 1520]
>>> print(hp.ang2pix([1, 2, 4, 8, 16], np.pi/2, 0))
[   4   12   72  336 1440]
>>> print(hp.ang2pix([1, 2, 4, 8, 16], 0, 0, lonlat=True))
[   4   12   72  336 1440]