mocpy.TimeMOC

class mocpy.TimeMOC(interval_set=None)[source]

Multi-order time coverage class. Experimental

Attributes:
LARK_PARSER_STR
consistency

Get a percentage of fill between the min and max time the moc is defined.

max_order

Depth of the smallest HEALPix cells found in the MOC instance.

max_time

Get the Time time of the tmoc last observation

min_time

Get the Time time of the tmoc first observation

total_duration

Get the total duration covered by the temporal moc

Methods

add_neighbours(self) Add all the pixels at max order in the neighbourhood of the moc
complement(self) Returns the complement of the MOC instance.
contains(self, times[, keep_inside, delta_t]) Get a mask array (e.g.
degrade_to_order(self, new_order) Degrades the MOC instance to a new, less precise, MOC.
difference(self, another_moc[, delta_t]) Difference between self and moc.
empty(self) Checks whether the MOC is empty.
from_fits(filename) Loads a MOC from a FITS file.
from_json(json_moc) Creates a MOC from a dictionary of HEALPix cell arrays indexed by their depth.
from_str(value) Create a MOC from a str.
from_time_ranges(min_times, max_times[, delta_t]) Create a TimeMOC from a range defined by two astropy.time.Time
from_times(times[, delta_t]) Create a TimeMOC from a astropy.time.Time
intersection(self, another_moc[, delta_t]) Intersection between self and moc.
order_to_time_resolution(order) Convert an TimeMoc order to its equivalent time
plot(self[, title, view]) Plot the TimeMoc in a time window.
remove_neighbours(self) Remove all the pixels at max order located at the bound of the moc
serialize(self[, format, optional_kw_dict]) Serializes the MOC into a specific format.
time_resolution_to_order(delta_time) Convert a time resolution to a TimeMoc order.
union(self, another_moc[, delta_t]) Union between self and moc.
write(self, path[, format, overwrite, …]) Writes the MOC to a file.
refine_to_order  
add_neighbours(self)[source]

Add all the pixels at max order in the neighbourhood of the moc

complement(self)

Returns the complement of the MOC instance.

Returns:
result : MOC or TimeMOC

The resulting MOC.

consistency

Get a percentage of fill between the min and max time the moc is defined.

A value near 0 shows a sparse temporal moc (i.e. the moc does not cover a lot of time and covers very distant times. A value near 1 means that the moc covers a lot of time without big pauses.

Returns:
result : float

fill percentage (between 0 and 1.)

contains(self, times, keep_inside=True, delta_t=<TimeDelta object: scale='tdb' format='sec' value=1800.0>)[source]

Get a mask array (e.g. a numpy boolean array) of times being inside (or outside) the TMOC instance.

Parameters:
times : astropy.time.Time

astropy times to check whether they are contained in the TMOC or not.

keep_inside : bool, 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.

delta_t : astropy.time.TimeDelta, optional

the duration of one observation. It is set to 30 min by default. This data is used to compute the more efficient TimeMOC order to represent the observations (Best order = the less precise order which is able to discriminate two observations separated by delta_t).

Returns:
array : darray

A mask boolean array

degrade_to_order(self, new_order)

Degrades 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_order : int

Maximum depth of the output degraded MOC.

Returns:
moc : MOC or TimeMOC

The degraded MOC.

difference(self, another_moc, delta_t=<TimeDelta object: scale='tdb' format='sec' value=1800.0>)[source]

Difference between self and moc. delta_t gives the possibility to the user to set a time resolution for performing the tmoc diff

Parameters:
another_moc : AbstractMOC

the MOC/TimeMoc to substract from self

delta_t : TimeDelta, optional

the duration of one observation. It is set to 30 min by default. This data is used to compute the more efficient TimeMoc order to represent the observations. (Best order = the less precise order which is able to discriminate two observations separated by delta_t)

Returns:
result : MOC or TimeMoc

MOC object whose interval set corresponds to : self - moc

empty(self)

Checks whether the MOC is empty.

A MOC is empty when its list of HEALPix cell ranges is empty.

Returns:
result: bool

True if the MOC instance is empty.

classmethod from_fits(filename)

Loads 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:
filename : str

The path to the FITS file.

Returns:
result : MOC or TimeMOC

The resulting MOC.

classmethod from_json(json_moc)

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

Parameters:
json_moc : dict(str

A dictionary of HEALPix cell arrays indexed by their depth.

Returns:
moc : MOC 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:
value : str

The MOC as a string following the grammar rules.

Returns:
moc : MOC or TimeMOC

The resulting MOC

Examples

>>> from mocpy import MOC
>>> moc = MOC.from_str("2/2-25,28,29 4/0 6/")
classmethod from_time_ranges(min_times, max_times, delta_t=<TimeDelta object: scale='tdb' format='sec' value=1800.0>)[source]

Create a TimeMOC from a range defined by two astropy.time.Time

Parameters:
min_times : astropy.time.Time

astropy times defining the left part of the intervals

max_times : astropy.time.Time

astropy times defining the right part of the intervals

delta_t : astropy.time.TimeDelta, optional

the duration of one observation. It is set to 30 min by default. This data is used to compute the more efficient TimeMOC order to represent the observations (Best order = the less precise order which is able to discriminate two observations separated by delta_t).

Returns:
time_moc : TimeMOC
classmethod from_times(times, delta_t=<TimeDelta object: scale='tdb' format='sec' value=1800.0>)[source]

Create a TimeMOC from a astropy.time.Time

Parameters:
times : astropy.time.Time

astropy observation times

delta_t : astropy.time.TimeDelta, optional

the duration of one observation. It is set to 30 min by default. This data is used to compute the more efficient TimeMOC order to represent the observations (Best order = the less precise order which is able to discriminate two observations separated by delta_t).

Returns:
time_moc : TimeMOC
intersection(self, another_moc, delta_t=<TimeDelta object: scale='tdb' format='sec' value=1800.0>)[source]

Intersection between self and moc. delta_t gives the possibility to the user to set a time resolution for performing the tmoc intersection

Parameters:
another_moc : AbstractMOC

the MOC/TimeMOC used for performing the intersection with self

delta_t : TimeDelta, optional

the duration of one observation. It is set to 30 min by default. This data is used to compute the more efficient TimeMoc order to represent the observations. (Best order = the less precise order which is able to discriminate two observations separated by delta_t)

Returns:
result : MOC or TimeMOC

MOC object whose interval set corresponds to : self & moc

max_order

Depth of the smallest HEALPix cells found in the MOC instance.

max_time

Get the Time time of the tmoc last observation

Returns:
max_time : Time

time of the last observation

min_time

Get the Time time of the tmoc first observation

Returns:
min_time : astropy.time.Time

time of the first observation

static order_to_time_resolution(order)[source]

Convert an TimeMoc order to its equivalent time

Parameters:
order : int

order to convert

Returns:
delta_t : TimeDelta

time equivalent to order

plot(self, title='TimeMoc', view=(None, None))[source]

Plot the TimeMoc in a time window.

This method uses interactive matplotlib. The user can move its mouse through the plot to see the time (at the mouse position).

Parameters:
title : str, optional

The title of the plot. Set to ‘TimeMoc’ by default.

view : (Time, Time), optional

Define the view window in which the observations are plotted. Set to (None, None) by default (i.e. all the observation time window is rendered).

remove_neighbours(self)[source]

Remove all the pixels at max order located at the bound of the moc

serialize(self, format='fits', optional_kw_dict=None)

Serializes the MOC into a specific format.

Possible formats are FITS, JSON and STRING

Parameters:
format : str

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

optional_kw_dict : dict

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

Returns:
result : astropy.io.fits.HDUList or JSON dictionary

The result of the serialization.

static time_resolution_to_order(delta_time)[source]

Convert a time resolution to a TimeMoc order.

Parameters:
delta_time : TimeDelta

time to convert

Returns:
order : int

The less precise order which is able to discriminate two observations separated by delta_time.

total_duration

Get the total duration covered by the temporal moc

Returns:
duration : TimeDelta

total duration of all the observation times of the tmoc total duration of all the observation times of the tmoc

union(self, another_moc, delta_t=<TimeDelta object: scale='tdb' format='sec' value=1800.0>)[source]

Union between self and moc. delta_t gives the possibility to the user to set a time resolution for performing the tmoc union

Parameters:
another_moc : AbstractMOC

the MOC/TimeMoc to bind to self

delta_t : TimeDelta, optional

the duration of one observation. It is set to 30 min by default. This data is used to compute the more efficient TimeMoc order to represent the observations. (Best order = the less precise order which is able to discriminate two observations separated by delta_t)

Returns:
result : MOC or TimeMoc

MOC object whose interval set corresponds to : self | moc

write(self, path, format='fits', overwrite=False, optional_kw_dict=None)

Writes the MOC to a file.

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

Parameters:
path : str, optional

The path to the file to save the MOC in.

format : str, 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”.

overwrite : bool, optional

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

optional_kw_dict : optional

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