LazyMask¶
- class spectral_cube.LazyMask(function, cube=None, data=None, wcs=None)[source]¶
Bases:
MaskBase
A boolean mask defined by the evaluation of a function on a fixed dataset.
This is conceptually identical to a fixed boolean mask as in
BooleanArrayMask
but defers the evaluation of the mask until it is needed.- Parameters
- functioncallable
The function to apply to
data
. This method should accept a numpy array, which will be a subset of the data array passed to __init__. It should return a boolean array, whereTrue
values indicate that which pixels are valid/unaffected by masking.- dataarray-like
The array to evaluate
function
on. This should support Numpy-like slicing syntax.- wcs
WCS
The WCS of the input data, which is used to define the coordinates for which the boolean mask is defined.
Attributes Summary
Methods Summary
any
()exclude
([data, wcs, view])Return a boolean array indicating which values should be excluded.
include
([data, wcs, view])Return a boolean array indicating which values should be included.
quicklook
(view[, wcs, filename, use_aplpy, ...])View a 2D slice of the mask, specified by view.
view
([view])Compatibility tool: if a numpy.ma.ufunc is run on the mask, it will try to grab a view of the mask, which needs to appear to numpy as a true array.
with_spectral_unit
(unit[, ...])Get a LazyMask copy with a WCS in the modified unit
Attributes Documentation
- dtype¶
- ndim¶
- shape¶
- size¶
Methods Documentation
- any()¶
- exclude(data=None, wcs=None, view=(), **kwargs)¶
Return a boolean array indicating which values should be excluded.
If
view
is passed, only the sliced mask will be returned, which avoids having to load the whole mask in memory. Otherwise, the whole mask is returned in-memory.kwargs are passed to _validate_wcs
- include(data=None, wcs=None, view=(), **kwargs)¶
Return a boolean array indicating which values should be included.
If
view
is passed, only the sliced mask will be returned, which avoids having to load the whole mask in memory. Otherwise, the whole mask is returned in-memory.kwargs are passed to _validate_wcs
- quicklook(view, wcs=None, filename=None, use_aplpy=True, aplpy_kwargs={})¶
View a 2D slice of the mask, specified by view.
- Parameters
- viewtuple
Slicing to apply to the mask. Must return a 2D slice.
- wcsastropy.wcs.WCS, optional
WCS object to use in plotting the mask slice.
- filenamestr, optional
Filename of the output image. Enables saving of the plot.
- use_aplpybool, optional
Try plotting with the aplpy package
- aplpy_kwargsdict, optional
kwargs passed to
FITSFigure
.
- view(view=())¶
Compatibility tool: if a numpy.ma.ufunc is run on the mask, it will try to grab a view of the mask, which needs to appear to numpy as a true array. This can be important for, e.g., plotting.
Numpy’s convention is that masked=True means “masked out”
Note
I don’t know if there are broader concerns or consequences from including this ‘view’ tool here.
- with_spectral_unit(unit, velocity_convention=None, rest_value=None)[source]¶
Get a LazyMask copy with a WCS in the modified unit
- Parameters
- unitu.Unit
Any valid spectral unit: velocity, (wave)length, or frequency. Only vacuum units are supported.
- velocity_conventionu.doppler_relativistic, u.doppler_radio, or u.doppler_optical
The velocity convention to use for the output velocity axis. Required if the output type is velocity.
- rest_valueu.Quantity
A rest wavelength or frequency with appropriate units. Required if output type is velocity. The cube’s WCS should include this already if the input type is velocity, but the WCS’s rest wavelength/frequency can be overridden with this parameter.