hipscat.pixel_math#

Utilities for performing fun math on healpix pixels

Submodules#

Package Contents#

Classes#

HealpixPixel

A HEALPix pixel, represented by an order and pixel number in NESTED ordering scheme

Functions#

get_healpix_pixel(...)

Function to convert argument of either HealpixPixel or a tuple of (order, pixel) to a

compute_hipscat_id(→ numpy.ndarray)

Compute the hipscat ID field.

hipscat_id_to_healpix(→ numpy.ndarray)

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

check_margin_bounds(r_asc, dec, pixel_order, pixel, ...)

Check a set of points in ra and dec space to see if they fall within the

compute_pixel_map(histogram[, highest_order, ...])

Generate alignment from high order pixels to those of equal or lower order

empty_histogram(highest_order)

Use numpy to create an histogram array with the right shape, filled with zeros.

generate_alignment(histogram[, highest_order, ...])

Generate alignment from high order pixels to those of equal or lower order

generate_constant_pixel_map(histogram, ...)

Special case of creating a destination pixel map where you want to

generate_destination_pixel_map(histogram, pixel_map)

Generate mapping from destination pixel to all the constituent pixels.

generate_histogram(data, highest_order[, ra_column, ...])

Generate a histogram of counts for objects found in data

get_edge(d_order, pix, edge)

Get all the pixels at order kk+dk bordering pixel pix.

get_margin(order, pix, d_order)

Get all the pixels at order order+dk bordering pixel pix.

pixel_is_polar(order, pix)

Checks if a healpixel is a polar pixel.

Attributes#

HealpixInputTypes

class HealpixPixel[source]#

A HEALPix pixel, represented by an order and pixel number in NESTED ordering scheme

see https://lambda.gsfc.nasa.gov/toolbox/pixelcoords.html for more information

property dir: int#

Directory number for the pixel.

This is necessary for file systems that limit to 10,000 subdirectories. The directory name will take the HiPS standard form of:

<catalog_base_dir>/Norder=<pixel_order>/Dir=<directory number>

Where the directory number is calculated using integer division as:

(pixel_number/10000)*10000
order: int#
pixel: int#
__post_init__() None[source]#

Initialize a HEALPix pixel :param order: HEALPix order :param pixel: HEALPix pixel number in NESTED ordering scheme

__str__() str[source]#

Return str(self).

__repr__()[source]#

Return repr(self).

__getitem__(key: int) int[source]#
convert_to_lower_order(delta_order: int) HealpixPixel[source]#

Returns the HEALPix pixel that contains the pixel at a lower order

Parameters:

delta_order – the difference in order to be subtracted from the current order to generate the pixel at a lower order. Must be non-negative

Returns:

A new HealpixPixel at the current order - delta_order which contains the current pixel

Raises:

ValueError – If delta_order is greater than the current order, a pixel cannot be generated at a negative order. Or if delta_order is negative

convert_to_higher_order(delta_order: int) List[HealpixPixel][source]#

Returns a list of HEALPix pixels making up the current pixel at a higher order

Parameters:

delta_order – the difference in order to be added to the current order to generate the pixels at a higher order. Must be non-negative

Returns:

A new HealpixPixel at the current order - delta_order which contains the current pixel

Raises:

ValueError – If delta_order + current order is greater than the maximum HEALPix order, pixels cannot be generated. Or if delta_order is negative

HealpixInputTypes[source]#
get_healpix_pixel(pixel: HealpixInputTypes) hipscat.pixel_math.healpix_pixel.HealpixPixel[source]#

Function to convert argument of either HealpixPixel or a tuple of (order, pixel) to a HealpixPixel

Parameters:

pixel – an object to be converted to a HealpixPixel object

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.

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

check_margin_bounds(r_asc, dec, pixel_order, pixel, margin_threshold, step=100, chunk_size=1000)[source]#
Check a set of points in ra and dec space to see if they fall within the

margin_threshold of a given healpixel.

Parameters:
  • r_asc (np.array) – one dimmensional array representing the ra of a given set of points.

  • dec (np.array) – one dimmensional array representing the dec of a given set of points.

  • pixel_order (int) – the order of the healpixel to be checked against.

  • pixel (int) – the healpixel to be checked against.

  • margin_threshold (double) – the size of the margin cache in arcseconds.

  • step (int) – the amount of samples we’ll get from the healpixel boundary to test the datapoints against. the higher the step the more granular the point checking and therefore the more accurate, however using a smaller step value might be helpful for super large datasets to save compute time if accuracy is less important.

  • chunk_size (int) – the size of the chunk of data points we process on a single thread at any given time. We only process a certain amount of the data at once as trying to check all data points at once can lead to memory issues.

Returns:

numpy.array of boolean values corresponding to the ra and dec coordinates checked against whether a given point is within any of the provided polygons.

compute_pixel_map(histogram, highest_order=10, lowest_order=0, threshold=1000000)[source]#

Generate alignment from high order pixels to those of equal or lower order

We may initially find healpix pixels at order 10, but after aggregating up to the pixel threshold, some final pixels are order 4 or 7. This method provides a map from destination healpix pixels at a lower order to the origin pixels at the highest order. This may be used as an input into later partitioning map reduce steps.

Parameters:
  • histogram (np.array) – one-dimensional numpy array of long integers where the value at each index corresponds to the number of objects found at the healpix pixel.

  • highest_order (int) – the highest healpix order (e.g. 0-10)

  • lowest_order (int) – the lowest healpix order (e.g. 1-5). specifying a lowest order constrains the partitioning to prevent spatially large pixels.

  • threshold (int) – the maximum number of objects allowed in a single pixel

Returns:

dictionary that maps the HealpixPixel (a pixel at destination order) to a tuple of origin pixel information.

The tuple contains the following:

  • 0 - the total number of rows found in this destination pixel

  • 1 - the set of indexes in histogram for the pixels at the original healpix order.

Raises:

ValueError – if the histogram is the wrong size, or some initial histogram bins exceed threshold.

empty_histogram(highest_order)[source]#

Use numpy to create an histogram array with the right shape, filled with zeros.

Parameters:

highest_order (int) – the highest healpix order (e.g. 0-10)

Returns:

one-dimensional numpy array of long integers, where the length is equal to the number of pixels in a healpix map of target order, and all values are set to 0.

generate_alignment(histogram, highest_order=10, lowest_order=0, threshold=1000000)[source]#

Generate alignment from high order pixels to those of equal or lower order

We may initially find healpix pixels at order 10, but after aggregating up to the pixel threshold, some final pixels are order 4 or 7. This method provides a map from pixels at order 10 to their destination pixel. This may be used as an input into later partitioning map reduce steps.

Parameters:
  • histogram (np.array) – one-dimensional numpy array of long integers where the value at each index corresponds to the number of objects found at the healpix pixel.

  • highest_order (int) – the highest healpix order (e.g. 5-10)

  • lowest_order (int) – the lowest healpix order (e.g. 1-5). specifying a lowest order constrains the partitioning to prevent spatially large pixels.

  • threshold (int) – the maximum number of objects allowed in a single pixel

Returns:

one-dimensional numpy array of integer 3-tuples, where the value at each index corresponds to the destination pixel at order less than or equal to the highest_order.

The tuple contains three integers:

  • order of the destination pixel

  • pixel number at the above order

  • the number of objects in the pixel

Raises:

ValueError – if the histogram is the wrong size, or some initial histogram bins exceed threshold.

generate_constant_pixel_map(histogram, constant_healpix_order)[source]#

Special case of creating a destination pixel map where you want to preserve some constant healpix order across the sky.

NB:
  • this method filters out pixels with no counts

  • no upper thresholds are applied, and single pixel may contain many rows

Parameters:
  • histogram (np.array) – one-dimensional numpy array of long integers where the value at each index corresponds to the number of objects found at the healpix pixel.

  • constant_healpix_order (int) – the desired healpix order (e.g. 0-10)

Returns:

dictionary that maps non-empty bins as HealpixPixel to a tuple of origin pixel information.

The tuple contains the following:

  • 0 - the total number of rows found in this destination pixel, same as the origin bin

  • 1 - the set of indexes in histogram for the pixels at the original healpix order, which will be a list containing only the origin pixel.

Raises:

ValueError – if the histogram is the wrong size

generate_destination_pixel_map(histogram, pixel_map)[source]#

Generate mapping from destination pixel to all the constituent pixels.

Parameters:
  • histogram (np.array) – one-dimensional numpy array of long integers where the value at each index corresponds to the number of objects found at the healpix pixel.

  • pixel_map (np.array) – one-dimensional numpy array of integer 3-tuples. See generate_alignment() for more details on this format.

Returns:

dictionary that maps the integer 3-tuple of a pixel at destination order to the set of indexes in histogram for the pixels at the original healpix order

generate_histogram(data: pandas.DataFrame, highest_order, ra_column='ra', dec_column='dec')[source]#

Generate a histogram of counts for objects found in data

Parameters:
  • data (pd.DataFrame) – tabular object data

  • highest_order (int) – the highest healpix order (e.g. 0-10)

  • ra_column (str) – where in the input to find the celestial coordinate, right ascension

  • dec_column (str) – where in the input to find the celestial coordinate, declination

Returns:

one-dimensional numpy array of long integers where the value at each index corresponds to the number of objects found at the healpix pixel.

Raises:

ValueError – if the ra_column or dec_column cannot be found in the input file.

get_edge(d_order, pix, edge)[source]#

Get all the pixels at order kk+dk bordering pixel pix. See hipscat/pixel_math/README.md for more info.

Parameters:
  • dk (int) – the change in k that we wish to find the margins for.

  • pix (int) – the healpix pixel to find margin pixels of.

  • edge (int) –

    the edge we want to find for the given pixel. (0-7)

    • 0: NE edge

    • 1: E corner

    • 2: SE edge

    • 3: S corner

    • 4: SW edge

    • 5: W corner

    • 6: NW edge

    • 7: N corner

Returns:

one-dimensional numpy array of long integers, filled with the healpix pixels at order kk+dk that border the given edge of pix.

get_margin(order, pix, d_order)[source]#

Get all the pixels at order order+dk bordering pixel pix. See hipscat/pixel_math/README.md for more info.

Parameters:
  • order (int) – the healpix order of pix.

  • pix (int) – the healpix pixel to find margin pixels of.

  • d_order (int) – the change in k that we wish to find the margins for. Must be greater than kk.

Returns:

one-dimensional numpy array of long integers, filled with the healpix pixels at order kk+dk that border pix.

pixel_is_polar(order, pix)[source]#

Checks if a healpixel is a polar pixel.

Parameters:
  • order (int) – the healpix order of the pixel to be checked.

  • pix (int) – the id of a healpixel to be checked. must be in the nested id scheme.

Returns:

tuple of a boolean and a string, the boolean indicating whether the pixel polar and the string denoting which pole it is (‘North’ or ‘South’)). string is empty in the case that the pixel isn’t polar.