mocpy.FrequencyMOC

class mocpy.FrequencyMOC(store_index)[source]

Multi-order frequency coverage class. Experimental.

Is a Frequency Coverage (F-MOC).

Args:

create_key : Object ensure __init__ is called by super-class/class-methods only store_index : index of the FrequencyMOC 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.

complement()

Return the complement of the MOC instance.

Returns:
resultMOC or TimeMOC

The resulting MOC.

contains(frequencies, keep_inside=True)[source]

Test is a frequency – or list of frequencies – is comprised in this FrequencyMOC.

Parameters:
frequenciesastropy.units.Quantity

astropy quantity (converted into Hz) to check whether they are contained in the F-MOC or not.

keep_insidebool, optional

True by default. If so the filtered table contains only observations that are located the MOC. If keep_inside is False, the filtered table contains all observations lying outside the MOC.

Returns:
arraynumpy.ndarray

A mask boolean array

Examples

>>> from mocpy import FrequencyMOC
>>> import astropy.units as u
>>> # Let's create a FMOC at order 10 for frequencies comprised between
>>> # 1Hz and 10Hz.
>>> fmoc = FrequencyMOC.from_frequency_ranges(10, 1 * u.Hz, 10 * u.Hz)
>>> # We can now test wether fmoc contains a list of frequencies between
>>> # 1Hz and 15Hz.
>>> fmoc.contains(range(1, 15, 1) * u.Hz)
array([ True,  True,  True,  True,  True,  True,  True,  True,  True,
       False, False, False, False, False])
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 FrequencyMOC instance to a new, less precise, FrequencyMOC.

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

Parameters:
new_orderint

Maximum depth of the output degraded F-MOC.

Returns:
FrequencyMOC

The degraded F-MOC.

Examples

>>> from mocpy import FrequencyMOC
>>> import astropy.units as u
>>> fmoc = FrequencyMOC.from_frequencies(40, 1 * u.Hz)
>>> fmoc
40/807453851648
>>> fmoc.degrade_to_order(10)
10/752
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.

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

flatten()

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

classmethod from_depth59_ranges(order, ranges)[source]

Create a FrequencyMOC from a set of FrequencyMOC ranges at order 59.

Parameters:
orderint, The resolution of the FrequencyMOC
rangesnumpy.ndarray or list

a N x 2 numpy array or list representing the set of depth 61 ranges.

Returns:
FrequencyMOC

The F-MOC

Examples

>>> from mocpy import FrequencyMOC
>>> FrequencyMOC.from_depth59_ranges(40, [[0, 10000000]])
36/0
38/4
40/
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_frequencies(order, frequencies)[source]

Create a FrequencyMOC from a astropy.units.Quantity that are internally converted in Hertz.

Parameters:
orderint

The resolution of the FrequencyMOC: see relative_precision_to_order

frequenciesastropy.units.Quantity

Quantity converted internally in Hertz

Returns:
FrequencyMOC

Examples

>>> from mocpy import FrequencyMOC
>>> import astropy.units as u
>>> FrequencyMOC.from_frequencies(42, [1e-6, 1e-3, 1] * u.Hz)
42/2544289697882 2887042656632 3229815406592
classmethod from_frequency_ranges(order, min_freq, max_freq)[source]

Create a FrequencyMOC from a range defined by two astropy.units.Quantity converted in Hertz.

Parameters:
orderint

The resolution of the FrequencyMOC: see relative_precision_to_order

min_freqastropy.units.Quantity

astropy quantity converted in Hertz and defining the left part of the intervals

max_freqastropy.units.Quantity

astropy quantity converted in Hertz and defining the right part of the intervals

Returns:
FrequencyMOC

Examples

>>> from mocpy import FrequencyMOC
>>> import astropy.units as u
>>> FrequencyMOC.from_frequency_ranges(10, [10, 40] * u.Hz, [20, 60] * u.Hz)
8/195
9/389 392 397-398
10/798
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_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 FrequencyMOC 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 F-MOC will be serialized before being saved. Possible formats are “ascii” or “json”. By default, format is set to “ascii”.

Examples

>>> from mocpy import FrequencyMOC
>>> FrequencyMOC.from_string("4/4")
4/4
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.

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

Load the FrequencyMOC 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 F-MOC from.

formatstr, optional {‘ascii’, ‘fits’, ‘json’}

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

property max_freq

Get the astropy.units.Quantity largest frequency of the F-MOC.

Returns:
max_freqastropy.units.Quantity

frequency of the last observation

Examples

>>> from mocpy import FrequencyMOC
>>> import astropy.units as u
>>> fmoc = FrequencyMOC.from_frequencies(10, [1, 10] * u.Hz)
>>> # this order is pretty low, thus the returned max frequency
>>> # corresponds to the high limit of the cell containing 10 Hz
>>> # at order 10
>>> print(fmoc.max_freq)
11.0 Hz
property max_index

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

property max_order

Depth/order of the F-MOC.

Returns:
np.uint8

Examples

>>> from mocpy import FrequencyMOC
>>> fmoc = FrequencyMOC.from_json({8: [12, 14, 16], 22: [120, 121, 122]})
>>> fmoc.max_order
22
property min_freq

Get the Quantity frequency of the F-MOC smallest frequency.

Returns:
min_freqastropy.units.Quantity

frequency of the first observation

Examples

>>> from mocpy import FrequencyMOC
>>> import astropy.units as u
>>> fmoc = FrequencyMOC.from_frequencies(10, [1, 10] * u.Hz)
>>> print(fmoc.min_freq)
1.0 Hz
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 FrequencyMOC
>>> FrequencyMOC.n_cells(0)
2
classmethod new_empty(max_depth)[source]

Create a new empty FrequencyMOC of given depth.

Parameters:
max_depthint

The resolution of the FrequencyMOC

Returns:
FrequencyMOC

Examples

>>> from mocpy import FrequencyMOC
>>> FrequencyMOC.new_empty(5)
5/
static order_to_relative_precision(order)[source]

Convert a FrequencyMOC order to a relative precision range.

Parameters:
orderint

order to convert

Returns:
numpy.ndarray

lower and upper relative precisions

Notes

In FMOCs, the precision of a cell depends on its position along the electromagnetic axis. The array returned by order_to_relative_precision corresponds to the lower and upper relative precisions. These must be multiplied by the value of the observed frequency to obtain the absolute upper and lower precisions.

In the code example, we see that at order 10:

\[F = 10_{- 0.6}^{+ 1}~\mathrm{Hz}\]
\[F = 1e3_{- 6e1}^{+ 1e2}~\mathrm{Hz}\]

At order 20 these precisions become:

\[F = 10_{- 6e-4}^{+ 1e-3}~\mathrm{Hz}\]
\[F = 1e3_{- 6e-2}^{+ 1e-1}~\mathrm{Hz}\]

Examples

>>> from mocpy import FrequencyMOC
>>> FrequencyMOC.order_to_relative_precision(10)
array([0.0625, 0.125 ])
>>> FrequencyMOC.order_to_relative_precision(20)
array([6.10351562e-05, 1.22070312e-04])
plot_frequencies(ax, color='blue', frequency_unit='Hz')[source]

Plot a frequency moc.

This method applies a matplotlib.collections.PatchCollection to an existing matplotlib.axes._axes.Axes object.

Parameters:
axmatplotlib.axes._axes.Axes
colorstr, default ‘blue’

any format supported by matplotlib for colors, see matplotlib.colors

length_unitstr or astropy.units.core.Unit, optional

any string or astropy.unit of physical type ‘frequency’, see astropy.units.physical.get_physical_type Defaults to Hertz ‘Hz’

Examples

>>> from mocpy import FrequencyMOC
>>> import matplotlib.pyplot as plt
>>> import astropy.units as u
>>> fmoc = FrequencyMOC.from_frequencies(10, [1, 0.1, 0.01, 0.001] * u.Hz)
>>> fig, ax = plt.subplots(figsize=(15, 1))
>>> fmoc.plot_frequencies(ax, color="pink", frequency_unit="1 / ks")
plot_wavelengths(ax, color='blue', length_unit='m')[source]

Plot a FrequencyMOC with a conversion to wavelengths.

This method applies a matplotlib.collections.PatchCollection to an existing matplotlib.axes._axes.Axes object.

Parameters:
axmatplotlib.axes._axes.Axes
colorstr, default ‘blue’

any format supported by matplotlib for colors, see matplotlib.colors.

length_unitstr or astropy.units.core.Unit, default ‘m’

any string or astropy.unit of physical type ‘length’, see astropy.units.get_physical_type Defaults to meters ‘m’

Examples

>>> from mocpy import FrequencyMOC
>>> import matplotlib.pyplot as plt
>>> import astropy.units as u
>>> fmoc = FrequencyMOC.from_frequencies(10, [1, 0.1, 0.01, 0.001] * u.Hz)
>>> fig, ax = plt.subplots(figsize=(15, 1))
>>> fmoc.plot_wavelengths(ax, color="lightblue", length_unit=u.nm)
static relative_precision_to_order(frequency_precision)[source]

Convert a relative frequency precision into a FrequencyMOC order.

Parameters:
frequency_precisionfloat

precision to be converted in FrequencyMOC order

Returns:
orderint

The less precise order fulfilling the given frequency_precision.

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.

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_depth59_ranges

Return the list of ranges this FrequencyMOC contains, at the maximum precision.

Examples

>>> from mocpy import FrequencyMOC
>>> import astropy.units as u
>>> fmoc = FrequencyMOC.from_frequency_ranges(59, 1 * u.Hz, 1.4 * u.Hz)
>>> print(fmoc.to_depth59_ranges)
[[423338364972826624 425139804823774822]]
to_hz_ranges()[source]

Return the Hertz ranges this FrequencyMOC contains, in Hertz.

Examples

>>> from mocpy import FrequencyMOC
>>> import astropy.units as u
>>> fmoc = FrequencyMOC.from_frequency_ranges(10, [1, 0.1, 0.01] * u.Hz, [1.5, 0.5, 0.05] * u.Hz)
>>> print(fmoc.to_hz_ranges())
[[0.00976562 0.05078125]
 [0.09375    0.5       ]
 [1.         1.5       ]]
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_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)

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’.