Classes
HealpixBase
- class mhealpy.HealpixBase(uniq=None, order=None, nside=None, base=None, scheme='ring')[source]
Bases:
object
Basic operations related to HEALPix pixelization, for which the map contents information is not needed. This class is conceptually very similar the the Healpix_Base class of Healpix_cxx.
Single resolution maps are fully defined by specifying their order (or NSIDE) and ordering scheme (“RING” or “NESTED”).
Multi-resolution maps follow an explicit “NUNIQ” scheme, with each pixel identfied by a _uniq_ number. No specific is needed nor guaranteed.
Warning
The initialization input is not validated by default. Consider calling is_mesh_valid() after initialization, otherwise results might be unexpected.
- Parameters
uniq (array) – Explicit numbering of each pixel in an “NUNIQ” scheme.
order (int) – Order of HEALPix map.
nside (int) – Alternatively, you can specify the NSIDE parameter.
scheme (str) – Healpix scheme. Either ‘RING’, ‘NESTED’ or ‘NUNIQ’
base (HealpixBase) – Alternatively, you can copy the properties of another HealpixBase object
- classmethod adaptive_moc_mesh(max_nside, split_fun)[source]
Return a MOC mesh with an adaptive resolution determined by an arbitrary function.
- Parameters
max_nside (int) – Maximum HEALPix nside to consider
split_fun (function) – This method should return
True
if a pixelorder (should be split into pixel of a higher) –
otherwise. (and False) –
integers (It takes two) –
start (inclusive) and
stop
(exclusive) –a (which correspond to a single pixel in nested rangeset format for) –
max_nside. (map of nside) –
- Returns
HealpixBase
- classmethod moc_from_pixels(nside, pixels, nest=False)[source]
Return a MOC mesh where a list of pixels are kept at a given nside, and every other pixel is appropiately downsampled.
Also see the more generic
adaptive_moc()
andadaptive_moc_mesh()
.- Parameters
nside (int) – Maximum healpix NSIDE (that is, the NSIDE for the pixel list)
pixels (array) – Pixels that must be kept at the finest pixelation
nest (bool) – Whether the pixels are a ‘NESTED’ or ‘RING’ scheme
- conformable(other)[source]
For single-resolution maps, return
True
if both maps have the same nside and scheme.For MOC maps, return True if both maps have the same list of UNIQ pixels (including the ordering)
- property npix
Get number of pixels.
For multi-resolutions maps, this corresponds to the number of utilized UNIQ pixels.
- Returns
int
- property order
Get map order
- Returns
int
- property nside
Get map NSIDE
- Returns
int
- property scheme
Return HEALPix scheme
- Returns
Either ‘NESTED’, ‘RING’ or ‘NUNIQ’
- Return type
str
- property is_nested
Return true if scheme is NESTED or NUNIQ
- Return
bool
- property is_ring
Return true if scheme is RING
- Return
bool
- property is_moc
Return true if this is a Multi-Dimensional Coverage (MOC) map (multi-resolution)
- Returns
bool
- pix_rangesets(nside=None, argsort=False)[source]
Get the equivalent range of child pixels in nested scheme for a map of equal or higher nside
- Parameters
nside (int or None) – Nside of output range sets. If None, the map nside will be used. nside = mhealpy.MAX_NSIDE returns the cached result
argsort (bool) – Also return also the indices that would sort the array.
- Returns
- With columns named ‘start’ (inclusive) and
’stop’ (exclusive)
- Return type
recarray
- pix_order_list()[source]
Get a list of lists containing all pixels sorted by order
- Returns
- (
pix_per_order
,nest_pix_per_order
) Each list has a size equal to the map order. Each element is a list of all pixels whose order matches the index of the list position. The first output contains the index of the pixels, while the second contains their coresponding pixel number in a nested scheme.
- (
- Return type
(list, list)
- pix2range(nside, pix)[source]
Get the equivalent range of child pixels in nested scheme for a map of equal or higher nside
- Parameters
nside (int) – Nside of output range sets
pix (int or array) – Pixel numbers
- Returns
- Start pixel (inclusive) and
stop pixel (exclusive)
- Return type
(int or array, int or array)
- pixarea(pix=0)[source]
Return area of pixel in steradians
- Parameters
pix (int or array) – Pixel number. Only relevant for MOC maps
- Returns
float or array
- pix2ang(pix, lonlat=False)[source]
Return the coordinates of the center of a pixel
- Parameters
pix (int or array) –
- Returns
(float or array, float or array)
- pix2vec(pix)[source]
Return a vector corresponding to the center of a pixel
- Parameters
pix (int or array) –
- Returns
Size (3,N)
- Return type
array
- ang2pix(theta, phi, lonlat=False)[source]
Get the pixel (as used in []) that contains a given coordinate
- Parameters
theta (float or array) – Zenith angle
phi (float or arrray) – Azimuth angle
- Returns
int or array
- vec2pix(x, y, z)[source]
Get the pixel (as used in []) that contains a given coordinate
- Parameters
theta (float or array) – Zenith angle
phi (float or arrray) – Azimuth angle
- Returns
int or array
- pix2uniq(pix)[source]
Get the UNIQ representation of a given pixel index.
- Parameters
pix (int) – Pixel number in the current scheme (as used for [])
- property uniq
Get an array with the NUNIQ numbers for all pixels
- nest2pix(pix)[source]
Get the corresponding pixel in the current grid for a pixel in NESTED scheme. For MOC map, return the pixel that contains it.
- Parameters
pix (int or array) – Pixel number in NESTED scheme. Must correspond to a map of the same order as the current.
- Returns
int or array
- get_interp_weights(theta, phi=None, lonlat=False)[source]
Return the 4 closest pixels on the two rings above and below the location and corresponding weights. Weights are provided for bilinear interpolation along latitude and longitude
- Parameters
theta (float or array) – Zenith angle (rad)
phi (float or array) – Azimuth angle (rad)
- Returns
- (pixels, weights), each with of (4,) if the input is scalar,
if (4,N) where N is size of theta and phi. For MOC maps, these pixel numbers might repeate.
- Return type
tuple
- get_all_neighbours(theta, phi=None, lonlat=False)[source]
Return the 8 nearest pixels. For MOC maps, these might repeat, as this is equivalent to raterizing the maps to the highest order, getting the neighbohrs, and then finding the pixels tha contain them.
- Parameters
theta (float or int or array) – Zenith angle (rad). If phi is
None
, these are assummed to be pixels numbers. For MOC maps, these are assumed to be pixel numbers in NESTED scheme for the equivalent single-resolution map of the highest order.phi (float or array or None) – Azimuth angle (rad)
- Returns
- pixel number of the SW, W, NW, N, NE, E, SE and S neighbours,
shape is (8,) if input is scalar, otherwise shape is (8, N) if input is of length N. If a neighbor does not exist (it can be the case for W, N, E and S) the corresponding pixel number will be -1.
- Return type
array
- is_mesh_valid()[source]
Return
True
if the map pixelization is valid. For single resolution this simply checks that the size is a valid NSIDE value. For MOC maps, it checks that every point in the sphere is covered by one and only one pixel.- Returns
True
- query_polygon(vertices, inclusive=False, fact=4)[source]
Returns the pixels whose centers lie within the convex polygon defined by the vertices array (if inclusive is False), or which overlap with this polygon (if inclusive is True).
- Parameters
vertices (float) – Vertex array containing the vertices of the polygon, shape (N, 3).
inclusive (bool) – f False, return the exact set of pixels whose pixels centers lie within the region; if True, return all pixels that overlap with the region.
fact (int) – Only used when inclusive=True. The overlapping test will be done at the resolution fact*nside. For NESTED ordering, fact must be a power of 2, less than 2**30, else it can be any positive integer. Default: 4.
- Returns
The pixels which lie within the given polygon.
- Return type
int array
- query_disc(vec, radius, inclusive=False, fact=4)[source]
- Parameters
vec (float, sequence of 3 elements) – The coordinates of unit vector defining the disk center.
radius (float) – The radius (in radians) of the disk
inclusive (bool) – f False, return the exact set of pixels whose pixels centers lie within the region; if True, return all pixels that overlap with the region.
fact (int) – Only used when inclusive=True. The overlapping test will be done at the resolution fact*nside. For NESTED ordering, fact must be a power of 2, less than 2**30, else it can be any positive integer. Default: 4.
- Returns
The pixels which lie within the given disc.
- Return type
int array
- query_strip(theta1, theta2, inclusive=False)[source]
Returns pixels whose centers lie within the colatitude range defined by theta1 and theta2 (if inclusive is False), or which overlap with this region (if inclusive is True). If theta1<theta2, the region between both angles is considered, otherwise the regions 0<theta<theta2 and theta1<theta<pi.
- Parameters
theta (float) – First colatitude (radians)
phi (float) – Second colatitude (radians)
inclusive (bool) – f False, return the exact set of pixels whose pixels centers lie within the region; if True, return all pixels that overlap with the region.
- Returns
The pixels which lie within the given strip.
- Return type
int array
- boundaries(pix, step=1)[source]
Returns an array containing vectors to the boundary of the nominated pixel.
The returned array has shape (3, 4*step), the elements of which are the x,y,z positions on the unit sphere of the pixel boundary. In order to get vector positions for just the corners, specify step=1.
- plot_grid(ax='mollview', ax_kw={}, step=32, coord='C', **kwargs)[source]
Plot the pixel boundaries of a Healpix grid
- Parameters
ax (WCSAxes or str) – Astropy’s WCSAxes to plot the map. Either an existing instance or the name of a registered projection.
ax_kw (dict) – Extra arguments if a new axes needs to be created.
step (int) – How many points per pixel side
coord (str) – Instrinsic coordinates of the map. Either ‘G’ (Galactic), ‘E’ (Ecliptic) , ‘C’ (Celestial = Equatorial) or any other coordinate frame recognized by astropy.
**kwargs – Passed to matplotlib.pyplot.plot()
- Returns
- The first return value
corresponds to the output
pyplot.plot()
for one of the pixels. The second is the astropy WCSAxes object used.
- Return type
matplotlib.lines.Line2D list, WCSAxes
HealpixMap
- class mhealpy.HealpixMap(data=None, uniq=None, order=None, nside=None, scheme='ring', base=None, density=False, dtype=None)[source]
Bases:
mhealpy.containers.healpix_base.HealpixBase
Object-oriented healpy wrapper with support for multi-resolutions maps (known as multi-order coverage map, or MOC).
You can instantiate a map by providing either:
Size (through
order
ornside
), and ascheme
(‘RING’ or ‘NESTED’). This will initialize an empty map.A list of UNIQ pixels. This will initialize a MOC map. Providing the values for each pixel is optional, zero-initialized by default.
An array (in
data
) and an a scheme (‘RING’ or ‘NESTED’). This will initialize the contents of the single-resolution map.A HealpixBase object. The data will be zero-initialized.
Warning
The initialization input is not validated by default. Consider calling is_mesh_valid() after initialization, otherwise results might be unexpected.
Regardless of the underlaying grid, you can operate on maps using
*
,/
,+
,-
,**
,==
andabs
. For binary operations the result always corresponds to the finest grid, so there is no loss of information. If any of the operands is a MOC, the result is a MOC with an appropiate updated grid.. If both operands have the same NSIDE, the scheme of the result corresponds to the left operand. If you want to preserve the grid for a specific operand, use*=
, /=, etc.The result of most binary operations have the same density parameter as the left-most operand, except the following cases: 1) the product of a density-like maps with a histogram-like results in histogram-like map. 2) the ratio between two histogram-like maps is a density-like map. Operations with scalars leave the map type unchanged.
Warning
Information might degrade if you use in-place operators (e.g.
*=
, /=)The maps are array-like, that is, the can be casted into a regular numpy array (as used by healpy), are iterable (over the pixel values) and can be used with built-in function such as
sum
andmax
.You can also access the value of pixels using regular numpy indexing with
[]
. For MOC maps, no specific pixel ordering is guaranted. For a given pixel numberipix
in the current grid, you can get the corresponding UNIQ pixel number usingm.pix2uniq(ipix)
.- Parameters
data (array) – Values to initialize map. Zero-initialized it not provided. The map NSIDE is deduced from the array size, unless
uniq
is specified in which case this is considered a multi-resolution map.uniq (array or HealpixBase) – List of NUNIQ pixel number to initialize a MOC map.
order (int) – Order of HEALPix map.
nside (int) – Alternatively, you can specify the NSIDE parameter.
scheme (str) – Healpix scheme. Either ‘RING’, ‘NESTED’ or ‘NUNIQ’
base (HealpixBase) – Specify the grid using a HealpixBase object
density (bool) – Whether the value of each pixel should be treated as counts in a histogram (
False
) or as the value of a [density] function evaluated at the center of the pixel (True
). This affect operations involving the splitting of a pixel.dtype (array) – Numpy data type. Will be ignored if data is provided.
- classmethod read_map(filename, field=None, uniq_field=0, hdu=1, density=False)[source]
Read a HEALPix map from a FITS file.
- Parameters
filename (Path) – Path to file
field (int) – Column where the map contents are. Default: 0 for single-resolution maps, 1 for MOC maps.
uniq_field (int) – Column where the UNIQ pixel numbers are. For MOC maps only.
hdu (int) – The header number to look at. Starts at 0.
density (bool) – Whether this is a histogram-like or a density-like map.
- Returns
HealpixMap
- write_map(filename, extra_maps=None, column_names=None, extra_header=None, overwrite=False, coordsys='C')[source]
Write map to disc.
- Parameters
filename (Path) – Path to output file
extra_maps (HealpixMap or array) – Save more maps in the same file as extra columns. Must be conformable.
column_names (str or array) – Name of colums. Must have the same length as the number for maps. Defaults to ‘CONTENTSn’, where
n
is the map number (ommited for a single map). For MOC maps, the pixel information is always stored in the first column, called ‘UNIQ’.coordsys (str) – Coordinate system. Celestial (‘C’), Galactic (‘G’) or Ecliptic (‘E’)
extra_header (iterable) – Iterable of (keyword, value, [comment]) tuples
overwrite (bool) – If True, overwrite the output file if it exists. Raises an OSError if False and the output file exists.
- classmethod adaptive_moc_mesh(max_nside, split_fun, density=False, dtype=None)[source]
Return a zero-initialized MOC map, with an adaptive resolution determined by an arbitrary function.
- Parameters
max_nside (int) – Maximum HEALPix nside to consider
split_fun (function) – This method should return
True
if a pixelorder (should be split into pixel of a higher) –
otherwise. (and False) –
integers (It takes two) –
start (inclusive) and
stop
(exclusive) –a (which correspond to a single pixel in nested rangeset format for) –
max_nside. (map of nside) –
density (bool) – Will be pass to HealpixMap initialization.
dtype (dtype) – Data type
- Returns
HealpixMap
- classmethod moc_from_pixels(nside, pixels, nest=False, density=False, dtype=None)[source]
Return a zero-initialize MOC map where a list of pixels are kept at a given nside, and every other pixel is appropiately downsampled.
Also see the more generic
adaptive_moc_mesh()
.- Parameters
nside (int) – Maximum healpix NSIDE (that is, the NSIDE for the pixel order list)
pixels (array) – Pixels that must be kept at the finest pixelation
nest (bool) – Whether the pixels are a ‘NESTED’ or ‘RING’ scheme
density (bool) – Wheather the map is density-like or histogram-like
dtype – Daty type
- classmethod moc_histogram(nside, samples, max_value, nest=False, weights=None)[source]
Generate an adaptive MOC map by histogramming samples.
If the number of samples is greater than the number of pixels in a map of the input
nside
, consider generating a single-resolution map and then use to_moc().Also see the more generic
adaptive_moc_mesh()
.- Parameters
nside (int) – Healpix NSIDE of the samples and maximum NSIDE of the output map
samples (int array) – List of pixels representing the samples. e.g. the output of healpy.ang2pix().
max_value – maximum number of samples (or sum of weights) per pixel. Note that due to limitations of the input
nside
, the output could contain pixels with a value largen than thisnest (bool) – Whether the samples are in NESTED or RING scheme
weights (array) – Optionally weight the samples. Both must have the same size.
- Returns
HealpixMap
- to_moc(max_value)[source]
Convert a single-resolution map into a MOC based on the maximum value a given pixel the latter should have.
Note
The maximum nside of the MOC map is the same as the nside of the single-resolution map, so the output map could contain pixels with a value greater than this.
If the map is already a MOC map, it will recompute the grid accordingly by combining uniq pixels. Uniq pixels are never split.
Also see the more generic
adaptive_moc_mesh()
.- Parameters
max_value – Maximum value per pixel of the MOC. Whether the map is histogram-like or density-like is taken into account.
- Returns
HealpixMap
- density(density=None, update=True)[source]
Switch between a density-like map and a histogram-like map.
- Parameters
density (bool or None) – Whether the value of each pixel should be treated as counts in a histogram (
False
) or as the value of a [density] function evaluated at the center of the pixel (True
). This affect operations involving the splitting of a pixel.None
will leave this paramter unchanged.update (bool) – If True, the values of the map will be updated accordingly. Otherwise only the density parameter is changed.
Note
The
updat=True
the pixels values are divided/multiplied by their effective number of pixels rather than their solid angle area. In order to achieve this scale the map by 4pi/m.npix.- Returns
The current density
- Return type
bool
- property data
Get the raw data in the form of an array.
- rasterize(nside, scheme)[source]
Convert to map of a given NSIDE and scheme
- Parameters
nside (int) – HEALPix NSIDE
scheme (str) – RING or NESTED
- Returns
HealpixMap
- plot(ax='mollview', ax_kw={}, rasterize=True, coord='C', cbar=True, **kwargs)[source]
Plot map. This is a wrapper for matplotlib.pyplot.imshow
Plots of multi-resolution maps are equivalent to plotting the equivalent rasterized single-resolution map –i.e. values are weighted based on pixel area.
- Parameters
ax (WCSAxes or str) – Astropy’s WCSAxes to plot the map. Either an existing instance or the name of a registered projection.
ax_kw (dict) – Extra arguments if a new axes needs to be created.
rasterize (bool) – If True, the resulting image is equivalent to having called
rasterize()
before plotting. This only affects multi-resolution histogram-like maps.coord (str) – Instrinsic coordinates of the map. Either ‘G’ (Galactic), ‘E’ (Ecliptic) , ‘C’ (Celestial = Equatorial) or any other coordinate frame recognized by astropy.
cbar (bool) – Whether to plot the colorbar.
**kwargs – Passed to matplotlib.pyplot.imshow
- Returns
- The first return value
corresponds to the output
imgshow
. The second is the astropy WCSAxes object used.
- Return type
AxesImage, WCSAxes
- get_interp_val(theta, phi, lonlat=False)[source]
Return the bi-linear interpolation value of a map using 4 nearest neighbours.
For MOC maps, this is equivalent to raterizing the map first to the highest order.
- Parameters
theta (float or array) – Zenith angle (rad)
phi (float or array) – Azimuth angle (rad)
- Returns
scalar or array
- ang2pix(theta, phi, lonlat=False)
Get the pixel (as used in []) that contains a given coordinate
- Parameters
theta (float or array) – Zenith angle
phi (float or arrray) – Azimuth angle
- Returns
int or array
- boundaries(pix, step=1)
Returns an array containing vectors to the boundary of the nominated pixel.
The returned array has shape (3, 4*step), the elements of which are the x,y,z positions on the unit sphere of the pixel boundary. In order to get vector positions for just the corners, specify step=1.
- conformable(other)
For single-resolution maps, return
True
if both maps have the same nside and scheme.For MOC maps, return True if both maps have the same list of UNIQ pixels (including the ordering)
- get_all_neighbours(theta, phi=None, lonlat=False)
Return the 8 nearest pixels. For MOC maps, these might repeat, as this is equivalent to raterizing the maps to the highest order, getting the neighbohrs, and then finding the pixels tha contain them.
- Parameters
theta (float or int or array) – Zenith angle (rad). If phi is
None
, these are assummed to be pixels numbers. For MOC maps, these are assumed to be pixel numbers in NESTED scheme for the equivalent single-resolution map of the highest order.phi (float or array or None) – Azimuth angle (rad)
- Returns
- pixel number of the SW, W, NW, N, NE, E, SE and S neighbours,
shape is (8,) if input is scalar, otherwise shape is (8, N) if input is of length N. If a neighbor does not exist (it can be the case for W, N, E and S) the corresponding pixel number will be -1.
- Return type
array
- get_interp_weights(theta, phi=None, lonlat=False)
Return the 4 closest pixels on the two rings above and below the location and corresponding weights. Weights are provided for bilinear interpolation along latitude and longitude
- Parameters
theta (float or array) – Zenith angle (rad)
phi (float or array) – Azimuth angle (rad)
- Returns
- (pixels, weights), each with of (4,) if the input is scalar,
if (4,N) where N is size of theta and phi. For MOC maps, these pixel numbers might repeate.
- Return type
tuple
- is_mesh_valid()
Return
True
if the map pixelization is valid. For single resolution this simply checks that the size is a valid NSIDE value. For MOC maps, it checks that every point in the sphere is covered by one and only one pixel.- Returns
True
- property is_moc
Return true if this is a Multi-Dimensional Coverage (MOC) map (multi-resolution)
- Returns
bool
- property is_nested
Return true if scheme is NESTED or NUNIQ
- Return
bool
- property is_ring
Return true if scheme is RING
- Return
bool
- nest2pix(pix)
Get the corresponding pixel in the current grid for a pixel in NESTED scheme. For MOC map, return the pixel that contains it.
- Parameters
pix (int or array) – Pixel number in NESTED scheme. Must correspond to a map of the same order as the current.
- Returns
int or array
- property npix
Get number of pixels.
For multi-resolutions maps, this corresponds to the number of utilized UNIQ pixels.
- Returns
int
- property nside
Get map NSIDE
- Returns
int
- property order
Get map order
- Returns
int
- pix2ang(pix, lonlat=False)
Return the coordinates of the center of a pixel
- Parameters
pix (int or array) –
- Returns
(float or array, float or array)
- pix2range(nside, pix)
Get the equivalent range of child pixels in nested scheme for a map of equal or higher nside
- Parameters
nside (int) – Nside of output range sets
pix (int or array) – Pixel numbers
- Returns
- Start pixel (inclusive) and
stop pixel (exclusive)
- Return type
(int or array, int or array)
- pix2uniq(pix)
Get the UNIQ representation of a given pixel index.
- Parameters
pix (int) – Pixel number in the current scheme (as used for [])
- pix2vec(pix)
Return a vector corresponding to the center of a pixel
- Parameters
pix (int or array) –
- Returns
Size (3,N)
- Return type
array
- pix_order_list()
Get a list of lists containing all pixels sorted by order
- Returns
- (
pix_per_order
,nest_pix_per_order
) Each list has a size equal to the map order. Each element is a list of all pixels whose order matches the index of the list position. The first output contains the index of the pixels, while the second contains their coresponding pixel number in a nested scheme.
- (
- Return type
(list, list)
- pix_rangesets(nside=None, argsort=False)
Get the equivalent range of child pixels in nested scheme for a map of equal or higher nside
- Parameters
nside (int or None) – Nside of output range sets. If None, the map nside will be used. nside = mhealpy.MAX_NSIDE returns the cached result
argsort (bool) – Also return also the indices that would sort the array.
- Returns
- With columns named ‘start’ (inclusive) and
’stop’ (exclusive)
- Return type
recarray
- pixarea(pix=0)
Return area of pixel in steradians
- Parameters
pix (int or array) – Pixel number. Only relevant for MOC maps
- Returns
float or array
- plot_grid(ax='mollview', ax_kw={}, step=32, coord='C', **kwargs)
Plot the pixel boundaries of a Healpix grid
- Parameters
ax (WCSAxes or str) – Astropy’s WCSAxes to plot the map. Either an existing instance or the name of a registered projection.
ax_kw (dict) – Extra arguments if a new axes needs to be created.
step (int) – How many points per pixel side
coord (str) – Instrinsic coordinates of the map. Either ‘G’ (Galactic), ‘E’ (Ecliptic) , ‘C’ (Celestial = Equatorial) or any other coordinate frame recognized by astropy.
**kwargs – Passed to matplotlib.pyplot.plot()
- Returns
- The first return value
corresponds to the output
pyplot.plot()
for one of the pixels. The second is the astropy WCSAxes object used.
- Return type
matplotlib.lines.Line2D list, WCSAxes
- query_disc(vec, radius, inclusive=False, fact=4)
- Parameters
vec (float, sequence of 3 elements) – The coordinates of unit vector defining the disk center.
radius (float) – The radius (in radians) of the disk
inclusive (bool) – f False, return the exact set of pixels whose pixels centers lie within the region; if True, return all pixels that overlap with the region.
fact (int) – Only used when inclusive=True. The overlapping test will be done at the resolution fact*nside. For NESTED ordering, fact must be a power of 2, less than 2**30, else it can be any positive integer. Default: 4.
- Returns
The pixels which lie within the given disc.
- Return type
int array
- query_polygon(vertices, inclusive=False, fact=4)
Returns the pixels whose centers lie within the convex polygon defined by the vertices array (if inclusive is False), or which overlap with this polygon (if inclusive is True).
- Parameters
vertices (float) – Vertex array containing the vertices of the polygon, shape (N, 3).
inclusive (bool) – f False, return the exact set of pixels whose pixels centers lie within the region; if True, return all pixels that overlap with the region.
fact (int) – Only used when inclusive=True. The overlapping test will be done at the resolution fact*nside. For NESTED ordering, fact must be a power of 2, less than 2**30, else it can be any positive integer. Default: 4.
- Returns
The pixels which lie within the given polygon.
- Return type
int array
- query_strip(theta1, theta2, inclusive=False)
Returns pixels whose centers lie within the colatitude range defined by theta1 and theta2 (if inclusive is False), or which overlap with this region (if inclusive is True). If theta1<theta2, the region between both angles is considered, otherwise the regions 0<theta<theta2 and theta1<theta<pi.
- Parameters
theta (float) – First colatitude (radians)
phi (float) – Second colatitude (radians)
inclusive (bool) – f False, return the exact set of pixels whose pixels centers lie within the region; if True, return all pixels that overlap with the region.
- Returns
The pixels which lie within the given strip.
- Return type
int array
- property scheme
Return HEALPix scheme
- Returns
Either ‘NESTED’, ‘RING’ or ‘NUNIQ’
- Return type
str
- property uniq
Get an array with the NUNIQ numbers for all pixels
- vec2pix(x, y, z)
Get the pixel (as used in []) that contains a given coordinate
- Parameters
theta (float or array) – Zenith angle
phi (float or arrray) – Azimuth angle
- Returns
int or array