mocpy.MOC

class mocpy.MOC(store_index)[source]

Multi-order spatial coverage class.

A MOC describes the coverage of an arbitrary region on the unit sphere. MOCs are usually used for describing the global coverage of catalog/image surveys such as GALEX or SDSS. A MOC corresponds to a list of HEALPix cells at different depths. This class gives you the possibility to:

  1. Define MOC objects:

  • From a FITS file that stores HEALPix cells (see load(path, 'fits')).

  • Directly from a list of HEALPix cells expressed either as a numpy structural array (see from_healpix_cells) or a simple python dictionary (see from_json).

  • From a list of sky coordinates (see from_skycoords, from_lonlat).

  • From a convex/concave polygon (see from_polygon).

  • From a cone (will be implemented in a next version).

  1. Perform fast logical operations between MOC objects:

  1. Plot the MOC objects:

  • Draw the MOC with its HEALPix cells (see fill)

  • Draw the perimeter of a MOC (see border)

  1. Get the sky coordinates defining the border(s) of MOC objects (see get_boundaries).

  2. Serialize MOC objects to astropy.io.fits.HDUList or JSON dictionary and save it to a file.

Is a Spatial Coverage (S-MOC).

Args:

create_key: Object ensure __init__ is called by super-class/class-methods only store_index: index of the S-MOC in the rust-side storage

add_neighbours()

Extend the MOC instance so that it includes the HEALPix cells touching its border.

The depth of the HEALPix cells added at the border is equal to the maximum depth of the MOC instance.

Returns:
mocMOC

self extended by one degree of neighbours.

barycenter()[source]

Return the Barycenter of the MOC.

border(ax, wcs, **kw_mpl_pathpatch)[source]

Draw the MOC border.s on a matplotlib axis.

This performs the projection of the sky coordinates defining the perimeter of the MOC to the pixel image coordinate system. You are able to specify various styling kwargs for matplotlib.patches.PathPatch (see the list of valid keywords).

Parameters:
axmatplotlib.axes.Axes

Matplotlib axis.

wcsastropy.wcs.WCS

WCS defining the World system <-> Image system projection.

kw_mpl_pathpatch

Plotting arguments for matplotlib.patches.PathPatch

Examples

>>> from mocpy import MOC
>>> from astropy.coordinates import Latitude, Longitude
>>> import astropy.units as u
>>> import matplotlib.pyplot as plt
>>> # Create a MOC
>>> lon = Longitude([5, -5, -5, 5], u.deg)
>>> lat = Latitude([5, 5, -5, -5], u.deg)
>>> moc = MOC.from_polygon(lon, lat)
>>> # Plot the MOC using matplotlib
>>> fig = plt.figure(figsize=(10, 10))
>>> wcs = moc.wcs(fig)
>>> ax = fig.add_subplot(projection=wcs)
>>> moc.border(ax, wcs, color='blue')
complement()

Return the complement of the MOC instance.

Returns:
resultMOC or TimeMOC

The resulting MOC.

contains(lon, lat, keep_inside=True)[source]

Test wether a MOC contains –or not– the given points. Returns a boolean mask array.

Deprecated since version 0.11.1: contains is replaced by contains_lonlat for naming consistency. Please consider switching.

Parameters:
lonastropy.coordinates.Longitude or its supertype astropy.units.Quantity

Right ascension array in deg

latastropy.coordinates.Latitude or its supertype astropy.units.Quantity

Declination array in deg

keep_insidebool, optional

True by default. If so the mask describes coordinates lying inside the MOC. If keep_inside is false, contains will return the mask of the coordinates lying outside the MOC.

Returns:
arrayndarray

A mask boolean array

contains_lonlat(lon, lat, keep_inside=True)[source]

Test wether a MOC contains (or not) the given points. Returns a boolean mask array.

Parameters:
lonastropy.coordinates.Longitude or its supertype astropy.units.Quantity

Right ascension array in deg

latastropy.coordinates.Latitude or its supertype astropy.units.Quantity

Declination array in deg

keep_insidebool, optional

True by default. If so the mask describes coordinates lying inside the MOC. If keep_inside is false, contains will return the mask of the coordinates lying outside the MOC.

Returns:
arrayndarray

A mask boolean array

Raises:
ValueErrorIf lon and lat have mismatched shapes.

Examples

>>> from mocpy import MOC
>>> import numpy as np
>>> from astropy.coordinates import Angle
>>> import astropy.units as u
>>> # create lists of coordinates
>>> lon = Angle(np.array([[1, 2, 3], [-2, -40, -5]]), unit=u.deg)
>>> lat = Angle(np.array([[20, 25, 10], [-60, 80, 0]]), unit=u.deg)
>>> # create a polygonal moc from these
>>> moc = MOC.from_polygon(lon=lon, lat=lat, max_depth=12)
>>> moc.contains_lonlat(lon=lon, lat=lat) # returns all true
array([[ True,  True,  True],
       [ True,  True,  True]])
contains_skycoords(skycoords, keep_inside=True)[source]

Return a boolean mask array of the positions lying inside (or outside) the MOC instance.

Parameters:
skycoordsastropy.coordinates.SkyCoord

The sky coordinates that will be tested.

keep_insidebool, optional

True by default. If so the mask describes coordinates lying inside the MOC. If keep_inside is false, contains will return the mask of the coordinates lying outside the MOC.

Returns:
arrayndarray

A mask boolean array

See also

contains_lonlat
contracted()

Return the MOC contracted by removing the internal border made of cells at the MOC maximum depth.

The only difference with respect to remove_neighbours is that contracted returns a new MOC instead of modifying the existing one.

Returns:
mocMOC

The extended MOC

degrade_to_order(new_order)[source]

Degrade the MOC instance to a new, less precise, MOC.

The maximum depth (i.e. the depth of the smallest HEALPix cells that can be found in the MOC) of the degraded MOC is set to new_order.

Parameters:
new_orderint

Maximum depth of the output degraded MOC.

Returns:
mocMOC

The degraded MOC.

difference(another_moc, *args)

Difference between the MOC instance and other MOCs.

Parameters:
another_mocMOC

The MOC used that will be substracted to self.

argsMOC

Other additional MOCs to perform the difference with.

Returns:
resultMOC or TimeMOC

The resulting MOC.

display_preview(y_size=300)[source]

Display a preview of the MOC (calling internally the to_rgba method).

Parameters:
y_sizethe number of pixels along the y-axis, default value is 300
empty()

(e.g. a numpy boolean array).

Returns:
result: bool

True if the MOC instance is empty.

extended()

Return the MOC extended by the external border made of cells at the MOC maximum depth.

The only difference with respect to add_neighbours is that extended returns a new MOC instead of modifying the existing one.

Returns:
mocMOC

The extended MOC

fill(ax, wcs, **kw_mpl_pathpatch)[source]

Draw the MOC on a matplotlib axis.

This performs the projection of the cells from the world coordinate system to the pixel image coordinate system. You are able to specify various styling kwargs for matplotlib.patches.PathPatch (see the list of valid keywords).

Parameters:
axmatplotlib.axes.Axes

Matplotlib axis.

wcsastropy.wcs.WCS

WCS defining the World system <-> Image system projection.

kw_mpl_pathpatch

Plotting arguments for matplotlib.patches.PathPatch.

Examples

>>> from mocpy import MOC
>>> import astropy.units as u
>>> import matplotlib.pyplot as plt
>>> # Create a MOC
>>> moc = MOC.from_ring(external_radius=2*u.deg,
...                     internal_radius=1*u.deg,
...                     lat=0*u.rad, lon=0*u.rad,
...                     max_depth=13,
...                    )
>>> # Plot the MOC using matplotlib
>>> fig = plt.figure(figsize=(10, 10))
>>> wcs = moc.wcs(fig)
>>> ax = fig.add_subplot(projection=wcs)
>>> moc.fill(ax, wcs, color='blue')
flatten()

Return the list of indices of all cells in the MOC at the MOC depth.

classmethod from_cone(lon, lat, radius, max_depth, delta_depth=2)[source]

Create a MOC from a cone.

The cone is centered around the (lon, lat) position with a radius expressed by radius.

Parameters:
lonastropy.coordinates.Longitude or its supertype astropy.units.Quantity

The longitude of the center of the cone.

latastropy.coordinates.Latitude or its supertype astropy.units.Quantity

The latitude of the center of the cone.

radiusastropy.coordinates.Angle

The radius angle of the cone.

max_depthint

Maximum HEALPix cell resolution.

delta_depthint, optional

To control the approximation, you can choose to perform the computations at a deeper depth using the depth_delta parameter. The depth at which the computations will be made will therefore be equal to max_depth + depth_delta.

Returns:
resultMOC

The resulting MOC

Examples

>>> from mocpy import MOC
>>> import astropy.units as u
>>> from astropy.coordinates import Angle, Longitude, Latitude
>>> moc = MOC.from_cone(
...  lon=Longitude(0 * u.deg),
...  lat=Latitude(0 * u.deg),
...  radius=Angle(10, u.deg),
...  max_depth=10
... )
classmethod from_depth29_ranges(max_depth, ranges)[source]

Create a MOC from a set of ranges of HEALPix Nested indices at order 29.

For each range, the lower bound is inclusive and the upper bound is exclusive. The range [0, 12*4^29[ represents the full-sky.

Parameters:
max_depthint, The resolution of the MOC
rangesndarray, optional

a N x 2 numpy array representing the set of depth 29 HEALPix nested ranges. defaults to np.zeros((0, 2), dtype=np.uint64)

Returns:
mocMOC

The MOC

classmethod from_elliptical_cone(lon, lat, a, b, pa, max_depth, delta_depth=2)[source]

Create a MOC from an elliptical cone.

The ellipse is centered around the (lon, lat) position. a (resp. b) corresponds to the semi-major axis magnitude (resp. semi-minor axis magnitude). pa is expressed as a Angle and defines the position angle of the elliptical cone.

Parameters:
lonastropy.coordinates.Longitude or its supertype astropy.units.Quantity

The longitude of the center of the elliptical cone.

latastropy.coordinates.Latitude or its supertype astropy.units.Quantity

The latitude of the center of the elliptical cone.

aastropy.coordinates.Angle

The semi-major axis angle of the elliptical cone.

bastropy.coordinates.Angle

The semi-minor axis angle of the elliptical cone.

paastropy.coordinates.Angle

The position angle (i.e. the angle between the north and the semi-major axis, east-of-north).

max_depthint

Maximum HEALPix cell resolution.

delta_depthint, optional

To control the approximation, you can choose to perform the computations at a deeper depth using the depth_delta parameter. The depth at which the computations will be made will therefore be equal to depth + depth_delta.

Returns:
resultMOC

The resulting MOC

Examples

>>> from mocpy import MOC
>>> import astropy.units as u
>>> from astropy.coordinates import Angle, Longitude, Latitude
>>> moc = MOC.from_elliptical_cone(
...  lon=Longitude(0 * u.deg),
...  lat=Latitude(0 * u.deg),
...  a=Angle(10, u.deg),
...  b=Angle(5, u.deg),
...  pa=Angle(0, u.deg),
...  max_depth=10
... )
classmethod from_fits(path_or_url, timeout=1000)

Load a MOC from a FITS file.

The specified FITS file must store the MOC (i.e. the list of HEALPix cells it contains) in a binary HDU table.

Parameters:
pathstr

The path to the FITS file.

timeoutfloat

Timeout for the query, defaults to 1000s

Returns:
resultMOC or TimeMOC

The resulting MOC.

classmethod from_fits_image(hdu, max_norder, mask=None)[source]

Create a MOC from an image stored as a FITS file.

Parameters:
hduHDU object

HDU containing the data of the image

max_norderint

The moc resolution.

masknumpy.ndarray, optional

A boolean array of the same size of the image where pixels having the value 1 are part of the final MOC and pixels having the value 0 are not.

Returns:
mocMOC

The resulting MOC.

classmethod from_fits_images(path_l, max_norder, hdu_index=0)[source]

Load a MOC from a set of FITS file images.

Assumes the data of the image is stored in the first HDU of the FITS file. Please call from_fits_image for passing another hdu than the first one.

Parameters:
path_l[str]

A list of path where the fits image are located.

max_norderint

The MOC resolution.

hdu_indexint

Index of the the HDU containing the image in each FITS file (default = 0)

Returns:
mocMOC

The union of all the MOCs created from the paths found in path_l.

classmethod from_healpix_cells(ipix, depth, max_depth)[source]

Create a MOC from a set of HEALPix cells at various depths.

Parameters:
ipixnumpy.ndarray

HEALPix cell indices in the NESTED notation. dtype must be np.uint64

depthnumpy.ndarray

Depth of the HEALPix cells. Must be of the same size of ipix. dtype must be np.uint8. Corresponds to the level of an HEALPix cell in astropy.healpix.

max_depthint, The resolution of the MOC (degrades on the fly input cells if necessary)
Returns:
mocMOC

The MOC

Raises:
IndexError

When ipix and depth do not have the same shape

classmethod from_ivorn(ivorn, nside=256)[source]

Create a MOC object from a given ivorn.

Parameters:
ivornstr
nsideint, optional

256 by default

Returns:
resultMOC

The resulting MOC.

classmethod from_json(json_moc)

Create a MOC from a dictionary of HEALPix cell arrays indexed by their depth.

Parameters:
json_mocdict(str[int]

A dictionary of HEALPix cell arrays indexed by their depth.

Returns:
mocMOC or TimeMOC

the MOC.

classmethod from_lonlat(lon, lat, max_norder)[source]

Create a MOC from astropy lon, lat astropy.units.Quantity.

Parameters:
lonastropy.coordinates.Longitude or its supertype astropy.units.Quantity

The longitudes of the sky coordinates belonging to the MOC.

latastropy.coordinates.Latitude or its supertype astropy.units.Quantity

The latitudes of the sky coordinates belonging to the MOC.

max_norderint

The depth of the smallest HEALPix cells contained in the MOC.

Returns:
resultMOC

The resulting MOC

classmethod from_multiordermap_fits_file(path, cumul_from=0.0, cumul_to=1.0, asc=False, strict=True, no_split=True, reverse_decent=False)[source]

Create a MOC from a mutli-order map FITS file.

HEALPix cells are first sorted by their values. The MOC contains the cells from which the cumulative value is between cumul_from and cumul_to. Cells being on the fence are recursively splitted and added until the depth of the cells is equal to max_norder.

For compatibility with Aladin, use no_split=False and reverse_decent=True

Remark: using no_split=False, the way the cells overlapping with the low and high thresholds are split is somewhat arbitrary.

Parameters:
pathstr or pathlib.Path

The path to the file to save the MOC in.

cumul_fromfloat

Cumulative value from which cells will be added to the MOC

cumul_tofloat

Cumulative value to which cells will be added to the MOC

asc: boolean

the cumulative value is computed from lower to highest densities instead of from highest to lowest

strict: boolean

(sub-)cells overlapping the cumul_from or cumul_to values are not added

no_split: boolean

cells overlapping the cumul_from or cumul_to values are not recursively split

reverse_decent: boolean

perform the recursive decent from the highest cell number to the lowest (to be compatible with Aladin)

Returns:
resultMOC

The resulting MOC

classmethod from_polygon(lon, lat, complement=False, max_depth=10)[source]

Create a MOC from a polygon.

The polygon is given as lon and lat astropy.units.Quantity that define the vertices of the polygon. Concave, convex and self-intersecting polygons are accepted.

Parameters:
lonastropy.coordinates.Longitude or its supertype astropy.units.Quantity

The longitudes defining the polygon. Can describe convex and concave polygons but not self-intersecting ones.

latastropy.coordinates.Latitude or its supertype astropy.units.Quantity

The latitudes defining the polygon. Can describe convex and concave polygons but not self-intersecting ones.

complementreturn the complement of the polygon. Set to False by default.
max_depthint, optional

The resolution of the MOC. Set to 10 by default.

Returns:
resultMOC

The resulting MOC

classmethod from_polygon_skycoord(skycoord, max_depth=10)[source]

Create a MOC from a polygon.

The polygon is given as an astropy.coordinates.SkyCoord that contains the vertices of the polygon. Concave, convex and self-intersecting polygons are accepted.

Parameters:
skycoordastropy.coordinates.SkyCoord

The sky coordinates defining the vertices of a polygon. It can describe a convex or concave polygon but not a self-intersecting one.

max_depthint, optional

The resolution of the MOC. Set to 10 by default.

Returns:
resultMOC

The resulting MOC

classmethod from_ring(lon, lat, internal_radius, external_radius, max_depth, delta_depth=2)[source]

Create a MOC from a ring.

The cone is centered around the (lon, lat) position with an internal radius expressed by internal_radius and an external radius expressed by external_radius.

Parameters:
lonastropy.coordinates.Longitude or its supertype astropy.units.Quantity

The longitude of the center of the ring.

latastropy.coordinates.Latitude or its supertype astropy.units.Quantity

The latitude of the center of the ring.

internal_radiusastropy.coordinates.Angle

The internal radius angle of the ring.

external_radiusastropy.coordinates.Angle

The external radius angle of the ring.

max_depthint

Maximum HEALPix cell resolution.

delta_depthint, optional

To control the approximation, you can choose to perform the computations at a deeper depth using the depth_delta parameter. The depth at which the computations will be made will therefore be equal to max_depth + depth_delta.

Returns:
resultMOC

The resulting MOC

Examples

>>> from mocpy import MOC
>>> import astropy.units as u
>>> from astropy.coordinates import Angle, Longitude, Latitude
>>> moc = MOC.from_ring(
...  lon=Longitude(0 * u.deg),
...  lat=Latitude(0 * u.deg),
...  internal_radius=Angle(5, u.deg),
...  external_radius=Angle(10, u.deg),
...  max_depth=10
... )
classmethod from_skycoords(skycoords, max_norder)[source]

Create a MOC from an astropy.coordinates.SkyCoord.

Parameters:
skycoordsastropy.coordinates.SkyCoord

The sky coordinates that will belong to the MOC.

max_norderint

The depth of the smallest HEALPix cells contained in the MOC.

Returns:
resultMOC

The resulting MOC

classmethod from_stmoc_time_fold(tmoc, stmoc)[source]

Build a new S-MOC from the fold operation of the given ST-MOC by the given T-MOC.

Parameters:
tmocTimeMoc
stmocSTMoc
classmethod from_str(value)

Create a MOC from a str.

This grammar is expressed is the MOC IVOA specification at section 2.3.2.

Parameters:
valuestr

The MOC as a string following the grammar rules.

Returns:
mocMOC or TimeMOC

The resulting MOC

See also

from_string

a duplicate of this method, with added fold option

Examples

>>> from mocpy import MOC
>>> moc = MOC.from_str("2/2-25 28 29 4/0 6/")
classmethod from_string(value, format='ascii')[source]

Deserialize the Spatial MOC from the given string.

Format can be ‘ascii’ or ‘json’, though the json format is not officially supported by the IVOA.

WARNING: the serialization must be strict, i.e. must not contain overlapping elements

Parameters:
formatstr, optional

The format in which the MOC will be serialized before being saved. Possible formats are “ascii” or “json”. By default, format is set to “ascii”.

classmethod from_url(url)[source]

Create a MOC object from a given url.

Parameters:
urlstr

The url of a FITS file storing a MOC.

Returns:
resultMOC

The resulting MOC.

classmethod from_valued_healpix_cells(uniq, values, max_depth=None, values_are_densities=False, cumul_from=0.0, cumul_to=1.0, asc=False, strict=True, no_split=True, reverse_decent=False)[source]

Create a MOC from a list of uniq associated with values.

HEALPix cells are first sorted by their values. The MOC contains the cells from which the cumulative value is between cumul_from and cumul_to. Cells being on the fence are recursively splitted and added until the depth of the cells is equal to max_norder.

For compatibility with Aladin, use no_split=False and reverse_decent=True

Remark: using no_split=False, the way the cells overlapping with the low and high thresholds are split is somewhat arbitrary.

Parameters:
uniqnumpy.ndarray

HEALPix cell indices written in uniq. dtype must be np.uint64

valuesnumpy.ndarray

Value associated with each uniq cells. dtype must be np.float64

max_depthint

The max depth of the MOC, should be at least as large as the depth corresponding of the smallest HEALPix cell found in uniq. Warnings: 1 - the depth of the returned MOC will be at least as deep as the smallest HEALPix cell found in uniq. 2 - contrary to MOCPy before v0.12, the user has to degrade the MOC if max_depth < smallest HEALPix cell depth.

values_are_densities: tell whether the values depend on the cell area or not
cumul_fromfloat

Cumulative value from which cells will be added to the MOC

cumul_tofloat

Cumulative value to which cells will be added to the MOC

asc: boolean

the cumulative value is computed from lower to highest densities instead of from highest to lowest

strict: boolean

(sub-)cells overlapping the cumul_from or cumul_to values are not added

no_split: boolean

cells overlapping the cumul_from or cumul_to values are not recursively split

reverse_decent: boolean

perform the recursive decent from the highest cell number to the lowest (to be compatible with Aladin)

Returns:
resultMOC

The resulting MOC

classmethod from_vizier_table(table_id, nside=256)[source]

Create a MOC object from a VizieR table.

Info: This method is already implemented in astroquery.cds. You can ask to get a mocpy.moc.MOC object from a vizier catalog ID.

Parameters:
table_idstr

table index

nsideint, optional

256 by default

Returns:
resultMOC

The resulting MOC.

get_boundaries(order=None)[source]

Return the sky coordinates defining the border(s) of the MOC.

The border(s) are expressed as a list of SkyCoord. Each SkyCoord refers to the coordinates of one border of the MOC (i.e. either a border of a connected MOC part or a border of a hole located in a connected MOC part). This function is currently not stable: encoding a vertex of a HEALPix cell (N, E, S, W) should not depend on the position of the vertex but rather on the uniq value (+ 2 bits to encode the direction of the vertex).

Parameters:
orderint

The depth of the MOC before computing its boundaries. A shallow depth leads to a faster computation. By default the maximum depth of the MOC is taken.

Returns:
coords: [SkyCoord]

A list of SkyCoord each describing one border.

Raises:
DeprecationWarning

This method is not stable and not tested! A future more stable algorithm will be implemented!

intersection(another_moc, *args)

Intersection between the MOC instance and other MOCs.

Parameters:
another_mocMOC

The MOC used for performing the intersection with self.

argsMOC

Other additional MOCs to perform the intersection with.

Returns:
resultMOC/TimeMOC

The resulting MOC.

largest_distance_from_coo_to_vertices(coo)[source]

Return the largest distance between the given coordinates and vertices of the MOC cells.

classmethod load(path, format='fits')[source]

Load the Spatial MOC from a file.

Format can be ‘fits’, ‘ascii’, or ‘json’, though the json format is not officially supported by the IVOA.

Parameters:
pathstr or pathlib.Path

The path to the file to load the MOC from.

formatstr, optional

The format from which the MOC is loaded. Possible formats are “fits”, “ascii” or “json”. By default, format is set to “fits”.

property max_index

Return the largest index (at the deepest absolute resolution) the MOC contains.

property max_order

Depth/order of the S-MOC.

property min_index

Return the smallest index (at the deepest absolute resolution) the MOC contains.

classmethod n_cells(depth)[source]

Get the number of cells for a given depth.

Parameters:
depthint

The depth. It is comprised between 0 and MAX_ORDER

Returns:
int

The number of cells at the given order

Examples

>>> from mocpy import MOC
>>> MOC.n_cells(0)
12
classmethod new_empty(max_depth)[source]

Create a new empty MOC of given depth.

Parameters:
max_depthint, The resolution of the MOC
Returns:
mocMOC

The MOC

static order_to_spatial_resolution(order)[source]

Convert a depth to its equivalent spatial resolution.

Parameters:
orderint

Spatial depth.

Returns:
spatial_resolutionAngle

Spatial resolution.

plot(title='MOC', frame=None)[source]

Plot the MOC object using a mollweide projection.

Deprecated: New fill and border methods produce more reliable results and allow you to specify additional matplotlib style parameters.

Parameters:
titlestr

The title of the plot

frameastropy.coordinates.BaseCoordinateFrame, optional

Describes the coordinate system the plot will be (ICRS, Galactic are the only coordinate systems supported).

query_simbad(max_rows=10000, timeout=1000)[source]

Query a view of SIMBAD data for SIMBAD objects in the coverage of the MOC instance.

Parameters:
max_rowsint, optional

maximum number of row returned

timeoutfloat, optional

timeout before aborting the query, default to 1000s

query_vizier_table(table_id: str, max_rows=10000, timeout=1000)[source]

Query a VizieR table for sources in the coverage of the MOC instance.

Parameters:
table_idstr

corresponds to a VizieR table id

max_rowsint, optional

maximum number of row returned

timeoutfloat, optional

timeout before aborting the query, default to 1000s

remove_neighbours()

Remove from the MOC instance the HEALPix cells located at its border.

The depth of the HEALPix cells removed is equal to the maximum depth of the MOC instance.

Returns:
mocMOC

self minus its HEALPix cells located at its border.

save(path, format='fits', overwrite=False, pre_v2=False, fold=0, fits_keywords=None)

Write the MOC to a file.

Format can be ‘fits’, ‘ascii’, or ‘json’, though the json format is not officially supported by the IVOA.

Parameters:
pathstr or pathlib.Path

The path to the file to save the MOC in.

formatstr, optional

The format in which the MOC is saved. Possible formats are “fits”, “ascii” or “json”. By default, format is set to “fits”.

overwritebool, optional

If the file already exists and you want to overwrite it, then set the overwrite keyword. Default to False.

fold: int

if >0, print ascii or json output with a maximum line width

fits_keywords: dict, optional

Additional keywords to add to the FITS header.

serialize(format='fits', optional_kw_dict=None, pre_v2=False)

Serialize the MOC into a specific format.

Possible formats are FITS, JSON and STRING

Parameters:
formatstr

‘fits’ by default. The other possible choice is ‘json’ or ‘str’.

optional_kw_dictdict

Optional keywords arguments added to the FITS header. Only used if format equals to ‘fits’.

Returns:
resultastropy.io.fits.HDUList or JSON dictionary

The result of the serialization.

property sky_fraction

Sky fraction covered by the MOC.

Examples

>>> from mocpy import MOC
>>> MOC.from_string("0/0-11").sky_fraction # all sky
1.0
static spatial_resolution_to_order(spatial_resolution)[source]

Convert a spatial resolution to a MOC order.

Parameters:
spatial_resolutionAngle

Spatial resolution

Returns:
orderint

The order corresponding to the spatial resolution

split(include_indirect_neighbours=False)[source]

Return the disjoint MOCs this MOC contains.

Parameters:
include_indirect_neighboursbool

if false, only consider cells having a common edge as been part of a same MOC if true, also consider cells having a common vertex as been part of the same MOC

split_count(include_indirect_neighbours=False)[source]

Return the number of disjoint MOCs the given MOC contains.

Parameters:
include_indirect_neighboursbool

if false, only consider cells having a common edge as been part of a same MOC if true, also consider cells having a common vertex as been part of the same MOC

symmetric_difference(another_moc, *args)

Symmetric difference (XOR) between the MOC instance and other MOCs.

Parameters:
another_mocMOC

The MOC used that will be substracted to self.

argsMOC

Other additional MOCs to perform the difference with.

Returns:
resultMOC or TimeMOC

The resulting MOC.

property to_depth29_ranges

Return the list of order 29 HEALPix nested ranges this MOC contains.

to_rgba(y_size=300)[source]

Create a matplotlib compatible RGBA preview of the given MOC.

Parameters:
y_sizethe number of pixels along the y-axis
Returns:
arrayA (2 * y_size, y_size, 4) array of 0-255 int values.
to_string(format='ascii', fold=0)

Write the MOC into a string.

Format can be ‘ascii’ or ‘json’, though the json format is not officially supported by the IVOA.

Parameters:
formatstr, optional

The format in which the MOC will be serialized before being saved. Possible formats are “ascii” or “json”. By default, format is set to “ascii”.

fold: int

if >0, print ascii or json output with a maximum line width

union(another_moc, *args)

Union between the MOC instance and other MOCs.

Parameters:
another_mocMOC

The MOC used for performing the union with self.

argsMOC

Other additional MOCs to perform the union with.

Returns:
resultMOC/TimeMOC

The resulting MOC.

property uniq_gen

Return a np.array of the generic uniq indices of the cell in the MOC.

Warning

This is not defined in the MOC standard and is not HEALPix scpecific.

Notes

  • It consists on the regular index with a sentinel bit placed at the immediate left of the index’s MSB. At a given depth, the sentinel bit is always put o the same bit.

  • Because the uniq HEALPix encoding is not adapted for non-HEALPIx indices.

property uniq_hpx

Return a np.array of the uniq HEALPIx indices of the cell in the MOC.

property uniq_zorder

Return a np.array of the zorder uniq indices of the cell in the MOC.

Warning

This is not defined in the MOC standard and is not HEALPix specific.

Notes

  • It consists on a regular index shifted on the left so that indices at all level have the same MSB. Plus a sentinel bit placed at the immediate right of the LSB.

  • Because the uniq HEALPix encoding is not adapted for non-HEALPIx indices AND because the natural ordering of such indices follow the same order as the zorder indices (which is very useful for streaming processing, e.g. when dealing with multi-order maps)

wcs(fig, coordsys='icrs', projection='AIT', rotation=None)[source]

Get a wcs that can be given to matplotlib to plot the MOC.

Parameters:
figfigure

The matplotlib figure used for plotting the MOC.

coordsysstr, optional

Coordinate system. Default to “icrs”. Must be in [“icrs”, “galactic”].

projectionstr, optional

World base -> Image base projection type. See http://docs.astropy.org/en/stable/wcs/#supported-projections for the projections currently supported in astropy. Default to Aitoff.

rotationAngle, optional

The angle of rotation. Default to no rotation.

Returns:
wcsWCS

The WCS that can be passed to mocpy.MOC.fill/border.

Examples

>>> from mocpy import MOC
>>> import matplotlib.pyplot as plt
>>> moc = MOC.from_str("2/2-25 28 29 4/0 6/")
>>> fig = plt.figure()
>>> moc.wcs(fig)
WCS Keywords

Number of WCS axes: 2
CTYPE : 'RA---AIT'  'DEC--AIT'
CRVAL : 92.29966711743452  54.33295312309193
CRPIX : 320.0  240.0
PC1_1 PC1_2  : 1.0  -0.0
PC2_1 PC2_2  : 0.0  1.0
CDELT : -0.27794934051515896  0.27794934051515896
NAXIS : 0  0
write(path, format='fits', overwrite=False, optional_kw_dict=None, pre_v2=False)

Write the MOC to a file.

Format can be ‘fits’ or ‘json’, though only the fits format is officially supported by the IVOA.

Parameters:
pathstr

The path to the file to save the MOC in.

formatstr, optional

The format in which the MOC will be serialized before being saved. Possible formats are “fits” or “json”. By default, format is set to “fits”.

overwritebool, optional

If the file already exists and you want to overwrite it, then set the overwrite keyword. Default to False.

optional_kw_dictoptional

Optional keywords arguments added to the FITS header. Only used if format equals to ‘fits’.