healsparse modules

HealSparseMap

class healsparse.healSparseMap.HealSparseMap(cov_map=None, cov_index_map=None, sparse_map=None, nside_sparse=None, healpix_map=None, nside_coverage=None, primary=None, sentinel=None, nest=True, metadata=None, _is_view=False)[source]

Bases: object

Class to define a HealSparseMap

__add__(other)[source]

Add a constant.

Cannot be used with recarray maps.

__and__(other)[source]

Perform a bitwise and with a constant.

Cannot be used with recarray maps.

__getitem__(key)[source]

Get part of a healpix map.

__iadd__(other)[source]

Add a constant, in place.

Cannot be used with recarray maps.

__iand__(other)[source]

Perform a bitwise and with a constant, in place.

Cannot be used with recarray maps.

__imul__(other)[source]

Multiply a constant, in place.

Cannot be used with recarray maps.

__init__(cov_map=None, cov_index_map=None, sparse_map=None, nside_sparse=None, healpix_map=None, nside_coverage=None, primary=None, sentinel=None, nest=True, metadata=None, _is_view=False)[source]

Instantiate a HealSparseMap.

Can be created with cov_index_map, sparse_map, and nside_sparse; or with healpix_map, nside_coverage. Also see HealSparseMap.read(), HealSparseMap.make_empty(), HealSparseMap.make_empty_like().

Parameters
  • cov_map (HealSparseCoverage, optional) – Coverage map object

  • cov_index_map (np.ndarray, optional) – Coverage index map, will be deprecated

  • sparse_map (np.ndarray, optional) – Sparse map

  • nside_sparse (int, optional) – Healpix nside for sparse map

  • healpix_map (np.ndarray, optional) – Input healpix map to convert to a sparse map

  • nside_coverage (int, optional) – Healpix nside for coverage map

  • primary (str, optional) – Primary key for recarray, required if dtype has fields.

  • sentinel (int or float, optional) – Sentinel value. Default is hp.UNSEEN for floating-point types, and minimum int for int types.

  • nest (bool, optional) – If input healpix map is in nest format. Default is True.

  • metadata (dict-like, optional) – Map metadata that can be stored in FITS header format.

  • _is_view (bool, optional) – This healSparse map is a view into another healsparse map. Not all features will be available. (Internal usage)

Returns

healSparseMap

Return type

HealSparseMap

__ior__(other)[source]

Perform a bitwise or with a constant, in place.

Cannot be used with recarray maps.

__ipow__(other)[source]

Divide a constant, in place.

Cannot be used with recarray maps.

__isub__(other)[source]

Subtract a constant, in place.

Cannot be used with recarray maps.

__itruediv__(other)[source]

Divide a constant, in place.

Cannot be used with recarray maps.

__ixor__(other)[source]

Perform a bitwise xor with a constant, in place.

Cannot be used with recarray maps.

__mul__(other)[source]

Multiply a constant.

Cannot be used with recarray maps.

__or__(other)[source]

Perform a bitwise or with a constant.

Cannot be used with recarray maps.

__pow__(other)[source]

Raise the map to a power.

Cannot be used with recarray maps.

__repr__()[source]

Return repr(self).

__setitem__(key, value)[source]

Set part of a healpix map

__str__()[source]

Return str(self).

__sub__(other)[source]

Subtract a constant.

Cannot be used with recarray maps.

__truediv__(other)[source]

Divide a constant.

Cannot be used with recarray maps.

__weakref__

list of weak references to the object (if defined)

__xor__(other)[source]

Perform a bitwise xor with a constant.

Cannot be used with recarray maps.

apply_mask(mask_map, mask_bits=None, mask_bit_arr=None, in_place=True)[source]

Apply an integer mask to the map. All pixels in the integer mask that have any bits in mask_bits set will be zeroed in the output map. The default is that this operation will be done in place, but it may be set to return a copy with a masked map.

Parameters
  • mask_map (HealSparseMap) – Integer mask to apply to the map.

  • mask_bits (int, optional) – Bits to be treated as bad in the mask_map. Default is None (all non-zero pixels are masked)

  • mask_bit_arr (list or np.ndarray, optional) – Array of bit values, used if mask_map is a wide_mask_map.

  • in_place (bool, optional) – Apply operation in place. Default is True

Returns

masked_map – self if in_place is True, a new copy otherwise

Return type

HealSparseMap

astype(dtype, sentinel=None)[source]

Convert sparse map to a different numpy datatype, including sentinel values. If sentinel is not specified the default for the converted datatype is used (healpy.UNSEEN for float, and -MAXINT for ints).

Parameters
  • dtype (numpy.dtype) – Valid numpy dtype for a single array.

  • sentinel (int or float, optional) – Converted map sentinel value.

Returns

sparse_map – New map with new data type.

Return type

HealSparseMap

check_bits_pix(pixels, bits, nest=True)[source]

Check the bits at the map for a set of pixels.

Parameters
  • pixel (np.ndarray) – Integer array of healpix pixels.

  • nest (bool, optional) – Are the pixels in nest scheme? Default is True.

  • bits (list) – List of bits to check

Returns

bit_flags – Array of np.bool flags on whether any of the input bits were set

Return type

np.ndarray

check_bits_pos(theta_or_ra, phi_or_dec, bits, lonlat=True)[source]

Check the bits at the map for an array of positions. Positions may be theta/phi co-latitude and longitude in radians, or longitude and latitude in degrees.

Parameters
  • theta_or_ra (float, array-like) – Angular coordinates of points on a sphere.

  • phi_or_dec (float, array-like) – Angular coordinates of points on a sphere.

  • lonlat (bool, optional) – If True, input angles are longitude and latitude in degrees. Otherwise, they are co-latitude and longitude in radians.

  • bits (list) – List of bits to check

Returns

bit_flags – Array of np.bool flags on whether any of the input bits were set

Return type

np.ndarray

clear_bits_pix(pixels, bits, nest=True)[source]

Clear bits of a wide_mask map.

Parameters
  • pixels (np.ndarray) – Integer array of sparse_map pixel values

  • bits (list) – List of bits to clear

static convert_healpix_map(healpix_map, nside_coverage, nest=True, sentinel=- 1.6375e+30)[source]

Convert a healpix map to a healsparsemap.

Parameters
  • healpix_map (np.ndarray) – Numpy array that describes a healpix map.

  • nside_coverage (int) – Nside for the coverage map to construct

  • nest (bool, optional) – Is the input map in nest format? Default is True.

  • sentinel (float, optional) – Sentinel value for null values in the sparse_map. Default is hp.UNSEEN

Returns

  • cov_map (HealSparseCoverage) – Coverage map with pixel indices

  • sparse_map (np.ndarray) – Sparse map of input values.

property coverage_map

Get the fractional area covered by the sparse map in the resolution of the coverage map

Returns

cov_map – Float array of fractional coverage of each pixel

Return type

np.ndarray

property coverage_mask

Get the boolean mask of the coverage map.

Returns

cov_mask – Boolean array of coverage mask.

Return type

np.ndarray

degrade(nside_out, reduction='mean', weights=None)[source]

Reduce the resolution, i.e., increase the pixel size of a given sparse map.

Parameters
  • nside_out (int) – Output Nside resolution parameter.

  • reduction (str) – Reduction method (mean, median, std, max, min, and, or, sum, prod, wmean).

  • weights (healSparseMap) – If the reduction is wmean this is the map with the weights to use. It should have the same characteristics as the original map.

Returns

healSparseMap – New map, at the desired resolution.

Return type

HealSparseMap

property dtype

get the dtype of the map

fracdet_map(nside)[source]

Get the fractional area covered by the sparse map at an arbitrary resolution. This output fracdet_map counts the fraction of “valid” sub-pixels (those that are not equal to the sentinel value) at the desired nside resolution.

Note: You should not compute the fracdet_map of an existing fracdet_map. To get a fracdet_map at a lower resolution, use the degrade method with the default “mean” reduction.

Parameters

nside (int) – Healpix nside for fracdet map. Must not be greater than sparse resolution or less than coverage resolution.

Returns

fracdet_map – Fractional coverage map.

Return type

HealSparseMap

generate_healpix_map(nside=None, reduction='mean', key=None, nest=True)[source]

Generate the associated healpix map

if nside is specified, then reduce to that nside

Parameters
  • nside (int) – Output nside resolution parameter (should be a multiple of 2). If not specified the output resolution will be equal to the parent’s sparsemap nside_sparse

  • reduction (str) – If a change in resolution is requested, this controls the method to reduce the map computing the “mean”, “median”, “std”, “max”, “min”, “sum” or “prod” (product) of the neighboring pixels to compute the “degraded” map.

  • key (str) – If the parent HealSparseMap contains recarrays, key selects the field that will be transformed into a HEALPix map.

  • nest (bool, optional) – Output healpix map should be in nest format?

Returns

hp_map – Output HEALPix map with the requested resolution.

Return type

np.ndarray

get_single(key, sentinel=None, copy=False)[source]

Get a single healpix map out of a recarray map, with the ability to override a sentinel value.

Parameters
  • key (str) – Field for the recarray

  • sentinel (int or float or None, optional) – Override the default sentinel value. Default is None (use default)

get_values_pix(pixels, nest=True, valid_mask=False)[source]

Get the map value for a set of pixelx.

Parameters
  • pixel (np.ndarray) – Integer array of healpix pixels.

  • nest (bool, optional) – Are the pixels in nest scheme? Default is True.

  • valid_mask (bool, optional) – Return mask of True/False instead of values

Returns

values – Array of values/validity from the map.

Return type

np.ndarray

get_values_pos(theta_or_ra, phi_or_dec, lonlat=True, valid_mask=False)[source]

Get the map value for the position. Positions may be theta/phi co-latitude and longitude in radians, or longitude and latitude in degrees.

Parameters
  • theta_or_ra (float, array-like) – Angular coordinates of points on a sphere.

  • phi_or_dec (float, array-like) – Angular coordinates of points on a sphere.

  • lonlat (bool, optional) – If True, input angles are longitude and latitude in degrees. Otherwise, they are co-latitude and longitude in radians.

  • valid_mask (bool, optional) – Return mask of True/False instead of values

Returns

values – Array of values/validity from the map.

Return type

np.ndarray

property is_integer_map

Check that the map is an integer map

Returns

is_integer_map

Return type

bool

property is_rec_array

Check that the map is a recArray map.

Returns

is_rec_array

Return type

bool

property is_unsigned_map

Check that the map is an unsigned integer map

Returns

is_unsigned_map

Return type

bool

property is_wide_mask_map

Check that the map is a wide mask

Returns

is_wide_mask_map

Return type

bool

classmethod make_empty(nside_coverage, nside_sparse, dtype, primary=None, sentinel=None, wide_mask_maxbits=None, metadata=None, cov_pixels=None)[source]

Make an empty map with nothing in it.

Parameters
  • nside_coverage (int) – Nside for the coverage map

  • nside_sparse (int) – Nside for the sparse map

  • dtype (str or list or np.dtype) – Datatype, any format accepted by numpy.

  • primary (str, optional) – Primary key for recarray, required if dtype has fields.

  • sentinel (int or float, optional) – Sentinel value. Default is hp.UNSEEN for floating-point types, and minimum int for int types.

  • wide_mask_maxbits (int, optional) – Create a “wide bit mask” map, with this many bits.

  • metadata (dict-like, optional) – Map metadata that can be stored in FITS header format.

  • cov_pixels (np.ndarray or list) – List of integer coverage pixels to pre-allocate

Returns

healSparseMap – HealSparseMap filled with sentinel values.

Return type

HealSparseMap

classmethod make_empty_like(sparsemap, nside_coverage=None, nside_sparse=None, dtype=None, primary=None, sentinel=None, wide_mask_maxbits=None, metadata=None, cov_pixels=None)[source]

Make an empty map with the same parameters as an existing map.

Parameters
  • sparsemap (HealSparseMap) – Sparse map to use as basis for new empty map.

  • nside_coverage (int, optional) – Coverage nside, default to sparsemap.nside_coverage

  • nside_sparse (int, optional) – Sparse map nside, default to sparsemap.nside_sparse

  • dtype (str or list or np.dtype, optional) – Datatype, any format accepted by numpy. Default is sparsemap.dtype

  • primary (str, optional) – Primary key for recarray. Default is sparsemap.primary

  • sentinel (int or float, optional) – Sentinel value. Default is sparsemap._sentinel

  • wide_mask_maxbits (int, optional) – Create a “wide bit mask” map, with this many bits.

  • metadata (dict-like, optional) – Map metadata that can be stored in FITS header format.

  • cov_pixels (np.ndarray or list) – List of integer coverage pixels to pre-allocate

Returns

healSparseMap – HealSparseMap filled with sentinel values.

Return type

HealSparseMap

property metadata

Return the metadata dict.

Returns

metadata

Return type

dict

property nside_coverage

Get the nside of the coverage map

Returns

nside_coverage

Return type

int

property nside_sparse

Get the nside of the sparse map

Returns

nside_sparse

Return type

int

property primary

Get the primary field

Returns

primary

Return type

str

classmethod read(filename, nside_coverage=None, pixels=None, header=False, degrade_nside=None, weightfile=None, reduction='mean')[source]

Read in a HealSparseMap.

Parameters
  • filename (str) – Name of the file to read. May be either a regular HEALPIX map or a HealSparseMap

  • nside_coverage (int, optional) – Nside of coverage map to generate if input file is healpix map.

  • pixels (list, optional) – List of coverage map pixels to read. Only used if input file is a HealSparseMap

  • header (bool, optional) – Return the fits header as well as map? Default is False.

  • degrade_nside (int, optional) – Degrade map to this nside on read. None means leave as-is.

  • weightfile (str, optional) – Floating-point map to supply weights for degrade wmean. Must be a HealSparseMap (weighted degrade not supported for healpix degrade-on-read).

  • reduction (str, optional) – Reduction method with degrade-on-read. (mean, median, std, max, min, and, or, sum, prod, wmean).

Returns

  • healSparseMap (HealSparseMap) – HealSparseMap from file, covered by pixels

  • header (fitsio.FITSHDR or astropy.io.fits (if header=True)) – Fits header for the map file.

set_bits_pix(pixels, bits, nest=True)[source]

Set bits of a wide_mask map.

Parameters
  • pixels (np.ndarray) – Integer array of sparse_map pixel values

  • bits (list) – List of bits to set

update_values_pix(pixels, values, nest=True, operation='replace')[source]

Update the values in the sparsemap for a list of pixels.

Parameters
  • pixels (np.ndarray) – Integer array of sparse_map pixel values

  • values (np.ndarray) – Value or Array of values. Must be same type as sparse_map

  • operation (str, optional) – Operation to use to update values. May be ‘replace’ (default), ‘or’, or ‘and’ (for bit masks)

property valid_pixels

Get an array of valid pixels in the sparse map.

Returns

valid_pixels

Return type

np.ndarray

valid_pixels_pos(lonlat=True, return_pixels=False)[source]

Get an array with the position of valid pixels in the sparse map.

Parameters
  • lonlat (bool, optional) – If True, input angles are longitude and latitude in degrees. Otherwise, they are co-latitude and longitude in radians.

  • return_pixels (bool, optional) – If true, return valid_pixels / co-lat / co-lon or valid_pixels / lat / lon instead of lat / lon

Returns

positions – By default it will return a tuple of the form (theta, phi) in radians unless lonlat = True, for which it will return (ra, dec) in degrees. If return_pixels = True, valid_pixels will be returned as first element in tuple.

Return type

tuple

property wide_mask_maxbits

Get the maximum number of bits stored in the wide mask.

Returns

wide_mask_maxbits – Maximum number of bits. 0 if not wide mask.

Return type

int

property wide_mask_width

Get the width of the wide mask

Returns

wide_mask_width – Width of wide mask array. 0 if not wide mask.

Return type

int

write(filename, clobber=False, nocompress=False)[source]

Write heal HealSparseMap to filename. Use the metadata property from the map to persist additional information in the fits header.

Parameters
  • filename (str) – Name of file to save

  • clobber (bool, optional) – Clobber existing file? Default is False.

  • nocompress (bool, optional) – If this is False, then integer maps will be compressed losslessly. Note that np.int64 maps cannot be compressed in the FITS standard.

HealSparseRandoms

healsparse.healSparseRandoms.make_uniform_randoms(sparse_map, n_random, rng=None)[source]

Make an array of uniform randoms.

Parameters
  • sparse_map (healsparse.HealSparseMap) – Sparse map object

  • n_random (int) – Number of randoms to generate

  • rng (np.random.RandomState, optional) – Pre-set Random number generator. Default is None.

Returns

  • ra_array (np.ndarray) – Float array of RAs (degrees)

  • dec_array (np.ndarray) – Float array of declinations (degrees)

healsparse.healSparseRandoms.make_uniform_randoms_fast(sparse_map, n_random, nside_randoms=8388608, rng=None)[source]

Make an array of uniform randoms.

Parameters
  • sparse_map (healsparse.HealSparseMap) – Sparse map object

  • n_random (int) – Number of randoms to generate

  • nside_randoms (int, optional) – Nside for pixel centers to select random points

  • rng (np.random.RandomState, optional) – Pre-set Random number generator. Default is None.

Returns

  • ra_array (np.ndarray) – Float array of RAs (degrees)

  • dec_array (np.ndarray) – Float array of declinations (degrees)

Operations

healsparse.operations.and_intersection(map_list)[source]

Bitwise or a list of HealSparseMaps as an intersection. Only pixels that are valid in all the input maps will have valid values in the output. Only works on integer maps.

Parameters

map_list (list of HealSparseMap) – Input list of maps to bitwise and

Returns

result – Bitwise and of maps

Return type

HealSparseMap

healsparse.operations.and_union(map_list)[source]

Bitwise and a list of HealSparseMaps as a union. Empty values will be treated as 0s in the bitwise and, and the output map will have a union of all the input map pixels. Only works in integer maps.

Parameters

map_list (list of HealSparseMap) – Input list of maps to bitwise and

Returns

result – Bitwise and of maps

Return type

HealSparseMap

healsparse.operations.max_intersection(map_list)[source]

Element-wise maximum of the intersection of a list of the HealSparseMaps.

Parameters

map_list (list of HealSparseMap) – Input list of maps to compute the maximum of

Returns

result – Element-wise maximum of maps

Return type

HealSparseMap

healsparse.operations.max_union(map_list)[source]

Element-wise maximum of the union of a list of HealSparseMaps.

Parameters

map_list (list of HealSparseMap) – Input list of maps to compute the maximum of

Returns

result – Element-wise maximum of maps

Return type

HealSparseMap

healsparse.operations.min_intersection(map_list)[source]

Element-wise minimum of the intersection of a list of HealSparseMaps.

Parameters

map_list (list of HealSparseMap) – Input list of maps to compute the minimum of

Returns

result – Element-wise minimum of maps

Return type

HealSparseMap

healsparse.operations.min_union(map_list)[source]

Element-wise minimum of the union of a list of HealSparseMaps.

Parameters

map_list (list of HealSparseMap) – Input list of maps to compute the minimum of

Returns

result – Element-wise minimum of maps

Return type

HealSparseMap

healsparse.operations.or_intersection(map_list)[source]

Bitwise or a list of HealSparseMaps as an intersection. Only pixels that are valid in all the input maps will have valid values in the output. Only works on integer maps.

Parameters

map_list (list of HealSparseMap) – Input list of maps to bitwise or

Returns

result – Bitwise or of maps

Return type

HealSparseMap

healsparse.operations.or_union(map_list)[source]

Bitwise or a list of HealSparseMaps as a union. Empty values will be treated as 0s in the bitwise or, and the output map will have a union of all the input map pixels. Only works in integer maps.

Parameters

map_list (list of HealSparseMap) – Input list of maps to bitwise or

Returns

result – Bitwise or of maps

Return type

HealSparseMap

healsparse.operations.product_intersection(map_list)[source]

Compute the product of a list of HealSparseMaps as an intersection. Only pixels that are valid in all the input maps will have valid values in the output.

Parameters

map_list (list of HealSparseMap) – Input list of maps to take the product

Returns

result – Product of maps

Return type

HealSparseMap

healsparse.operations.product_union(map_list)[source]

Compute the product of a list of HealSparseMaps as a union. Empty values will be treated as 1s in the product, and the output map will have a union of all the input map pixels.

Parameters

map_list (list of HealSparseMap) – Input list of maps to take the product

Returns

result – Product of maps

Return type

HealSparseMap

healsparse.operations.sum_intersection(map_list)[source]

Sum a list of HealSparseMaps as an intersection. Only pixels that are valid in all the input maps will have valid values in the output.

Parameters

map_list (list of HealSparseMap) – Input list of maps to sum

Returns

result – Summation of maps

Return type

HealSparseMap

healsparse.operations.sum_union(map_list)[source]

Sum a list of HealSparseMaps as a union. Empty values will be treated as 0s in the summation, and the output map will have a union of all the input map pixels.

Parameters

map_list (list of HealSparseMap) – Input list of maps to sum

Returns

result – Summation of maps

Return type

HealSparseMap

healsparse.operations.ufunc_intersection(map_list, func, filler_value=0)[source]

Apply numpy ufunc to the intersection of a list of HealSparseMaps.

Parameters
  • map_list (list of HealSparseMap) – Input list of maps where the operation is applied

  • func (np.ufunc) – Numpy universal function to apply

  • filler_value (int or float) – Starting value

Returns

result – Resulting map

Return type

HealSparseMap

healsparse.operations.ufunc_union(map_list, func, filler_value=0)[source]

Apply numpy ufunc to the union of a list of HealSparseMaps.

Parameters
  • map_list (list of HealSparseMaps) – Input list of maps where the operation is applied

  • func (np.ufunc) – Numpy universal function to apply

  • filler_value (int or float) – Starting value and filler for the union

Returns

result – Resulting map

Return type

HealSparseMap

healsparse.operations.xor_intersection(map_list)[source]

Bitwise xor a list of HealSparseMaps as an intersection. Only pixels that are valid in all the input maps will have valid values in the output. Only works on integer maps.

Parameters

map_list (list of HealSparseMap) – Input list of maps to bitwise xor

Returns

result – Bitwise xor of maps

Return type

HealSparseMap

healsparse.operations.xor_union(map_list)[source]

Bitwise xor a list of HealSparseMaps as a union. Empty values will be treated as 0s in the bitwise or, and the output map will have a union of all the input map pixels. Only works in integer maps.

Parameters

map_list (list of HealSparseMap) – Input list of maps to bitwise xor

Returns

result – Bitwise xor of maps

Return type

HealSparseMap

Geometry Library

class healsparse.geom.Circle(*, ra, dec, radius, value)[source]

Bases: healsparse.geom.GeomBase

__init__(*, ra, dec, radius, value)[source]
Parameters
  • ra (float) – ra in degrees

  • dec (float) – dec in degrees

  • radius (float) – radius in degrees

  • value (number) – Value for pixels in the map

__repr__()[source]

Return repr(self).

property dec

get the dec value

get_pixels(*, nside)[source]

get the pixels associated with this circle

Parameters

nside (int) – Nside for the pixels

property ra

get the ra value

property radius

get the radius value

class healsparse.geom.GeomBase[source]

Bases: object

base class for goemetric objects that can convert themselves to maps

__weakref__

list of weak references to the object (if defined)

get_map(*, nside_coverage, nside_sparse, dtype, wide_mask_maxbits=None)[source]

get a healsparse map corresponding to this geometric primitive

Parameters
  • nside_coverage (int) – nside of coverage map

  • nside_sparse (int) – nside of sparse map

  • dtype (np.dtype) – dtype of the output array

  • wide_mask_maxbits (int, optional) – Create a “wide bit mask” map, with this many bits.

Returns

Return type

HealSparseMap

get_map_like(sparseMap)[source]

Get a healsparse map corresponding to this geometric primitive, with the same parameters as an input sparseMap.

Parameters

sparseMap (healsparse.HealSparseMap) – Input map to match parameters

Returns

Return type

HealSparseMap

get_pixels(*, nside)[source]

get pixels for this map

Parameters

nside (int) – Nside for the pixels

property is_integer_value

Check if the value is an integer type

property value

get the value to be used for all pixels in the map

class healsparse.geom.Polygon(*, ra, dec, value)[source]

Bases: healsparse.geom.GeomBase

__init__(*, ra, dec, value)[source]

represent a polygon

both counter clockwise and clockwise order for polygon vertices works

Parameters
  • ra (array) – ra of vertices in degrees, size [nvert]

  • dec (array) – dec of vertices in degrees, size [nvert]

  • value (number) – Value for pixels in the map

__repr__()[source]

Return repr(self).

property dec

get the dec value

get_pixels(*, nside)[source]

get the pixels associated with this polygon

Parameters

nside (int) – Nside for the pixels

property ra

get the ra value

property vertices

get the dec value

healsparse.geom.realize_geom(geom, smap, type='or')[source]

Realize geometry objects in a map

Parameters
  • geom (geometric primitive or list thereof) – List of Geom objects, e.g. Circle, Polygon

  • smap (HealSparseMaps) – The map in which to realize the objects

  • type (string) – Way to combine the list of geometric objects. Default is to “or” them

HealSparseCoverage

class healsparse.healSparseCoverage.HealSparseCoverage(cov_index_map, nside_sparse)[source]

Bases: object

Class to define a HealSparseCoverage map

__init__(cov_index_map, nside_sparse)[source]

Instantiate a HealSparseCoverage map.

Returns

healSparseCoverage

Return type

HealSparseCoverage

__weakref__

list of weak references to the object (if defined)

append_pixels(sparse_map_size, new_cov_pix, check=True, copy=True)[source]

Append new pixels to the coverage map

Parameters
  • sparse_map_size (int) – Size of current sparse map

  • new_cov_pix (np.ndarray) – Array of new coverage pixels

property bit_shift

Get the bit_shift for the coverage map

Returns

bit_shift – Number of bits to shift from coarse to fine maps

Return type

int

cov_pixels(sparse_pixels)[source]

Get coverage pixel numbers (nest) from a set of sparse pixels.

Parameters

sparse_pixels (np.ndarray) – Array of sparse pixels

Returns

cov_pixels – Coverage pixel numbers (nest format)

Return type

np.ndarray

cov_pixels_from_index(index)[source]

Get the coverage pixels from the sparse map index.

Parameters

index (np.ndarray) – Array of indices in sparse map

Returns

cov_pixels – Coverage pixel numbers (nest format)

Return type

np.ndarray

property coverage_mask

Get the boolean mask of the coverage map.

Returns

cov_mask – Boolean array of coverage mask.

Return type

np.ndarray

initialize_pixels(cov_pix)[source]

Initialize pixels in the index map

Parameters

cov_pix (np.ndarray) – Array of coverage pixels

classmethod make_empty(nside_coverage, nside_sparse)[source]

Make an empty coverage map.

Parameters
  • nside_coverage (int) – Healpix nside for the coverage map

  • nside_sparse (int) – Healpix nside for the sparse map

Returns

healSparseCoverage – HealSparseCoverage from file

Return type

HealSparseCoverage

classmethod make_from_pixels(nside_coverage, nside_sparse, cov_pixels)[source]

Make an empty coverage map.

Parameters
  • nside_coverage (int) – Healpix nside for the coverage map

  • nside_sparse (int) – Healpix nside for the sparse map

  • cov_pixels (np.ndarray) – Array of coverage pixels

Returns

healSparseCoverage – HealSparseCoverage from file

Return type

HealSparseCoverage

property nfine_per_cov

Get the number of fine (sparse) pixels per coarse (coverage) pixel

Returns

nfine_per_cov – Number of fine (sparse) pixels per coverage pixel

Return type

int

property nside_coverage

Get the nside of the coverage map

Returns

nside_coverage

Return type

int

property nside_sparse

Get the nside of the associated sparse map

Returns

nside_sparse

Return type

int

classmethod read(filename_or_fits)[source]

Read in HealSparseCoverage.

Returns

healSparseCoverage – HealSparseCoverage from file

Return type

HealSparseCoverage

Concatenation

healsparse.cat_healsparse_files.cat_healsparse_files(file_list, outfile, check_overlap=False, clobber=False, in_memory=False, nside_coverage_out=None, or_overlap=False)[source]

Concatenate healsparse files together in a memory-efficient way.

Parameters
  • file_list (list of str) – List of filenames to concatenate

  • outfile (str) – Output filename

  • check_overlap (bool, optional) – Check that each file has a unique sparse map. This may be slower.

  • clobber (bool, optional) – Clobber existing outfile

  • in_memory (bool, optional) – Do operations in-memory (required unless fitsio is available).

  • nside_coverage_out (int, optional) – Output map with specific nside_coverage. Default is nside_coverage of first map in file_list.

  • or_overlap (bool, optional) – If True compute the or overlap of two integer maps when concatenating.