hpmoc.plot module
Plotting commands implemented using Astropy (rather than HEALPy) for greater
flexibility, robustness, and cross-platform compatibility. Meant to replace
the plotters in plotters
. These provide similar functionality to the
plotters provided by the ligo.skymap
package as well as healpy
; those
packages are not included as dependencies, but the high-performance rendering
tools included in hpmoc
are totally compatible with those plotting tools
through the PartialUniqSkymap.render
interface (in fact, the previous
version of these plotting scripts used healpy
for its included projection
axes and coordinate transforms).
Available projections are drawn from the FITS standard. You can specify your
own projections using the appropriate FITS headers plugged into an
astropy.wcs.WCS
world coordinate system instance as the projection
argument to a subplot or astropy.visualization.wcsaxes.WCSAxes
instance.
The projection code follows the “4-3” form defined in the
FITS Definition document, i.e. by specifying the coordinate type
(e.g. RA--
for right ascension) right-padded by enough dashes to fill
4 characters, followed by an additional -
and the 3-character projection
code (one of those listed below) to specify the projection. This information is
then stored in the CTYPEi
headers. Pixel-scaling through the CDELTi
headers has a projection-dependent normalization and further depends on pixel
resolution of the final image. Details on WCS header fields are available in
the FITS Definition, and the normalization factors can be extracted from
each projection’s definition in
Representations of celestial coordinates in FITS.
The complete list of available projections can be found in the
FITS Definition, with concrete transformation definitions given in
the Representations of Celestial Coordinates in FITS. The table of
available projections is reproduced below, with sections linking to the
corresponding section in Representations of Celestial Coordinates in FITS.
Note that not all valid WCS projections can be displayed by astropy at time of
writing; in particular, the HEALPIX HPX
projection does not work out of the
box, which is one of the (many) motivations for this plotting library.
Code |
φ_0 |
θ_0 |
Properties1 |
Projection name |
---|---|---|---|---|
AZP |
0◦ |
90◦ |
Sect. 5.1.1 |
Zenithal perspective |
SZP |
0◦ |
90◦ |
Sect. 5.1.2 |
Slant zenithal perspective |
TAN |
0◦ |
90◦ |
Sect. 5.1.3 |
Gnomonic |
STG |
0◦ |
90◦ |
Sect. 5.1.4 |
Stereographic |
SIN |
0◦ |
90◦ |
Sect. 5.1.5 |
Slant orthographic |
ARC |
0◦ |
90◦ |
Sect. 5.1.6 |
Zenithal equidistant |
ZPN |
0◦ |
90◦ |
Sect. 5.1.7 |
Zenithal polynomial |
ZEA |
0◦ |
90◦ |
Sect. 5.1.8 |
Zenithal equal-area |
AIR |
0◦ |
90◦ |
Sect. 5.1.9 |
Airy |
CYP |
0◦ |
0◦ |
Sect. 5.2.1 |
Cylindrical perspective |
CEA |
0◦ |
0◦ |
Sect. 5.2.2 |
Cylindrical equal area |
CAR |
0◦ |
0◦ |
Sect. 5.2.3 |
Plate carrée |
MER |
0◦ |
0◦ |
Sect. 5.2.4 |
Mercator |
SFL |
0◦ |
0◦ |
Sect. 5.3.1 |
Samson-Flamsteed |
PAR |
0◦ |
0◦ |
Sect. 5.3.2 |
Parabolic |
MOL |
0◦ |
0◦ |
Sect. 5.3.3 |
Mollweide |
AIT |
0◦ |
0◦ |
Sect. 5.3.4 |
Hammer-Aitoff |
COP |
0◦ |
θa |
Sect. 5.4.1 |
Conic perspective |
COE |
0◦ |
θa |
Sect. 5.4.2 |
Conic equal-area |
COD |
0◦ |
θa |
Sect. 5.4.3 |
Conic equidistant |
COO |
0◦ |
θa |
Sect. 5.4.4 |
Conic orthomorphic |
BON |
0◦ |
0◦ |
Sect. 5.5.1 |
Bonne’s equal area |
PCO |
0◦ |
0◦ |
Sect. 5.5.2 |
Polyconic |
TSC |
0◦ |
0◦ |
Sect. 5.6.1 |
Tangential spherical cube |
CSC |
0◦ |
0◦ |
Sect. 5.6.2 |
COBE quadrilateralized spherical cube |
QSC |
0◦ |
0◦ |
Sect. 5.6.3 |
Quadrilateralized spherical cube |
HPX |
0◦ |
0◦ |
Sect. 6 2 |
HEALPix grid |
See also
hpmoc.partial.PartialUniqSkymap
, astropy.wcs.WCS
, astropy.visualization.wcsaxes
, matplotlib.axes.Axes
, matplotlib.figure.Figure
, ligo.skymap
, healpy
- hpmoc.plot.dec_exclusions(ax)
Get declination values to exclude from ticks.
- hpmoc.plot.get_colormap(cmap, bad=None)
- hpmoc.plot.get_frame_class(projection: Union[str, astropy.wcs.WCS, astropy.io.fits.Header] = 'MOL', frame_class: Optional[Union[str, astropy.visualization.wcsaxes.frame.BaseFrame]] = None, vdelta: Optional[float] = None, hdelta: Optional[float] = None, rot: Optional[Union[Tuple[float, float, float], Tuple[float, float]]] = None, scatter: List[hpmoc.points.PointsTuple] = (), sigmas: List[float] = (1,)) astropy.visualization.wcsaxes.frame.BaseFrame
Get the frame class associated with a given projection. If already given a
frame_class
, returns it immediately, making this function idempotent. Will otherwise try to determine the best frame choice from the given arguments.- Parameters
projection (str, WCS, or Header, optional) – The projection to use. See
get_wcs
for details. Only used if astr
.frame_class (str or BaseFrame) – The frame class; optionally pass it to ensure that an argument is an instance of
BaseFrame
.vdelta (float, optional) – Manual vdelta override. Give up if set and return
RectangularFrame
.hdelta (float, optional) – Manual hdelta override. Give up if set and return
RectangularFrame
.rot – See the note in
get_wcs
on these arguments. If both rotation and window are set for ARC (zenithal equidistant), switch toRectangularFrame
.scatter – See the note in
get_wcs
on these arguments. If both rotation and window are set for ARC (zenithal equidistant), switch toRectangularFrame
.sigmas – See the note in
get_wcs
on these arguments. If both rotation and window are set for ARC (zenithal equidistant), switch toRectangularFrame
.
- Returns
frame_class – The frame class best-suited to this type of projection.
- Return type
astropy.visualization.wcsaxes.frame.BaseFrame
- Raises
IndexError – If the specified projection could not be found.
- hpmoc.plot.get_projection(projection, *args, **kwargs)
If
projection
is already anastropy.wcs.WCS
instance, return it; if it is anastropy.io.fits.Header
, returnWCS(projection)
. If it is a string, try to create a newWCS
usingget_wcs
with*args
and**kwargs
passed to that function; failing that, try to interpretprojection
as a raw FITS header and instantiate aWCS
therefrom.
- hpmoc.plot.get_ticks(include, exclude, ticks, delta, transform)
- hpmoc.plot.get_wcs(projection: str, width: Optional[int] = None, height: int = 180, hdelta: Optional[float] = None, vdelta: Optional[float] = None, rot: Optional[Union[Tuple[float, float, float], Tuple[float, float]]] = None, facing_sky: bool = True, scatter: List[hpmoc.points.PointsTuple] = (), sigmas: List[float] = (1,)) astropy.wcs.WCS
Get a
WCS
instance by name to match the given parameters.- Parameters
projection (str) –
The following projections are available by default. See
hpmoc.plot
documentation for instructions on constructing your own WCS headers for other plot styles, which can be passed in instead using aWCS
instance. You can also use this approach to plot the skymap over an existingWCS
taken from another fits file, making it easy to plot skymaps over other astrophysical data. This function will returnWCS
instances for the following projection types:Cylindrical/Pseudo-cylindrical (For all-sky plots)
MOL: Mollweide, mollview, Homolographic, Homalographic, Babinet, Elliptical
AIT: Hammer-Aitoff, Aitoff, Hammer, Aitov, Hammer equal-area
CAR: Carée, Plate carée, Caree, Plate caree, Cartesian, Tyre, cartview, Equidistant cylindrical
CEA: Cylindrical equal-area, Lambert equal area
SFL: Sanson-Flamsteed
PAR: Parabolic, Craster
Zenithal
TAN: gnomonic, gnomview, Central, zoom
SIN: Slant Orthographic, Orthographic, Globe, orthview
ARC: Zenithal Equidistant, azimuthal equidistant, azeqview, Postel, Equidistant, Globular
ZEA: Zenithal Equal-area, Azimuthal equivalent, polar azimuthal, Lambert azimuthal equivalent, Lambert azimuthal equal-area, Lambert polar azimuthal, Lambert
width (int, optional) – Width of the image in pixels. If not provided, default to the height for zenithal projections and twice the height for azimuthal projections.
height (int) – Height of the image in pixels.
hdelta (float, optional) – The CDELT1 value, if you wish to override the default. Note that the actual angular width of a pixel at the reference point depends on this value as well as the projection used. Ignored if
ax
is given.vdelta (float, optional) – The CDELT2 value, if you wish to override the default. Note that the actual angular height of a pixel at the reference point depends on this value as well as the projection used. Ignored if
ax
is given.rot (float, float, float) or (float, float), optional) – Euler angles for rotations about the Z, X, Z axes. These are immediately translated to
CRVAL1, CRVAL2, LONPOLE
in the returnedWCS
; that is, the first two angles specify the angle of the center of the image, while the last specifies an additional rotation of the reference pole about the Z-axis. Seehpmoc.plot
documentation for further references. Defaults to centering on RA = 180, dec = 0 (the convention used in most plots produced by LIGO). If only two angles are provided,LONPOLE
will be determined automatically such that the center of the skymap is given by the (RA, dec) given and the orientation of the plot is otherwise unchanged.facing_sky (bool) – Whether the projection is outward- or inward-facing. Equivalent to reversing the direction of the longitude.
scatter (List[PointsTuple]) – A list of collections of point sources that will be plotted. If only one collection containing one point source is provided and
rot=None
(the default), thenrot
will be set to the location of that single point source.sigmas (List[float]) – See
plot
for meaning. Ifsigmas
is non-empty;rot
is set by a single point source inscatter
as described above;hdelta
andvdelta
are bothNone
; and theprojection
is an alias for ARC (zenithal equidistant), thenhdelta
andvdelta
will be set so that the returned frame is 1.5x the size of the largest error region plotted in the smallest axis. This provides an easy way to visually emphasize singular points of interest.
- Returns
wcs – A
WCS
instance that can be used for rendering and plotting.- Return type
astropy.wcs.WCS
- Raises
IndexError – If the specified projection could not be found or if an argument is specified incorrectly.
- hpmoc.plot.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]]]], fig: Optional[Union[matplotlib.figure.Figure, matplotlib.gridspec.GridSpec, dict]] = None, projections: List[Union[str, astropy.wcs.WCS, astropy.io.fits.Header, List[Union[str, astropy.wcs.WCS, astropy.io.fits.Header, astropy.visualization.wcsaxes.WCSAxes]]]] = ('MOL',), scatters: Optional[List[List[hpmoc.points.PointsTuple]]] = None, subplot_height: float = 4, ncols: int = 2, hspace: float = 0.2, wspace: float = 0.2, wshrink: Union[float, List[float]] = 1.0, subplot_kwargs: Optional[List[Optional[List[dict]]]] = None, left=0.0, right=1.0, bottom=0.0, top=1.0, **kwargs) Tuple[matplotlib.gridspec.GridSpec, List[List[astropy.visualization.wcsaxes.WCSAxes]]]
Make a grid plot of multiple skymaps (optionally with scatterplots for each).
- Parameters
*skymaps ('hpmoc.PartialUniqSkymap', array, or (array, array)) – The skymaps to plot. Can be a
PartialUniqSkymap
, a single-resolution HEALPix skymap in NEST ordering only, or a tuple of (pixel values, NUNIQ indices) accepted as the first two arguments toPartialUniqSkymap
.fig (matplotlib.figure.Figure or dict, optional) – The figure to plot to. If not provided, a new figure will be created. If a dictonary is provided, it will be passed as keyword arguments to create a new figure. If a
GridSpec
is provided, then the figure to which it is attached will be used, and thatGridSpec
will be used to define the layout.projections (List[Union[str, WCS, fits.Header, List[Union[str, WCS, fits.Header, WCSAxes]]]], optional) – A list of projections (see the
projection
argument ofplot
) to use for each skymap inskymaps
. If multiple projections are specified, they will be plotted alongside each other; if this makes the figure too wide, change the number of columns in the grid withncols
. You can also pass a list of lists of axes of the type returned by this function (which you should do while also passing aGridSpec
asfig
), allowing you to plot multiple layers of data to the same grid plot.scatters (List[List[PointsTuple]], optional) – Scatterplots to use, one list for each skymap containing the sets of points to plot for that skymap. For any of the skymaps which are
PartialUniqSkymap
instances, you can default to plotting that instance’spoint_sources
by passingNone
instead of a list of point sources; this is also the default behavior ifscatters=None
.subplot_height (float, optional) – The height of each subplot (in inches). Width is automatically determined. Overall figure height will be this times
ncols
. Ignored iffig
is a pre-existing figure, or iffigsize
is specified as a keyword argument infig
(though this is not recommended usage and should be avoided unless you know what you’re doing).ncols (int, optional) – How many columns of skymaps to include in the grid. NB: All subplots for a single skymap are counted as a single column (in contrast to
matplotlib.gridspec.GridSpec
, which counts each subplot in a single column). The number of rows is determined automatically from the number of skymaps provided combined withncols
.hspace (float, optional) – How much vertical space (height) to reserve between subplots.
wspace (float, optional) – How much horizontal space (width) to reserve between subplots.
wshrink (float or list of floats, optional) – Scale the plot widths by this much. Useful if you are adding color bars to preserve spacing. If a list, each element corresponds to a projection.
subplot_kwargs (List[Optional[List[dict]]], optional) – Lists of keyword argument dictionaries that will be used for each subplot. The first index specifies the plotter, and the second index specifies the skymap from
skymaps
. This behavior allows you to easily specify lists of keyword arguments for specific projections (since often) only one of theprojections
requires skymap-specific parameters). NB: These subplot-specific keyword arguments take precedence over**kwargs
for their respective subplots. Projection-related keyword arguments are ignored ifprojections
is a list of lists of axes; see theplot
documentation for further details.left (float) – Bounds for the GridSpec within the plot. Move
top
down, for example, if you want more space for a super-title withplt.suptitle
.right (float) – Bounds for the GridSpec within the plot. Move
top
down, for example, if you want more space for a super-title withplt.suptitle
.bottom (float) – Bounds for the GridSpec within the plot. Move
top
down, for example, if you want more space for a super-title withplt.suptitle
.top (float) – Bounds for the GridSpec within the plot. Move
top
down, for example, if you want more space for a super-title withplt.suptitle
.**kwargs – Keyword arguments applied to all projections; again, see
plot
for usages and caveats.
- Returns
gs (matplotlib.gridspec.GridSpec) – The
GridSpec
defining the shape of the subplots. Access the plotted figure asgs.figure
. Reuse this layout by passinggs
to another invocation ofgridplot
.axs (List[List[‘astropy.visualization.wcsaxes.WCSAxes’]]) – The axes which were plotted. The first index (outer list) corresponds to the projection used, and the second index (inner lists) correspond to the skymap. Reuse the axes and projections by passing
projections=axs
to another invocation ofgridplot
.
- Raises
ValueError – If
scatters
is provided and is ill-formatted or not of the same length asskymaps
; ifsubplot_kwargs
is included and is not of the same length asprojections
, or if its elements are neitherNone
nor lists of the same length asskymaps
; iffig
is passed as aGridSpec
which is not compatible with the rest of the arguments given; or if one of theprojections
is specified as a string but cannot be found in this module.
See also
- hpmoc.plot.plot(skymap: 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]]]], *scatter: hpmoc.points.PointsTuple, vmin: Optional[float] = None, vmax: Optional[float] = None, cmap: Optional[Union[str, matplotlib.colors.Colormap]] = 'gist_heat_r', missing_color: Optional[Union[str, Tuple[int, int, int]]] = None, nan_color: Optional[Union[str, Tuple[int, int, int]]] = None, alpha: float = 1.0, sigmas: Iterable[float] = (1,), scatter_labels: Union[bool, dict] = True, ax: Optional[Union[astropy.visualization.wcsaxes.WCSAxes, astropy.visualization.wcsaxes.WCSAxesSubplot]] = None, projection: Union[str, astropy.wcs.WCS, astropy.io.fits.Header] = 'Mollweide', frame_class: Optional[Union[str, astropy.visualization.wcsaxes.frame.BaseFrame]] = None, width: Optional[int] = None, height: int = 180, hdelta: Optional[float] = None, vdelta: Optional[float] = None, rot: Optional[Union[Tuple[float, float, float], Tuple[float, float]]] = None, facing_sky: bool = True, fig: Optional[Union[matplotlib.figure.Figure, dict]] = None, subplot: Optional[Union[Tuple[int, int, int], matplotlib.gridspec.SubplotSpec]] = None, cr: Iterable[float] = (), cr_format: Callable[[float, float], str] = None, cr_filled: bool = False, cr_kwargs: Optional[dict] = None, cr_label_kwargs: Optional[dict] = None, pixels: Union[bool, dict] = False, pixels_format: Optional[Callable[[hpmoc.PartialUniqSkymap], Iterable[str]]] = None, cbar: Union[bool, dict] = False) Union[astropy.visualization.wcsaxes.WCSAxes, astropy.visualization.wcsaxes.WCSAxesSubplot]
- Parameters
skymap ('hpmoc.PartialUniqSkymap', array, or (array, array)) – The skymap to plot. Can be a
PartialUniqSkymap
, a single-resolution HEALPix skymap in NEST ordering only, or a tuple of (pixel values, NUNIQ indices) accepted as the first two arguments toPartialUniqSkymap
.scatter (PointsTuple) – Point-sources to plot as a scatter-map, with disks showing their error regions. Provide multiple (ideally with different colors) to plot many populations at once.
vmin (float, optional) – The smallest value in the color map used to plot
skymap
. SetNone
to have it calculated automatically.vmax (float, optional) – The largest value in the color map used to plot
skymap
. SetNone
to have it calculated automatically.cmap (str or matplotlib.colors.Colormap) – The color map to use to plot
skymap
. Note that the colors for the point sources inscatter
are set using thergba
parameter inPointsTuple
and will not be affected by this value. IfNone
, the skymap itself will not be plotted; this can be useful if overlaying multiple skymaps.missing_color (str or (int, int, int), optional) – The color to use for missing parts of the skymap. If not provided, they simply will not be shown.
nan_color (str or (int, int, int), optional) – The color to use for parts of the map that are included but equal to
np.nan
. If not provided, will be the same asmissing_color
(meaning transparent if that argument is not provided), i.e. they will be plotted as if the pixels were missing. Will not use thealpha
argument (to allow for highlighting NaN values more clearly; to use the samealpha
value, you can create an RGBA tuple including the desired color and your shared alpha value, e.g.nan_color=[*matplotlib.colors.to_rgb('pink'), alpha]
.alpha (float, optional) – The opacity of the plotted skymap image.
sigmas (Iterable[float], optional) – The size of the error region about each point source to plot in units of its error parameter sigma.
scatter_labels (bool, optional) – Whether to show labels for the scattered points. If
True
, display either their labels (if defined) or their indices within thePointsTuple.points
list. If given as adict
, will be interpreted as keyword arguments to pass toWCSAxes.text
, which can be used to control the appearance of the point source labels.ax (WCSAxes or WCSAxesSubplot, optional) – Axes to plot to. If provided, all other arguments pertaining to creating a
WCS
andWCSAxes
instance are ignored, and these axes are used instead.projection (str, WCS, or Header, optional) – Either provide the name of the projection (see
get_wcs
docstring for valid names) to create a newWCS
for this plot, or provide a ready-madeWCS
instance or FITSHeader
from which such an instance can be crafted. In the first case, you will need to specify other parameters needed to fully define the world coordinate system. In the latter cases, you might need to customizeframe
. Ignored ifax
is given.frame_class (BaseFrame or str, optional) – The frame type to use for this plot, e.g. a
RectangularFrame
(for a plate carée/Cartesian plot) or anEllipticalFrame
for a Mollweide plot. Selected automatically whenprojection
is specified by name, otherwise defaults toRectangularFrame
. You can also specify frame_class=’rectangular’ or frame_class=’elliptical’ to choose one of these two frames.width (int, optional) – Width of the image in pixels. If not provided, default to the height for zenithal projections and twice the height for azimuthal projections. Ignored if
ax
is given or ifprojection
is aWCS
instance.height (int, optional) – The height of the plot in pixels. Ignored if
ax
is given or ifprojection
is aWCS
instance.hdelta (float, optional) – The CDELT1 value, if you wish to override the default. Note that the actual angular width of a pixel at the reference point depends on this value as well as the projection used. Ignored if
ax
is given or ifprojection
is aWCS
instance.vdelta (float, optional) – The CDELT2 value, if you wish to override the default. Note that the actual angular height of a pixel at the reference point depends on this value as well as the projection used. Ignored if
ax
is given or ifprojection
is aWCS
instance.rot ((float, float, float) or (float, float), optional) – Euler angles for rotations about the Z, X, Z axes. These are immediately translated to
CRVAL1, CRVAL2, LONPOLE
in the returnedWCS
; that is, the first two angles specify the angle of the center of the image, while the last specifies an additional rotation of the reference pole about the Z-axis. Seehpmoc.plot
documentation for further references. Defaults to centering on RA = 180, dec = 0 (the convention used in most plots produced by LIGO). If only two angles are provided,LONPOLE
will be determined automatically such that the center of the skymap is given by the (RA, dec) given and the orientation of the plot is otherwise unchanged. Ignored ifax
is given or ifprojection
is aWCS
instance.facing_sky (bool, optional) – Whether the projection is outward- or inward-facing. Equivalent to reversing the direction of the longitude. Ignored if
ax
is given or ifprojection
is aWCS
instance.fig (matplotlib.figure.Figure or dict, optional) – The figure to plot to. If not provided, a new figure will be created. If a dictonary is provided, it will be passed as keyword arguments to create a new figure. Ignored if
ax
is given.subplot ((int, int, int) or SubplotSpec, optional) – If provided, initialize the plot as a subplot using the standard
Figure.subplot
matplotlib interface, returning aWCSAxesSubplot
instance rather than aWCSAxes
instance. Ignored ifax
is given.cr (Iterable[float], optional) – If provided, plot contour lines around the credible regions specified in this list. For example,
cr=[0.9]
will plot contours around the smallest region containing 90% of the skymap’s integrated value.cr_format (Callable[[float, float], str], optional) – A function taking the CR level (e.g.
0.9
for 90% CR) and the actual value of the skymap on that contour and returning a string to be used to label each contour. If not provided, contours will not be labeled. Ignored ifcr
is empty.cr_filled (bool, optional) – Whether to fill the contours using
contourf
rather thancontour
. Ignored ifcr
is empty.cr_kwargs (dict, optional) – Additional arguments to pass to either
contour
orcontourf
. Ignored ifcr
is empty.cr_label_kwargs (dict, optional) – Arguments to pass to
clabel
governing the display format of contour labels.pixels (bool or dict, optional) – Whether to plot pixel borders. You should probably keep this
False
if you are doing an all-sky plot, since border plotting is slow for a large number of pixels, and the boundaries will not be visible anyway unless each visible pixel’s size is comparable to the overall size of the plot window. If you just want to see information about e.g. the size of pixels across the whole sky, consider plottingskymap.orders(as_skymap=True)
to see a color map of pixel sizes instead. If given as adict
, will be interpreted as keyword arguments to pass toWCSAxes.plot
, which can be used to control the appearance of the pixel borders.pixels_format (Callable[[PartialUniqSkymap], Iterable[str]], optional) – A function that takes a skymap and returns a string for each pixel. Will be called on the selection of pixels overlapping with the visible area, so don’t worry about optimizing it, since just plotting the borders will be slow enough for a large number of visible pixels. The returned string will be plotted over the center of each pixel, which again is only useful if the pixels are large with respect to the size of the plot window.
cbar (bool or dict, optional) – If
True
, add a colorbar for the plotted skymap. If adict
is provided, pass it as the keyword arguments tomatplotlib.pyplot.colorbar
. Ignored ifcmap=None
.
- Returns
ax – The axes that were just plotted to.
- Return type
WCSAxes or WCSAxesSubplot
See also
hpmoc.partial.PartialUniqSkymap
,matplotlib.figure.Figure
,matplotlib.axes.Axes
,matplotlib.pyplot.colorbar
,matplotlib.gridspec.GridSpec
,matplotlib.gridspec.SubplotSpec
,astropy.io.fits.Header
,astropy.wcs.WCS
,astropy.visualization.wcsaxes.WCSAxes
- hpmoc.plot.ra_exclusions(ax)
Get RA values to exclude from ticks.
- hpmoc.plot.register_dec_exclusion(cls, func)
- hpmoc.plot.register_ra_exclusion(cls, func)