hipscat.pixel_math.hipscat_id#

Compute the hipscat index field

This index is defined as a 64-bit integer which has two parts:
  • healpix pixel (at order 19)

  • incrementing counter (within same healpix, for uniqueness)

|------------------------------------------|-------------------|
|<-----    healpixel at order 19    ------>|<--   counter   -->|

This provides us with an increasing index, that will not overlap between spatially partitioned data files.

See the README in this directory for more information on the fiddly pixel math.

Module Contents#

Functions#

compute_hipscat_id(→ numpy.ndarray)

Compute the hipscat ID field.

_compute_hipscat_id_from_mapped_pixels(mapped_pixels, ...)

hipscat_id_to_healpix(→ numpy.ndarray)

Convert some hipscat ids to the healpix pixel at the specified order

healpix_to_hipscat_id(→ numpy.uint64 | numpy.ndarray)

Convert a healpix pixel to a hipscat_id

Attributes#

HIPSCAT_ID_COLUMN

HIPSCAT_ID_HEALPIX_ORDER

HIPSCAT_ID_MAX

HIPSCAT_ID_COLUMN = '_hipscat_index'[source]#
HIPSCAT_ID_HEALPIX_ORDER = 19[source]#
HIPSCAT_ID_MAX[source]#
compute_hipscat_id(ra_values: List[float], dec_values: List[float]) numpy.ndarray[source]#

Compute the hipscat ID field.

Parameters:
  • ra_values (List[float]) – celestial coordinates, right ascension

  • dec_values (List[float]) – celestial coordinates, declination

Returns:

one-dimensional numpy array of int64s with hipscat IDs for the sky positions.

Raises:

ValueError – if the length of the input lists don’t match.

_compute_hipscat_id_from_mapped_pixels(mapped_pixels, offset_counter)[source]#
hipscat_id_to_healpix(ids: List[int], target_order: int = HIPSCAT_ID_HEALPIX_ORDER) numpy.ndarray[source]#

Convert some hipscat ids to the healpix pixel at the specified order This is just bit-shifting the counter away.

Parameters:
  • ids (List[int64]) – list of well-formatted hipscat ids

  • target_order (int64) – Defaults to HIPSCAT_ID_HEALPIX_ORDER. The order of the pixel to get from the hipscat ids.

Returns:

numpy array of target_order pixels from the hipscat id

healpix_to_hipscat_id(order: int | List[int], pixel: int | List[int], counter: int | List[int] = 0) numpy.uint64 | numpy.ndarray[source]#

Convert a healpix pixel to a hipscat_id

This maps the healpix pixel to the lowest pixel number within that pixel at order 19, then shifts and adds the given counter to get a hipscat_id.

Useful for operations such as filtering by hipscat_id.

Parameters:
  • order (int64 | List[int64]) – order of pixel to convert

  • pixel (int64 | List[int64]) – pixel number in nested ordering of pixel to convert

  • counter (int64 | List[int64]) – 0): counter value in converted hipscat id

Returns:

hipscat id or numpy array of hipscat ids