hpmoc.partial module

A partial HEALPix skymap class supporting multi-resolution HEALPix skymaps in NUNIQ ordering.

class hpmoc.partial.PartialUniqSkymap(s, u, copy=False, name=None, point_sources=None, meta=None, empty=None, compress=False, interp='nearest')

Bases: hpmoc.abstract.AbstractPartialUniqSkymap

A HEALPix skymap object residing in memory with NUNIQ ordering. Only a subset of the full sky. You can index into a PartialUniqSkymap with NUNIQ indices to get a skymap with the same shape (optionally padding missing values with a second index argument). You can also use index notation to set pixel values at the specified NUNIQ index locations.

ang_dist(ra, dec, degrees=True)

Get distances from each pixel in this skymap to the point at right-ascension ra and declination dec.

Parameters
  • ra (array-like or astropy.units.Quantity) – Right-ascension of the point

  • dec (array-like or astropy.units.Quantity) – Declination of the point

  • degrees (bool, optional) – If ra and dec are astropy.units.Quantity instances, they will be automatically converted. If they are unitless scalars, they will be interpreted as degrees if degrees=True, radians otherwise.

Returns

ang_dist – The distances of each pixel in this skymap to the point at ra, dec in degrees.

Return type

astropy.units.Quantity

Examples

We should find that the distance from any pixel to the North pole is equal to 90 minus the declination (within some small error):

>>> import numpy as np
>>> from astropy.units import deg
>>> skymap = PartialUniqSkymap(*([4+np.arange(12)]*2))
>>> _, dec = skymap.coords()
>>> dec
<Quantity [ 41.8103149,  41.8103149,  41.8103149,  41.8103149,   0.       ,
             0.       ,   0.       ,   0.       , -41.8103149, -41.8103149,
           -41.8103149, -41.8103149] deg>
>>> Δθ⃗ = skymap.ang_dist(32, 90)
>>> Δθ⃗
<Quantity [0.84106867, 0.84106867, 0.84106867, 0.84106867, 1.57079633,
           1.57079633, 1.57079633, 1.57079633, 2.30052398, 2.30052398,
           2.30052398, 2.30052398] rad>
>>> np.all(abs(Δθ⃗+dec-90*deg).to('deg').value<1e-13)
True

Likewise, the distance from any pixel to the South pole should be equal to 90 plus the declination:

>>> not np.around(skymap.ang_dist(359, -90)-dec-90*deg,
...               14).value.any()
True
apply(func: Callable, copy: bool = True, quantity: bool = True) __class__

Map a function to this skymap’s pixels and return a new skymap whose values are the returned values of func. Note that func(self.s) must therefore return values corresponding to the same pixels as the input skymap. Use this to, for example, get the logarithm of this skymap as self.apply(np.log). If copy == True, make a copy of the NUNIQ indices self.u; otherwise, share them with the new skymap. To operate directly on the skymap value (even if it is stored as an astropy.units.Quantity), use quantity = True. Note that this will strip units from the result.

area()

Area per-pixel for pixels in this skymap in astropy.unit.sr.

astype(dtype, copy=True, **kwargs)

Return a new PartialUniqSkymap with the data-type of s set to dtype. If copy=True, always make sure both u and s are copies of the original data in the new array. Otherwise, re-use u and (if possible given the provided dtype and **kwargs) s. copy and **kwargs are passed on to s.astype to make the conversion.

azeqview(*scatter, **kwargs)
cartview(*scatter, **kwargs)
compress(stype=None, utype=None, **kwargs)

Eliminate redundant pixels with utils.uniq_minimize and store indices u in the smallest integer size that represents all values.

Parameters
  • stype (type, optional) – If provided, store s as this type. Defaults to s.dtype.

  • utype (type, optional) – If provided, store u as this type. Defaults to the smallest int type required to store all values of u.

  • kwargs – Keyword arguments to pass to hpmoc.utils.uniq_minimize.

Returns

compressed – A compressed version of this skymap.

Return type

PartialUniqSkymap

coords()

Get the sky coordinates (right-ascension and declination, ICRS) corresponding to each pixel in the skymap.

Returns

ra_dec – 2D array whose first row is the right-ascension and second row is the declination (in degrees) of each pixel. You can get each of these individually with ra, dec = self.coords(). self.coords()[:, i] corresponds to RA, Dec for self.s[i].

Return type

astropy.units.Quantity

copy()

Return a copy of this instance with separate copies of its data. The copy can be edited without affecting the original.

fill(nside=None, pad=None, as_skymap=False)

Return a full-sky nested HEALPix skymap at NSIDE resolution nside.

Parameters
  • nside (int) – HEALPix NSIDE value of the output map. If not provided, use the highest NSIDE value in this skymap’s nside values to preserve detail.

  • pad (float, optional) – Fill in missing values with pad (if not provided, use healpy.UNSEEN). Preserves astropy.units.Unit of this skymap’s pixel values (if s is an astropy.units.Quantity).

  • as_skymap (bool, optional) – If True, return a PartialUniqSkymap instance with the new pixelization (instead of a bare array with implicit indexing).

Returns

s – The filled-in skymap, either as an array if as_skymap == False or as a new PartialUniqSkymap instance.

Return type

array or PartialUniqSkymap

fixed(nside=None)

Re-raster to a fixed NSIDE. Like fill but for partial skymaps.

Parameters

nside (int) – HEALPix NSIDE value of the output map. If not provided, use the highest NSIDE value in this skymap’s nside values to preserve detail.

gnomview(*scatter, **kwargs)
gridplot(*skymaps: Union[hpmoc.PartialUniqSkymap, nptyping.types._ndarray.NDArray[typing.Any, Any], Tuple[nptyping.types._ndarray.NDArray[typing.Any, Any], Optional[Union[nptyping.types._ndarray.NDArray[typing.Any, int], astropy.wcs.WCS, str]]]], **kwargs) Tuple[matplotlib.gridspec.GridSpec, List[List[astropy.visualization.wcsaxes.WCSAxes]]]

Plot this skymap and any others in skymaps. A thin wrapper around hpmoc.plot.gridplot.

Parameters
  • *skymaps (PartialUniqSkymap or map-like) – Skymaps to pass to gridplot. Can be anything accepted by hpmoc.plot.plot, which hpmoc.plot.gridplot will use to display them.

  • **kwargs – Keyword arguments for hpmoc.plot.gridplot.

intersection(u)

See utils.uniq_intersection.

max()

Maximum skymap value = self.s.max().

min()

Minimum skymap value = self.s.min().

mollview(*scatter, **kwargs)
multiplot(*, nest: bool = True, **kwargs)

Call plotters.multiplot with the default transform suitable for a PartialUniqSkymap.

Parameters
  • *skymapsₗ (List[Union[PartialUniqSkymap, array]]) – Skymaps to plot. Can be PartialUniqSkymap instances or full-sky single-resolution skymaps.

  • **kwargs – Keyword arguments passed to plotters.multiplot.

Returns

fig – A new matplotlib figure containing the specified subplots.

Return type

matplotlib.figure.Figure

See also

plot.multiplot, PartialUniqSkymap.plot, plot.plot

nside(as_skymap=False, copy=False, **kwargs)

Pixel NSIDE values. If as_skymap=True, return as a PartialUniqSkymap instance (with **kwargs passed to init).

orders(as_skymap=False, copy=True, **kwargs)

HEALPix order values. If as_skymap=True, return as a PartialUniqSkymap instance (with **kwargs passed to init).

orthview(*scatter, **kwargs)
plot(*scatter: hpmoc.points.PointsTuple, **kwargs) Union[astropy.visualization.wcsaxes.WCSAxes, astropy.visualization.wcsaxes.WCSAxesSubplot]

Plot this skymap. A thin wrapper around hpmoc.plot.plot that automatically includes scatter plots for this instance’s point_sources.

Parameters
  • scatter (PointsTuple) – Additional point sources to plot (in addition to those stored in self.point_sources).

  • kwargs – Keywork arguments to be passed to hpmoc.plot.plot.

Returns

ax – The axes that were just plotted to.

Return type

WCSAxes or WCSAxesSubplot

See also

hpmoc.plot.plot

point_sources: List[hpmoc.points.PointsTuple]
quantiles(quantiles: Iterable[float]) Tuple[Iterator, nptyping.types._ndarray.NDArray, float]

Get an iterator of downselected skymaps partitioned by quantiles. For example, get the smallest sky area containing 90% of the probability (or whatever other intensive quantity this skymap represents) with quantiles=[0.1, 1].

Parameters

quantiles (Iterable[float]) – Quantiles from which to select pixels. Must be in ascending order with values in the interval [0, 1]. These will form endpoints for partitions of the skymap. For example, [0.1, 0.9] will omit the lowest and highest value pixels, giving the intermediate pixels accounting for 80% of the integrated skymap. Note that quantiles returned by this function are non-intersecting and half-open on the right (as with python indices), with the exception of 1 for the last pixel; for example, [0, 1] will include all pixels, [0.5, 1] will include the highest density pixels accounting for 50% of the integrated skymap value, [0, 0.5, 1] will partition the skymap into non-intersecting sets of pixels accounting for the high- and low-density partitions of the skymap by integrated value, etc.

Returns

  • partitions (Iterator[PartialUniqSkymap]) – A generator containing non-overlapping skymaps falling into the specified quantiles. Will have length one less than the number of quantiles, which form the boundaries of the partitions. Always an iterable, even if only two quantiles are provided; you can unpack a single value with, e.g., x, = m.quantiles(...).

  • levels (astropy.units.Quantity or array) – Values of the skymap at each quantile. Useful for e.g. contour plots (though PartialUniqSkymap.plot will handle this automatically).

  • norm (astropy.units.Quantity or float) – The integrated value of the skymap. Useful for calculating the integral of each quantile.

Raises

ValueError – If quantiles has length less than 2; if its values are not in order and contained in the interval [0, 1]; if nside and skymap cannot be broadcast together; if any values in skymap are negative; or if the total integrated skymap equals zero, in which case quantiles are undefined.

read(*, strategy: str = 'basic', **kwargs)

Read a PartialUniqSkymap from file using astropy.table.Table.read. See IoRegistry attributes for details. When called as a bound method or passed a PartialUniqSkymap as the first argument, uses that skymap as a mask for the bound skymap and only loads those pixels.

Parameters
  • mask (PartialUniqSkymap, optional) – Only read pixels overlapping with this mask.

  • file (file or str) – The file object or filename to read from.

  • *args – Arguments to pass on to IoRegistry.

  • strategy (str, optional) – Name of the IoRegistry strategy to use for reads/writes.

  • **kwargs – Keyword arguments to pass on to a specific IoStrategy. See IoRegistry properties for details.

Returns

m – A new PartialUniqSkymap instance with the specified data.

Return type

PartialUniqSkymap

Examples

For multi-skymap fits files, you can load the full set of skymaps with hdul = astropy.io.fits.open(fitsfile) and then load the skymap of interest with PartialUniqSkymap.read(hdul[i]). With fits files, you can memory-map the resulting skymap to save memory (at the expense of speed) by adding memmap=True as a keyword argument.

See also

astropy.table.Table.read, astropy.table.Table.read.help, astropy.io.fits.open, IoRegistry, IoRegistry.basic

render(u⃗o, pad=None, valid=None, mask_missing=False)

Like reraster, but u⃗ᵒ does not need to be unique. Use this to e.g. render a skymap to a plot. Unlike reraster, will not return a PartialUniqSkymap; instead, simply returns the pixel values corresponding to u⃗ᵒ.

u⃗ᵒ can also be an astropy.wcs.WCS world coordinate system, in which case the returned array will contain the pixel values of this skymap in that coordinate system (with regions outside of the projection set to np.nan).

Parameters
  • u⃗ᵒ (array or astropy.wcs.WCS) – The pixels to fill. If an array, treated as UNIQ indices (duplicates allowed); if WCS, treated as a set of pixels to render to.

  • pad (float, optional) – A pad value to use for pixels not contained in the maps. Defaults to None, which will raise an error if suitable values cannot be found for every valid pixel in u⃗ᵒ (this does not apply to values outside a WCS projection, which will take on np.nan values).

  • valid (array, optional) – If provided, results will be scattered into an array of the same shape as valid, filling the indices where valid==True. The number of True values in valid must therefore equal the length of u⃗ᵒ. This argument only makes sense if u⃗ᵒ is an array of NUNIQ indices; if it is a WCS instance and valid is provided, an error is raised. Use valid to produce plots or to reuse indices produced by wcs2mask_and_uniq in several render invocations. See note on how mask_missing affects the result.

  • mask_missing (bool) – If mask_missing=True, return a np.ma.core.MaskedArray. Missing values are tolerated and are marked as True in the mask_missing. They will be set to pad or 0 in the data field. If valid is also provided, then the output will still be a np.ma.core.MaskedArray, but will be set to True wherever valid == False in addition to wherever pixels are missing (and will still take on masked values of np.nan in the invalid regions).

Returns

pixels – Pixel values at locations specified in u⃗ᵒ. If u⃗ᵒ is a WCS instance, then values outside of the projection will be set to np.nan.

Return type

array

reraster(u⃗o, pad=None, mask_missing=False, copy=True)

Return a new PartialUniqSkymap instance with the same pixel values rerasterized to match the output NUNIQ indices u⃗ᵒ. Fill in missing values in the output skymap with pad. If pad is not provided and this skymap does not cover the full region defined in u⃗ᵒ, raises a ValueError. Preserves astropy.units.Unit of this skymap’s pixel values (if s is an astropy.units.Quantity). If copy is False, use u⃗ᵒ as the indices of the new skymap; otherwise, use a copy.

s: Union[nptyping.types._ndarray.NDArray[typing.Any, Any], astropy.units.Quantity]
sort(copy=True)

Sort this skymap by UNIQ indices u (sorting s as well, of course). If copy=True, copy u and s and return a new PartialUniqSkymap; otherwise, sort them in-place and return this PartialUniqSkymap.

to(*args, **kwargs)

Convert the units of this skymap’s pixels (if they are stored as an astropy.units.Quantity instance).

Parameters
  • *args – Passed on to astropy.units.Quantity.to.

  • **kwargs – Passed on to astropy.units.Quantity.to.

Raises

TypeError – If self.s is not a Quantity.

to_table(name=None, uname='UNIQ')

Return a new astropy.table.Table whose UNIQ column is the NUNIQ indices self.u and PIXELS (or self.name, if set) column is the skymap pixel values s. Optionally override the pixel value column name and/or the NUNIQ column name with the name and uname arguments respectively.

u: nptyping.types._ndarray.NDArray[typing.Any, int]
property unit

self.s.unit, if defined; otherwise None.

unzip_atlas()

Return 12 sub-skymaps corresponding to the HEALPix base pixels.

unzip_orders()

Return a list of PartialUniqSkymap instances corresponding to the parts of this sky imaged and each HEALPix order. Length equals the maximum order of this skymap. Empty terms indicate that this skymap does not have pixels of the corresponding HEALPix order.

property value

Get a dimensionless view of this skymap (no effect if s is not an astropy.units.Quantity).

write(file: Union[IO, str], *args, strategy: str = 'basic', **kwargs)

Write a PartialUniqSkymap instance to file using the specified IoStrategy. See IoRegistry attributes for details.

Parameters
  • file (file or str) – The file object or filename to write to.

  • *args – Arguments to pass on to IoRegistry.

  • strategy (str, optional) – Name of the IoRegistry strategy to use for reads/writes.

  • **kwargs – Keyword arguments to pass on to a specific IoStrategy. See IoRegistry properties for details.

See also

astropy.table.Table.write, astropy.table.Table.write.help, astropy.io.fits.open, IoRegistry, IoRegistry.basic

hpmoc.partial.bool_to_uint8(s)

Convert a boolean value to np.uint8.

hpmoc.partial.diadic_dunder(pad=None, coarse=False, post=None)

Implement diadic dunder methods like __add__, __radd__, etc. for scalar, array, and PartialUniqSkymap arguments using uniq_diadic. pad and coarse are passed to uniq_diadic when performing operations between PartialUniqSkymap instances. If provided, run post on the result array before creating a new PartialUniqSkymap instance out of it (for example, to cast booleans to integers).