hats.inspection.visualize_catalog#
Generate a molleview map with the pixel densities of the catalog
NB: Testing validity of generated plots is currently not tested in our unit test suite.
Functions#
|
Read the object spatial distribution information from a healpix FITS file. |
|
Create a visual map of the density of input points of a catalog on-disk. |
|
Create a visual map of the pixel density of the catalog. |
|
Create a visual map of the pixel density of a list of pixels. |
|
Plots a moc |
|
Returns a MOC that matches the plot window defined by a WCS |
|
Culls a mapping of ipix to values to pixels that are inside the plot window defined by a WCS |
|
Merges any pixels too small in a map to a lower order, with the map values within a lower order pixel |
|
Modified from mocpy.moc.plot.culling_backfacing_cells.from_moc |
|
Compute HEALPix vertices. |
|
Plot a map of HEALPix pixels to values as a colormap across a projection of the sky |
|
Initializes matplotlib Figure and WCSAxes if they do not exist |
|
Perform the plotting of a healpix pixel map. |
Module Contents#
- _read_point_map(catalog_base_dir)[source]#
Read the object spatial distribution information from a healpix FITS file.
- Parameters:
catalog_base_dir – path to a catalog
- Returns:
one-dimensional numpy array of long integers where the value at each index corresponds to the number of objects found at the healpix pixel.
- plot_density(catalog: hats.catalog.Catalog, *, plot_title: str | None = None, order=None, unit=None, **kwargs)[source]#
Create a visual map of the density of input points of a catalog on-disk. :param catalog: :type catalog: hats.catalog.Catalog :param plot_title: Optional title for the plot :type plot_title: str :param order: Optionally reduce the display healpix order, and aggregate smaller tiles. :type order: int :param kwargs: Additional args to pass to plot_healpix_map
- plot_pixels(catalog: hats.catalog.healpix_dataset.healpix_dataset.HealpixDataset, plot_title: str | None = None, **kwargs)[source]#
Create a visual map of the pixel density of the catalog.
- Parameters:
catalog (hats.catalog.Catalog)
plot_title (str) – Optional title for the plot
kwargs – Additional args to pass to plot_healpix_map
- plot_pixel_list(pixels: list[hats.pixel_math.HealpixPixel], plot_title: str = '', projection='MOL', color_by_order=True, **kwargs)[source]#
Create a visual map of the pixel density of a list of pixels.
- Parameters:
pixels – healpix pixels (order and pixel number) to visualize
plot_title (str) – heading for the plot
projection (str) – The projection to use. Available projections listed at https://docs.astropy.org/en/stable/wcs/supported_projections.html
color_by_order (bool) – Whether to color the pixels by their order. True by default.
kwargs – Additional args to pass to plot_healpix_map
- plot_moc(moc: mocpy.MOC, *, projection: str = 'MOL', title: str = '', fov: astropy.units.Quantity | tuple[astropy.units.Quantity, astropy.units.Quantity] = None, center: astropy.coordinates.SkyCoord | None = None, wcs: astropy.wcs.WCS = None, frame_class: Type[astropy.visualization.wcsaxes.frame.BaseFrame] | None = None, ax: astropy.visualization.wcsaxes.WCSAxes | None = None, fig: matplotlib.figure.Figure | None = None, **kwargs) tuple[matplotlib.figure.Figure, astropy.visualization.wcsaxes.WCSAxes] [source]#
Plots a moc
By default, a new matplotlib figure and axis will be created, and the projection will be a Molleweide projection across the whole sky.
- Parameters:
moc (mocpy.MOC) – MOC to plot
projection (str) – The projection to use in the WCS. Available projections listed at https://docs.astropy.org/en/stable/wcs/supported_projections.html
title (str) – The title of the plot
fov (Quantity or Sequence[Quantity, Quantity] | None) – The Field of View of the WCS. Must be an astropy Quantity with an angular unit, or a tuple of quantities for different longitude and latitude FOVs (Default covers the full sky)
center (SkyCoord | None) – The center of the projection in the WCS (Default: SkyCoord(0, 0))
wcs (WCS | None) – The WCS to specify the projection of the plot. If used, all other WCS parameters are ignored and the parameters from the WCS object is used.
frame_class (Type[BaseFrame] | None) – The class of the frame for the WCSAxes to be initialized with. if the ax kwarg is used, this value is ignored (By Default uses EllipticalFrame for full sky projection. If FOV is set, RectangularFrame is used)
ax (WCSAxes | None) – The matplotlib axes to plot onto. If None, an axes will be created to be used. If specified, the axes must be an astropy WCSAxes, and the wcs parameter must be set with the WCS object used in the axes. (Default: None)
fig (Figure | None) – The matplotlib figure to add the axes to. If None, one will be created, unless ax is specified (Default: None)
**kwargs – Additional kwargs to pass to mocpy.MOC.fill
- Returns:
Tuple[Figure, WCSAxes] - The figure and axes used to plot the healpix map
- get_fov_moc_from_wcs(wcs: mocpy.WCS) mocpy.MOC | None [source]#
Returns a MOC that matches the plot window defined by a WCS
Modified from mocpy.moc.plot.utils.build_plotting_moc
- Parameters:
wcs (astropy.WCS) – The wcs object with the plot’s projection
- Returns:
The moc which defines the area of the sky that would be visible in a WCSAxes with the given WCS
- cull_to_fov(depth_ipix_d: dict[int, tuple[numpy.ndarray, numpy.ndarray]], wcs)[source]#
Culls a mapping of ipix to values to pixels that are inside the plot window defined by a WCS
Modified from mocpy.moc.plot.utils.build_plotting_moc
Any pixels too small are merged to a lower order, with the map values within a lower order pixel being sampled
- Parameters:
depth_ipix_d (Dict[int, Tuple[np.ndarray, np.ndarray]]) – Map of HEALPix order to a tuple of 2 arrays (the ipix array of pixel numbers in NESTED ordering, and the values of the pixels)
wcs (astropy.WCS) – The wcs object with the plot’s projection
- Returns:
A new map with the same datatype of depth_ipix_d, with any pixels outside the plot’s FOV removed, and any pixels too small merged with their map values subsampled.
- _merge_too_small_pixels(depth_ipix_d: dict[int, tuple[numpy.ndarray, numpy.ndarray]], wcs)[source]#
Merges any pixels too small in a map to a lower order, with the map values within a lower order pixel being sampled
- cull_from_pixel_map(depth_ipix_d: dict[int, tuple[numpy.ndarray, numpy.ndarray]], wcs, max_split_depth=7)[source]#
Modified from mocpy.moc.plot.culling_backfacing_cells.from_moc
Create a new MOC that do not contain the HEALPix cells that are backfacing the projection.
- Parameters:
depth_ipix_d (Dict[int, Tuple[np.ndarray, np.ndarray]]) – Map of HEALPix order to a tuple of 2 arrays (the ipix array of pixel numbers in NESTED ordering, and the values of the pixels)
wcs (astropy.WCS) – The wcs object with the plot’s projection
max_split_depth – the max depth to split backfacing cells to
- Returns:
A new map with the same datatype of depth_ipix_d, with backfacing cells split into higher order
- compute_healpix_vertices(depth, ipix, wcs, step=1)[source]#
Compute HEALPix vertices.
Modified from mocpy.moc.plot.fill.compute_healpix_vertices
- Parameters:
depth (int) – The depth of the HEALPix cells.
ipix (numpy.ndarray) – The HEALPix cell index given as a np.uint64 numpy array.
wcs (astropy.wcs.WCS) – A WCS projection
- Returns:
path_vertices, codes
- Return type:
(np.array, np.array)
- plot_healpix_map(healpix_map: numpy.ndarray, *, projection: str = 'MOL', title: str = '', cmap: str | matplotlib.colors.Colormap = 'viridis', norm: matplotlib.colors.Normalize | None = None, ipix: numpy.ndarray | None = None, depth: numpy.ndarray | None = None, cbar: bool = True, fov: astropy.units.Quantity | tuple[astropy.units.Quantity, astropy.units.Quantity] = None, center: astropy.coordinates.SkyCoord | None = None, wcs: astropy.wcs.WCS = None, frame_class: Type[astropy.visualization.wcsaxes.frame.BaseFrame] | None = None, ax: astropy.visualization.wcsaxes.WCSAxes | None = None, fig: matplotlib.figure.Figure | None = None, **kwargs)[source]#
Plot a map of HEALPix pixels to values as a colormap across a projection of the sky
Plots the given healpix pixels on a spherical projection defined by a WCS. Colors each pixel based on the corresponding value in a map.
The map can be across all healpix pixels at a given order, or specify a subset of healpix pixels with the ipix and depth parameters.
By default, a new matplotlib figure and axis will be created, and the projection will be a Molleweide projection across the whole sky.
Additional kwargs will be passed to the creation of a matplotlib PathCollection object, which is the artist that draws the tiles.
- Parameters:
healpix_map (np.ndarray) – Array of map values for the healpix tiles. If ipix and depth are not specified, the length of this array will be used to determine the healpix order, and will plot healpix pixels with pixel index corresponding to the array index in NESTED ordering. If ipix and depth are specified, all arrays must be of the same length, and the pixels specified by the ipix and depth arrays will be plotted with their values specified in the healpix_map array.
projection (str) – The projection to use in the WCS. Available projections listed at https://docs.astropy.org/en/stable/wcs/supported_projections.html
title (str) – The title of the plot
cmap (str | Colormap) – The matplotlib colormap to plot with
norm (Normalize | None) – The matplotlib normalization to plot with
ipix (np.ndarray | None) – Array of HEALPix NESTED pixel indices. Must be used with depth, and arrays must be the same length
depth (np.ndarray | None) – Array of HEALPix pixel orders. Must be used with ipix, and arrays must be the same length
cbar (bool) – If True, includes a color bar in the plot (Default: True)
fov (Quantity or Sequence[Quantity, Quantity] | None) – The Field of View of the WCS. Must be an astropy Quantity with an angular unit, or a tuple of quantities for different longitude and latitude FOVs (Default covers the full sky)
center (SkyCoord | None) – The center of the projection in the WCS (Default: SkyCoord(0, 0))
wcs (WCS | None) – The WCS to specify the projection of the plot. If used, all other WCS parameters are ignored and the parameters from the WCS object is used.
frame_class (Type[BaseFrame] | None) – The class of the frame for the WCSAxes to be initialized with. if the ax kwarg is used, this value is ignored (By Default uses EllipticalFrame for full sky projection. If FOV is set, RectangularFrame is used)
ax (WCSAxes | None) – The matplotlib axes to plot onto. If None, an axes will be created to be used. If specified, the axes must be an astropy WCSAxes, and the wcs parameter must be set with the WCS object used in the axes. (Default: None)
fig (Figure | None) – The matplotlib figure to add the axes to. If None, one will be created, unless ax is specified (Default: None)
**kwargs – Additional kwargs to pass to creating the matplotlib PathCollection artist
- Returns:
Tuple[Figure, WCSAxes] - The figure and axes used to plot the healpix map
- initialize_wcs_axes(projection: str = 'MOL', fov: astropy.units.Quantity | tuple[astropy.units.Quantity, astropy.units.Quantity] = None, center: astropy.coordinates.SkyCoord | None = None, wcs: astropy.wcs.WCS = None, frame_class: Type[astropy.visualization.wcsaxes.frame.BaseFrame] | None = None, ax: astropy.visualization.wcsaxes.WCSAxes | None = None, fig: matplotlib.figure.Figure | None = None, **kwargs)[source]#
Initializes matplotlib Figure and WCSAxes if they do not exist
- Parameters:
projection (str) – The projection to use in the WCS. Available projections listed at https://docs.astropy.org/en/stable/wcs/supported_projections.html
fov (Quantity or Sequence[Quantity, Quantity] | None) – The Field of View of the WCS. Must be an astropy Quantity with an angular unit, or a tuple of quantities for different longitude and latitude FOVs (Default covers the full sky)
center (SkyCoord | None) – The center of the projection in the WCS (Default: SkyCoord(0, 0))
wcs (WCS | None) – The WCS to specify the projection of the plot. If used, all other WCS parameters are ignored and the parameters from the WCS object is used.
frame_class (Type[BaseFrame] | None) – The class of the frame for the WCSAxes to be initialized with. if the ax kwarg is used, this value is ignored (By Default uses EllipticalFrame for full sky projection. If FOV is set, RectangularFrame is used)
ax (WCSAxes | None) – The matplotlib axes to plot onto. If None, the current axes will be used. If the current axes is not the correct WCSAxes type, a new figure and axes will be created to be used. If specified, the axes must be an astropy WCSAxes, and the wcs parameter will be ignored and the wcs of the axes used. (Default: None)
fig (Figure | None) – The matplotlib figure to add the axes to. If None, the current figure will be used, unless ax is specified, in which case the figure of the ax will be used. If there is no current figure, one will be created. (Default: None)
kwargs – additional kwargs to pass to figure initialization