diff --git a/src/.pylintrc b/src/.pylintrc index b25fe570..1762c338 100644 --- a/src/.pylintrc +++ b/src/.pylintrc @@ -293,7 +293,7 @@ max-bool-expr=5 max-branches=20 # Maximum number of locals for function / method body. -max-locals=20 +max-locals=25 # Maximum number of parents for a class (see R0901). max-parents=7 diff --git a/src/hats/catalog/healpix_dataset/healpix_dataset.py b/src/hats/catalog/healpix_dataset/healpix_dataset.py index bfb94987..75c68c38 100644 --- a/src/hats/catalog/healpix_dataset/healpix_dataset.py +++ b/src/hats/catalog/healpix_dataset/healpix_dataset.py @@ -15,6 +15,7 @@ from hats.catalog.dataset import Dataset from hats.catalog.dataset.table_properties import TableProperties from hats.catalog.partition_info import PartitionInfo +from hats.inspection import plot_pixels from hats.io import file_io, paths from hats.io.file_io import read_parquet_metadata from hats.pixel_math import HealpixPixel @@ -225,3 +226,11 @@ def align( return align_with_mocs( self.pixel_tree, other_cat.pixel_tree, self.moc, other_cat.moc, alignment_type=alignment_type ) + + def plot_pixels(self, **kwargs): + """Create a visual map of the pixel density of the catalog. + + Args: + kwargs: Additional args to pass to `hipscat.inspection.visualize_catalog.plot_healpix_map` + """ + return plot_pixels(self, **kwargs) diff --git a/src/hats/inspection/visualize_catalog.py b/src/hats/inspection/visualize_catalog.py index 7f1e3765..90c424a9 100644 --- a/src/hats/inspection/visualize_catalog.py +++ b/src/hats/inspection/visualize_catalog.py @@ -3,17 +3,38 @@ NB: Testing validity of generated plots is currently not tested in our unit test suite. """ -from typing import List +from __future__ import annotations -import healpy as hp -import matplotlib.colors as mcolors +import warnings +from typing import TYPE_CHECKING, Dict, List, Tuple, Type + +import astropy.units as u +import astropy.wcs +import cdshealpix import numpy as np +from astropy.coordinates import ICRS, Angle, SkyCoord +from astropy.units import Quantity +from astropy.visualization.wcsaxes.frame import BaseFrame, EllipticalFrame +from astropy.wcs.utils import pixel_to_skycoord, skycoord_to_pixel from matplotlib import pyplot as plt +from matplotlib.axes import Axes +from matplotlib.collections import PathCollection +from matplotlib.colors import Colormap, Normalize +from matplotlib.figure import Figure +from matplotlib.path import Path +from mocpy import MOC, WCS +from mocpy.moc.plot.culling_backfacing_cells import backface_culling +from mocpy.moc.plot.fill import compute_healpix_vertices +from mocpy.moc.plot.utils import _set_wcs -import hats.pixel_math.healpix_shim as hs -from hats.catalog import Catalog from hats.io import file_io, paths from hats.pixel_math import HealpixPixel +from hats.pixel_tree.moc_filter import perform_filter_by_moc +from hats.pixel_tree.pixel_tree import PixelTree + +if TYPE_CHECKING: + from hats.catalog import Catalog + from hats.catalog.healpix_dataset.healpix_dataset import HealpixDataset def _read_point_map(catalog_base_dir): @@ -29,118 +50,377 @@ def _read_point_map(catalog_base_dir): return file_io.read_fits_image(map_file_pointer) -def plot_points(catalog: Catalog, projection="moll", **kwargs): +def plot_points(catalog: Catalog, **kwargs): """Create a visual map of the input points of an in-memory catalog. Args: catalog (`hats.catalog.Catalog`) Catalog to display - projection (str) The map projection to use. Valid values include: - - moll - Molleweide projection (default) - - gnom - Gnomonic projection - - cart - Cartesian projection - - orth - Orthographic projection + kwargs: Additional args to pass to `plot_healpix_map` """ if not catalog.on_disk: raise ValueError("on disk catalog required for point-wise visualization") point_map = _read_point_map(catalog.catalog_base_dir) - _plot_healpix_map(point_map, projection, f"Catalog point density map - {catalog.catalog_name}", **kwargs) + return plot_healpix_map(point_map, title=f"Catalog point density map - {catalog.catalog_name}", **kwargs) -def plot_pixels(catalog: Catalog, projection="moll", **kwargs): +def plot_pixels(catalog: HealpixDataset, **kwargs): """Create a visual map of the pixel density of the catalog. Args: catalog (`hats.catalog.Catalog`) Catalog to display - projection (str) The map projection to use. Valid values include: - - moll - Molleweide projection (default) - - gnom - Gnomonic projection - - cart - Cartesian projection - - orth - Orthographic projection + kwargs: Additional args to pass to `plot_healpix_map` """ pixels = catalog.get_healpix_pixels() - plot_pixel_list( + return plot_pixel_list( pixels=pixels, plot_title=f"Catalog pixel density map - {catalog.catalog_name}", - projection=projection, **kwargs, ) -def plot_pixel_list(pixels: List[HealpixPixel], plot_title: str = "", projection="moll", **kwargs): +def plot_pixel_list(pixels: List[HealpixPixel], plot_title: str = "", projection="MOL", **kwargs): """Create a visual map of the pixel density of a list of pixels. Args: pixels: healpix pixels (order and pixel number) to visualize plot_title (str): heading for the plot - projection (str) The map projection to use. Valid values include: - - moll - Molleweide projection (default) - - gnom - Gnomonic projection - - cart - Cartesian projection - - orth - Orthographic projection + projection (str): The projection to use. Available projections listed at + https://docs.astropy.org/en/stable/wcs/supported_projections.html + kwargs: Additional args to pass to `plot_healpix_map` """ - max_order = np.max(pixels).order - min_order = np.min(pixels).order - - if max_order == min_order: - colors = [plt.cm.viridis(0.0), plt.cm.viridis(0.1)] # pylint: disable=no-member - cmap = mcolors.LinearSegmentedColormap.from_list("my_colormap", colors, 1) - kwargs["cbar"] = False - plot_title = f"Norder {max_order} {plot_title}" - else: - num_colors = max_order - min_order + 1 - colors = plt.cm.viridis(np.linspace(0, 1, num_colors)) # pylint: disable=no-member - cmap = mcolors.LinearSegmentedColormap.from_list("my_colormap", colors, num_colors) - - order_map = np.full(hs.order2npix(max_order), hs.unseen_pixel()) - - for pixel in pixels: - explosion_factor = 4 ** (max_order - pixel.order) - exploded_pixels = [ - *range( - pixel.pixel * explosion_factor, - (pixel.pixel + 1) * explosion_factor, - ) - ] - order_map[exploded_pixels] = pixel.order - _plot_healpix_map(order_map, projection, plot_title, cmap=cmap, **kwargs) - - -def get_projection_method(projection): - """Get the healpy plotting method for a specified projection string + orders = np.array([p.order for p in pixels]) + ipix = np.array([p.pixel for p in pixels]) + order_map = orders.copy() + fig, ax = plot_healpix_map( + order_map, projection=projection, title=plot_title, ipix=ipix, depth=orders, cbar=False, **kwargs + ) + col = ax.collections[0] + col_array = col.get_array() + plt.colorbar( + col, + boundaries=np.arange(np.min(col_array) - 0.5, np.max(col_array) + 0.6, 1), + ticks=np.arange(np.min(col_array), np.max(col_array) + 1), + label="order", + ) + return fig, ax + + +def cull_to_fov(depth_ipix_d: Dict[int, Tuple[np.ndarray, np.ndarray]], wcs): + """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 + + Args: + 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. + """ + + depth_ipix_d = _merge_too_small_pixels(depth_ipix_d, wcs) + + # Get the MOC delimiting the FOV polygon + width_px = int(wcs.wcs.crpix[0] * 2.0) # Supposing the wcs is centered in the axis + height_px = int(wcs.wcs.crpix[1] * 2.0) + + # Compute the sky coordinate path delimiting the viewport. + # It consists of a closed polygon of (4 - 1)*4 = 12 vertices + x_px = np.linspace(0, width_px, 4) + y_px = np.linspace(0, height_px, 4) + + x, y = np.meshgrid(x_px, y_px) + + x_px = np.append(x[0, :-1], x[:-1, -1]) + x_px = np.append(x_px, x[-1, 1:][::-1]) + x_px = np.append(x_px, x[:-1, 0]) + + y_px = np.append(y[0, :-1], y[:-1, -1]) + y_px = np.append(y_px, y[-1, :-1]) + y_px = np.append(y_px, y[1:, 0][::-1]) + + # Disable the output of warnings when encoutering NaNs. + warnings.filterwarnings("ignore") + # Inverse projection from pixel coordinate space to the world coordinate space + viewport = pixel_to_skycoord(x_px, y_px, wcs) + # If one coordinate is a NaN we exit the function and do not go further + ra_deg, dec_deg = viewport.icrs.ra.deg, viewport.icrs.dec.deg + warnings.filterwarnings("default") + + if np.isnan(ra_deg).any() or np.isnan(dec_deg).any(): + return depth_ipix_d + + # Create a rough MOC (depth=3 is sufficient) from the viewport + moc_viewport = MOC.from_polygon_skycoord(viewport, max_depth=3) + + output_dict = {} + + # The moc to plot is the INPUT_MOC & MOC_VIEWPORT. For small FOVs this can reduce + # a lot the time to draw the MOC along with its borders. + for d, (ip, vals) in depth_ipix_d.items(): + ip_argsort = np.argsort(ip) + ip_sorted = ip[ip_argsort] + pixel_tree = PixelTree(np.vstack([ip_sorted, ip_sorted + 1]).T, order=d) + ip_viewport_mask = perform_filter_by_moc( + pixel_tree.to_depth29_ranges(), moc_viewport.to_depth29_ranges + ) + output_dict[d] = (ip_sorted[ip_viewport_mask], vals[ip_argsort][ip_viewport_mask]) + return output_dict + + +def _merge_too_small_pixels(depth_ipix_d: Dict[int, Tuple[np.ndarray, np.ndarray]], wcs): + """Merges any pixels too small in a map to a lower order, with the map values within a lower order pixel + being sampled""" + # Get the WCS cdelt giving the deg.px^(-1) resolution. + cdelt = wcs.wcs.cdelt + # Convert in rad.px^(-1) + cdelt = np.abs((2 * np.pi / 360) * cdelt[0]) + # Get the minimum depth such as the resolution of a cell is contained in 1px. + depth_res = int(np.floor(np.log2(np.sqrt(np.pi / 3) / cdelt))) + depth_res = max(depth_res, 0) + + max_depth = max(depth_ipix_d.keys()) + + # Combine healpix pixels smaller than 1px in the plot + if max_depth > depth_res: + warnings.warn( + "This plot contains HEALPix pixels smaller than a pixel of the plot. Some values may be lost" + ) + new_ipix_d = {} + for d, (ip, vals) in depth_ipix_d.items(): + if d <= depth_res: + new_ipix_d[d] = (ip, vals) + else: + ipix_depth_res = ip >> (2 * (d - depth_res)) + # Get the unique pixels at the maximum depth resolution + unique_ipix, unique_indices = np.unique(ipix_depth_res, return_index=True) + # Get the first values from the map for each lower order pixel + vals_depth_res = vals[unique_indices] + if depth_res not in new_ipix_d: + new_ipix_d[depth_res] = (unique_ipix, vals_depth_res) + else: + # combine with existing pixels if they exist + ipix_depth_res = np.concatenate([new_ipix_d[depth_res][0], unique_ipix]) + vals_depth_res = np.concatenate([new_ipix_d[depth_res][1], vals_depth_res]) + ip_argsort = np.argsort(ipix_depth_res) + new_ipix_d[depth_res] = (ipix_depth_res[ip_argsort], vals_depth_res[ip_argsort]) + depth_ipix_d = new_ipix_d + return depth_ipix_d + + +def cull_from_pixel_map(depth_ipix_d: Dict[int, Tuple[np.ndarray, np.ndarray]], wcs, max_split_depth=7): + """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. Args: - projection (str): The map projection to use. Valid values include: - - moll - Molleweide projection (default) - - gnom - Gnomonic projection - - cart - Cartesian projection - - orth - Orthographic projection + 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: - The healpy method that plots a HEALPix map with the specified projection + A new map with the same datatype of depth_ipix_d, with backfacing cells split into higher order """ - if projection == "moll": - projection_method = hp.mollview - elif projection == "gnom": - projection_method = hp.gnomview - elif projection == "cart": - projection_method = hp.cartview - elif projection == "orth": - projection_method = hp.orthview - else: - raise NotImplementedError(f"unknown projection: {projection}") - return projection_method - - -def _plot_healpix_map(healpix_map, projection, title, cmap="viridis", **kwargs): - """Perform the plotting of a healpix pixel map. + depths = list(depth_ipix_d.keys()) + min_depth = min(depths) + max_depth = max(depths) + ipixels, vals = depth_ipix_d[min_depth] + + # Split the cells located at the border of the projection + # until at least the max_split_depth + max_split_depth = max(max_split_depth, max_depth) + + ipix_d = {} + for depth in range(min_depth, max_split_depth + 1): + # for each depth, check if pixels are too large, or wrap around projection, and split into pixels at + # higher order + if depth < 3: + too_large_ipix = ipixels + too_large_vals = vals + else: + ipix_lon, ipix_lat = cdshealpix.vertices(ipixels, depth) + + ipix_lon = ipix_lon[:, [2, 3, 0, 1]] + ipix_lat = ipix_lat[:, [2, 3, 0, 1]] + ipix_vertices = SkyCoord(ipix_lon, ipix_lat, frame=ICRS()) + + # Projection on the given WCS + xp, yp = skycoord_to_pixel(coords=ipix_vertices, wcs=wcs) + _, _, frontface_id = backface_culling(xp, yp) + + # Get the pixels which are backfacing the projection + backfacing_ipix = ipixels[~frontface_id] # backfacing + backfacing_vals = vals[~frontface_id] + frontface_ipix = ipixels[frontface_id] + frontface_vals = vals[frontface_id] + + ipix_d.update({depth: (frontface_ipix, frontface_vals)}) + + too_large_ipix = backfacing_ipix + too_large_vals = backfacing_vals + + next_depth = depth + 1 + + # get next depth if there is one, or use empty array as default + ipixels = np.array([], dtype=ipixels.dtype) + vals = np.array([], dtype=vals.dtype) + + if next_depth in depth_ipix_d: + ipixels, vals = depth_ipix_d[next_depth] + + # split too large ipix into next order, with each child getting the same map value as parent + + too_large_child_ipix = np.repeat(too_large_ipix << 2, 4) + np.tile( + np.array([0, 1, 2, 3]), len(too_large_ipix) + ) + too_large_child_vals = np.repeat(too_large_vals, 4) + + ipixels = np.concatenate((ipixels, too_large_child_ipix)) + vals = np.concatenate((vals, too_large_child_vals)) + + return ipix_d + + +def plot_healpix_map( + healpix_map: np.ndarray, + *, + projection: str = "MOL", + title: str = "", + cmap: str | Colormap = "viridis", + norm: Normalize | None = None, + ipix: np.ndarray | None = None, + depth: np.ndarray | None = None, + cbar: bool = True, + fov: Quantity | Tuple[Quantity, Quantity] = None, + center: SkyCoord | None = None, + wcs: astropy.wcs.WCS = None, + frame_class: Type[BaseFrame] | None = None, + ax: Axes | None = None, + fig: Figure | None = None, + **kwargs, +): + """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. Args: - healpix_map: array containing the map - projection: projection type to display - title: title used in image plot - cmap: matplotlib colormap to use + 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 (Axes | None): The matplotlib axes to plot onto. If None, an axes will be created to be used. If + specified, the axes must be initialized with a WCS for the projection, and passed to the method + with the WCS parameter. (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, Axes] - The figure and axes used to plot the healpix map """ - projection_method = get_projection_method(projection) + if ipix is None or depth is None: + order = int(np.ceil(np.log2(len(healpix_map) / 12) / 2)) + ipix = np.arange(len(healpix_map)) + depth = np.full(len(healpix_map), fill_value=order) + if fig is None: + if ax is not None: + fig = ax.get_figure() + else: + fig = plt.figure(figsize=(10, 5)) + if frame_class is None and fov is None and wcs is None: + frame_class = EllipticalFrame + if fov is None: + fov = (320 * u.deg, 160 * u.deg) + if center is None: + center = SkyCoord(0, 0, unit="deg", frame="icrs") + if ax is None: + if wcs is None: + wcs = WCS( + fig, + fov=fov, + center=center, + coordsys="icrs", + rotation=Angle(0, u.degree), + projection=projection, + ).w + ax = fig.add_subplot(1, 1, 1, projection=wcs, frame_class=frame_class) + elif wcs is None: + raise ValueError( + "if ax is provided, wcs must also be provided with the projection used in initializing ax" + ) + _plot_healpix_value_map(ipix, depth, healpix_map, ax, wcs, cmap=cmap, norm=norm, cbar=cbar, **kwargs) + plt.grid() + plt.ylabel("Dec") + plt.xlabel("RA") + plt.title(title) + return fig, ax + + +def _plot_healpix_value_map(ipix, depth, values, ax, wcs, cmap="viridis", norm=None, cbar=True, **kwargs): + """Perform the plotting of a healpix pixel map.""" + + # create dict mapping depth to ipix and values + depth_ipix_d = {} + + for d in np.unique(depth): + mask = depth == d + depth_ipix_d[d] = (ipix[mask], values[mask]) + + # cull backfacing and out of fov cells + fov_culled_d = cull_to_fov(depth_ipix_d, wcs) + culled_d = cull_from_pixel_map(fov_culled_d, wcs) + + # Generate Paths for each pixel and add to ax + plt_paths = [] + cum_vals = [] + for d, (ip, vals) in culled_d.items(): + vertices, codes = compute_healpix_vertices(depth=d, ipix=ip, wcs=wcs) + for i in range(len(ip)): + plt_paths.append(Path(vertices[5 * i : 5 * (i + 1)], codes[5 * i : 5 * (i + 1)])) + cum_vals.append(vals) + col = PathCollection(plt_paths, cmap=cmap, norm=norm, **kwargs) + col.set_array(np.concatenate(cum_vals)) + ax.add_collection(col) + + # Add color bar + if cbar: + plt.colorbar(col) - projection_method(healpix_map, title=title, nest=True, cmap=cmap, **kwargs) - plt.plot() + # Set projection + _set_wcs(ax, wcs) diff --git a/src/hats/pixel_tree/pixel_alignment.py b/src/hats/pixel_tree/pixel_alignment.py index 3c779547..60fd6b4e 100644 --- a/src/hats/pixel_tree/pixel_alignment.py +++ b/src/hats/pixel_tree/pixel_alignment.py @@ -409,7 +409,7 @@ def filter_alignment_by_moc(alignment: PixelAlignment, moc: MOC) -> PixelAlignme PixelAlignment object with the filtered mapping and tree """ moc_ranges = moc.to_depth29_ranges - tree_29_ranges = alignment.pixel_tree.tree << (2 * (29 - alignment.pixel_tree.tree_order)) + tree_29_ranges = alignment.pixel_tree.to_depth29_ranges() tree_mask = perform_filter_by_moc(tree_29_ranges, moc_ranges) new_tree = PixelTree(alignment.pixel_tree.tree[tree_mask], alignment.pixel_tree.tree_order) return PixelAlignment( diff --git a/src/hats/pixel_tree/pixel_tree.py b/src/hats/pixel_tree/pixel_tree.py index 71aaa3bf..5b4039f7 100644 --- a/src/hats/pixel_tree/pixel_tree.py +++ b/src/hats/pixel_tree/pixel_tree.py @@ -92,6 +92,10 @@ def to_moc(self) -> MOC: """Returns the MOC object that covers the same pixels as the tree""" return MOC.from_healpix_cells(self.pixels.T[1], self.pixels.T[0], self.tree_order) + def to_depth29_ranges(self) -> np.ndarray: + """Returns the ranges of the pixels in the tree at depth 29""" + return self.tree << (2 * (29 - self.tree_order)) + @classmethod def from_healpix(cls, healpix_pixels: Sequence[HealpixInputTypes], tree_order=None) -> PixelTree: """Build a tree from a list of constituent healpix pixels diff --git a/tests/hats/inspection/test_visualize_catalog.py b/tests/hats/inspection/test_visualize_catalog.py index 0cac97ea..a57f90f5 100644 --- a/tests/hats/inspection/test_visualize_catalog.py +++ b/tests/hats/inspection/test_visualize_catalog.py @@ -1,92 +1,671 @@ -import healpy as hp +import astropy.units as u +import matplotlib.pyplot as plt +import numpy as np import pytest +from astropy.coordinates import Angle, SkyCoord +from astropy.visualization.wcsaxes.frame import EllipticalFrame, RectangularFrame +from matplotlib.colors import LogNorm, Normalize +from matplotlib.pyplot import get_cmap +from mocpy import MOC, WCS +from mocpy.moc.plot.culling_backfacing_cells import from_moc +from mocpy.moc.plot.fill import compute_healpix_vertices +from mocpy.moc.plot.utils import build_plotting_moc -from hats.catalog import Catalog -from hats.inspection import plot_pixel_list, plot_pixels, plot_points -from hats.loaders import read_hats +from hats.inspection import plot_pixels +from hats.inspection.visualize_catalog import cull_from_pixel_map, cull_to_fov, plot_healpix_map # pylint: disable=no-member -@pytest.mark.timeout(10) -def test_generate_projections(small_sky_dir, mocker): - """Basic test that map data can be generated""" +DEFAULT_CMAP_NAME = "viridis" +DEFAULT_FOV = (320 * u.deg, 160 * u.deg) +DEFAULT_CENTER = SkyCoord(0, 0, unit="deg", frame="icrs") +DEFAULT_COORDSYS = "icrs" +DEFAULT_ROTATION = Angle(0, u.degree) +DEFAULT_PROJECTION = "MOL" - cat = read_hats(small_sky_dir) - projection = "moll" - mocker.patch("healpy.mollview") - plot_pixels(cat, projection=projection) - hp.mollview.assert_called_once() - hp.mollview.reset_mock() - plot_points(cat, projection=projection) - hp.mollview.assert_called_once() +def test_plot_healpix_pixels(): + order = 3 + length = 10 + ipix = np.arange(length) + pix_map = np.arange(length) + depth = np.full(length, fill_value=order) + fig, ax = plot_healpix_map(pix_map, ipix=ipix, depth=depth) + assert len(ax.collections) == 1 + col = ax.collections[0] + paths = col.get_paths() + assert len(paths) == length + assert col.get_cmap() == get_cmap(DEFAULT_CMAP_NAME) + assert isinstance(col.norm, Normalize) + assert col.norm.vmin == min(pix_map) + assert col.norm.vmax == max(pix_map) + assert col.colorbar is not None + assert col.colorbar.cmap == get_cmap(DEFAULT_CMAP_NAME) + assert col.colorbar.norm == col.norm + wcs = WCS( + fig, + fov=DEFAULT_FOV, + center=DEFAULT_CENTER, + coordsys=DEFAULT_COORDSYS, + rotation=DEFAULT_ROTATION, + projection=DEFAULT_PROJECTION, + ).w + for path, ipix in zip(paths, ipix): + verts, codes = compute_healpix_vertices(order, np.array([ipix]), wcs) + np.testing.assert_array_equal(path.vertices, verts) + np.testing.assert_array_equal(path.codes, codes) + np.testing.assert_array_equal(col.get_array(), pix_map) + assert ax.frame_class == EllipticalFrame - projection = "gnom" - mocker.patch("healpy.gnomview") - plot_pixels(cat, projection=projection) - hp.gnomview.assert_called_once() - hp.gnomview.reset_mock() - plot_points(cat, projection=projection) - hp.gnomview.assert_called_once() - projection = "cart" - mocker.patch("healpy.cartview") - plot_pixels(cat, projection=projection) - hp.cartview.assert_called_once() - hp.cartview.reset_mock() - plot_points(cat, projection=projection) - hp.cartview.assert_called_once() +def test_plot_healpix_pixels_different_order(): + order = 6 + length = 1000 + ipix = np.arange(length) + pix_map = np.arange(length) + depth = np.full(length, fill_value=order) + fig, ax = plot_healpix_map(pix_map, ipix=ipix, depth=depth) + assert len(ax.collections) == 1 + col = ax.collections[0] + paths = col.get_paths() + assert len(paths) == length + wcs = WCS( + fig, + fov=DEFAULT_FOV, + center=DEFAULT_CENTER, + coordsys=DEFAULT_COORDSYS, + rotation=DEFAULT_ROTATION, + projection=DEFAULT_PROJECTION, + ).w - projection = "orth" - mocker.patch("healpy.orthview") - plot_pixels(cat, projection=projection) - hp.orthview.assert_called_once() - hp.orthview.reset_mock() - plot_points(cat, projection=projection) - hp.orthview.assert_called_once() + all_verts, all_codes = compute_healpix_vertices(order, ipix, wcs) + for i, (path, ipix) in enumerate(zip(paths, ipix)): + verts, codes = all_verts[i * 5 : (i + 1) * 5], all_codes[i * 5 : (i + 1) * 5] + np.testing.assert_array_equal(path.vertices, verts) + np.testing.assert_array_equal(path.codes, codes) + np.testing.assert_array_equal(col.get_array(), pix_map) -def test_generate_map_unknown_projection(small_sky_dir): - """Test for error with unknown projection type""" +def test_order_0_pixels_split_to_order_3(): + map_value = 0.5 + order_0_pix = 4 + ipix = np.array([order_0_pix]) + pix_map = np.array([map_value]) + depth = np.array([0]) + fig, ax = plot_healpix_map(pix_map, ipix=ipix, depth=depth) + assert len(ax.collections) == 1 + col = ax.collections[0] + paths = col.get_paths() + length = 4**3 + order3_ipix = np.arange(length * order_0_pix, length * (order_0_pix + 1)) + assert len(paths) == length + wcs = WCS( + fig, + fov=DEFAULT_FOV, + center=DEFAULT_CENTER, + coordsys=DEFAULT_COORDSYS, + rotation=DEFAULT_ROTATION, + projection=DEFAULT_PROJECTION, + ).w + all_verts, all_codes = compute_healpix_vertices(3, order3_ipix, wcs) + for i, (path, ipix) in enumerate(zip(paths, order3_ipix)): + verts, codes = all_verts[i * 5 : (i + 1) * 5], all_codes[i * 5 : (i + 1) * 5] + np.testing.assert_array_equal(path.vertices, verts) + np.testing.assert_array_equal(path.codes, codes) + np.testing.assert_array_equal(col.get_array(), np.full(length, fill_value=map_value)) - cat = read_hats(small_sky_dir) - with pytest.raises(NotImplementedError): - plot_pixels(cat, projection=None) - with pytest.raises(NotImplementedError): - plot_pixels(cat, projection="") +def test_edge_pixels_split_to_order_7(): + map_value = 0.5 + order_0_pix = 2 + ipix = np.array([order_0_pix]) + pix_map = np.array([map_value]) + depth = np.array([0]) + fig, ax = plot_healpix_map(pix_map, ipix=ipix, depth=depth) + assert len(ax.collections) == 1 + edge_pixels = {0: [order_0_pix]} + for iter_ord in range(1, 8): + edge_pixels[iter_ord] = [p * 4 + i for p in edge_pixels[iter_ord - 1] for i in (2, 3)] + non_edge_pixels = {} + pixels_ord3 = np.arange(4**3 * order_0_pix, 4**3 * (order_0_pix + 1)) + non_edge_pixels[3] = pixels_ord3[~np.isin(pixels_ord3, edge_pixels[3])] + for iter_ord in range(4, 8): + pixels_ord = np.concatenate([np.arange(4 * pix, 4 * (pix + 1)) for pix in edge_pixels[iter_ord - 1]]) + non_edge_pixels[iter_ord] = pixels_ord[~np.isin(pixels_ord, edge_pixels[iter_ord])] + col = ax.collections[0] + paths = col.get_paths() + length = sum(len(x) for x in non_edge_pixels.values()) + assert len(paths) == length + wcs = WCS( + fig, + fov=DEFAULT_FOV, + center=DEFAULT_CENTER, + coordsys=DEFAULT_COORDSYS, + rotation=DEFAULT_ROTATION, + projection=DEFAULT_PROJECTION, + ).w + ords = np.concatenate([np.full(len(x), fill_value=o) for o, x in non_edge_pixels.items()]) + pixels = np.concatenate([np.array(x) for _, x in non_edge_pixels.items()]) + for path, iter_ord, pix in zip(paths, ords, pixels): + verts, codes = compute_healpix_vertices(iter_ord, np.array([pix]), wcs) + np.testing.assert_array_equal(path.vertices, verts) + np.testing.assert_array_equal(path.codes, codes) + np.testing.assert_array_equal(col.get_array(), np.full(length, fill_value=map_value)) - with pytest.raises(NotImplementedError): - plot_pixels(cat, projection="projection") +def test_cull_from_pixel_map(): + order = 1 + ipix = np.arange(12 * 4**order) + pix_map = np.arange(12 * 4**order) + map_dict = {order: (ipix, pix_map)} + fig = plt.figure(figsize=(10, 5)) + wcs = WCS( + fig, + fov=DEFAULT_FOV, + center=DEFAULT_CENTER, + coordsys=DEFAULT_COORDSYS, + rotation=DEFAULT_ROTATION, + projection=DEFAULT_PROJECTION, + ).w + culled_dict = cull_from_pixel_map(map_dict, wcs) + mocpy_culled = from_moc({str(order): ipix}, wcs) + for iter_ord, (pixels, m) in culled_dict.items(): + np.testing.assert_array_equal(pixels, mocpy_culled[str(iter_ord)]) + map_indices = pixels >> (2 * (iter_ord - order)) + np.testing.assert_array_equal(m, pix_map[map_indices]) -def test_plot_pixel_map_order1(small_sky_order1_dir, mocker): - """Basic test that map data can be generated""" - mocker.patch("healpy.mollview") - cat = read_hats(small_sky_order1_dir) - plot_pixels(cat) - hp.mollview.assert_called_once() +def test_cull_to_fov(): + order = 4 + ipix = np.arange(12 * 4**order) + pix_map = np.arange(12 * 4**order) + map_dict = {order: (ipix, pix_map)} + fig = plt.figure(figsize=(10, 5)) + wcs = WCS( + fig, + fov=(20 * u.deg, 10 * u.deg), + center=DEFAULT_CENTER, + coordsys=DEFAULT_COORDSYS, + rotation=DEFAULT_ROTATION, + projection=DEFAULT_PROJECTION, + ).w + culled_dict = cull_to_fov(map_dict, wcs) + moc = MOC.from_healpix_cells(ipix, np.full(ipix.shape, fill_value=order), max_depth=order) + mocpy_culled = build_plotting_moc(moc, wcs) + for iter_ord, (pixels, m) in culled_dict.items(): + for p in pixels: + assert ( + len( + MOC.from_healpix_cells(np.array([p]), np.array([iter_ord]), max_depth=iter_ord) + .intersection(mocpy_culled) + .to_depth29_ranges + ) + > 0 + ) + ord_ipix = np.arange(12 * 4**iter_ord) + ord_non_pixels = ord_ipix[~np.isin(ord_ipix, pixels)] + for p in ord_non_pixels: + assert ( + len( + MOC.from_healpix_cells(np.array([p]), np.array([iter_ord]), max_depth=iter_ord) + .intersection(mocpy_culled) + .to_depth29_ranges + ) + == 0 + ) + map_indices = pixels >> (2 * (iter_ord - order)) + np.testing.assert_array_equal(m, pix_map[map_indices]) - hp.mollview.reset_mock() - plot_points(cat) - hp.mollview.assert_called_once() +def test_cull_to_fov_subsamples_high_order(): + order = 10 + ipix = np.arange(12 * 4**order) + pix_map = np.arange(12 * 4**order) + map_dict = {order: (ipix, pix_map)} + fig = plt.figure(figsize=(10, 5)) + wcs = WCS( + fig, + fov=DEFAULT_FOV, + center=DEFAULT_CENTER, + coordsys=DEFAULT_COORDSYS, + rotation=DEFAULT_ROTATION, + projection=DEFAULT_PROJECTION, + ).w + with pytest.warns(match="smaller"): + culled_dict = cull_to_fov(map_dict, wcs) + # Get the WCS cdelt giving the deg.px^(-1) resolution. + cdelt = wcs.wcs.cdelt + # Convert in rad.px^(-1) + cdelt = np.abs((2 * np.pi / 360) * cdelt[0]) + # Get the minimum depth such as the resolution of a cell is contained in 1px. + depth_res = int(np.floor(np.log2(np.sqrt(np.pi / 3) / cdelt))) + depth_res = max(depth_res, 0) + assert depth_res < order + for iter_ord, (pixels, m) in culled_dict.items(): + assert iter_ord == depth_res + assert np.all(np.isin(ipix >> (2 * (order - depth_res)), pixels)) + map_indices = pixels << (2 * (order - depth_res)) + np.testing.assert_array_equal(m, pix_map[map_indices]) -def test_visualize_in_memory_catalogs(catalog_info, catalog_pixels, mocker): - """Test behavior of visualization methods for non-on-disk catalogs and pixel data.""" - mocker.patch("healpy.mollview") - catalog = Catalog(catalog_info, catalog_pixels) - plot_pixels(catalog) - hp.mollview.assert_called_once() - hp.mollview.reset_mock() - plot_pixel_list(catalog_pixels, plot_title="My special catalog") - hp.mollview.assert_called_once() +def test_cull_to_fov_subsamples_multiple_orders(): + depth = np.array([0, 5, 8, 10]) + ipix = np.array([10, 5, 4, 2]) + pix_map = np.array([1, 2, 3, 4]) + map_dict = {depth[i]: (ipix[[i]], pix_map[[i]]) for i in range(len(depth))} + fig = plt.figure(figsize=(10, 5)) + wcs = WCS( + fig, + fov=DEFAULT_FOV, + center=DEFAULT_CENTER, + coordsys=DEFAULT_COORDSYS, + rotation=DEFAULT_ROTATION, + projection=DEFAULT_PROJECTION, + ).w + with pytest.warns(match="smaller"): + culled_dict = cull_to_fov(map_dict, wcs) + # Get the WCS cdelt giving the deg.px^(-1) resolution. + cdelt = wcs.wcs.cdelt + # Convert in rad.px^(-1) + cdelt = np.abs((2 * np.pi / 360) * cdelt[0]) + # Get the minimum depth such as the resolution of a cell is contained in 1px. + depth_res = int(np.floor(np.log2(np.sqrt(np.pi / 3) / cdelt))) + depth_res = max(depth_res, 0) + assert depth_res < max(depth) - hp.mollview.reset_mock() - with pytest.raises(ValueError, match="on disk catalog required"): - plot_points(catalog) - hp.mollview.assert_not_called() + assert list(culled_dict.keys()) == [0, 5, depth_res] + + assert culled_dict[0] == (np.array([10]), np.array([1])) + assert culled_dict[5] == (np.array([5]), np.array([2])) + small_pixels_map = pix_map[2:] + small_pixels_converted = ipix[2:] >> (2 * (depth[2:] - depth_res)) + small_pixels_argsort = np.argsort(small_pixels_converted) + assert np.all(culled_dict[depth_res][0] == small_pixels_converted[small_pixels_argsort]) + assert np.all(culled_dict[depth_res][1] == small_pixels_map[small_pixels_argsort]) + + +def test_plot_healpix_map(): + order = 1 + ipix = np.arange(12 * 4**order) + pix_map = np.arange(12 * 4**order) + fig, ax = plot_healpix_map(pix_map) + assert len(ax.collections) == 1 + col = ax.collections[0] + paths = col.get_paths() + wcs = WCS( + fig, + fov=DEFAULT_FOV, + center=DEFAULT_CENTER, + coordsys=DEFAULT_COORDSYS, + rotation=DEFAULT_ROTATION, + projection=DEFAULT_PROJECTION, + ).w + map_dict = {order: (ipix, pix_map)} + culled_dict = cull_from_pixel_map(map_dict, wcs) + all_vals = [] + start_i = 0 + for iter_ord, (pixels, pix_map) in culled_dict.items(): + all_verts, all_codes = compute_healpix_vertices(iter_ord, pixels, wcs) + for i, _ in enumerate(pixels): + verts, codes = all_verts[i * 5 : (i + 1) * 5], all_codes[i * 5 : (i + 1) * 5] + path = paths[start_i + i] + np.testing.assert_array_equal(path.vertices, verts) + np.testing.assert_array_equal(path.codes, codes) + all_vals.append(pix_map) + start_i += len(pixels) + assert start_i == len(paths) + np.testing.assert_array_equal(np.concatenate(all_vals), col.get_array()) + + +def test_plot_wcs_params(): + order = 3 + length = 10 + ipix = np.arange(length) + pix_map = np.arange(length) + depth = np.full(length, fill_value=order) + fig, ax = plot_healpix_map( + pix_map, + ipix=ipix, + depth=depth, + fov=(100 * u.deg, 50 * u.deg), + center=SkyCoord(10, 10, unit="deg", frame="icrs"), + projection="AIT", + ) + assert len(ax.collections) == 1 + col = ax.collections[0] + paths = col.get_paths() + assert len(paths) == length + wcs = WCS( + fig, + fov=(100 * u.deg, 50 * u.deg), + center=SkyCoord(10, 10, unit="deg", frame="icrs"), + coordsys=DEFAULT_COORDSYS, + rotation=DEFAULT_ROTATION, + projection="AIT", + ).w + for path, ipix in zip(paths, np.arange(len(pix_map))): + verts, codes = compute_healpix_vertices(order, np.array([ipix]), wcs) + np.testing.assert_array_equal(path.vertices, verts) + np.testing.assert_array_equal(path.codes, codes) + np.testing.assert_array_equal(col.get_array(), pix_map) + assert ax.get_transform("icrs") is not None + assert ax.frame_class == RectangularFrame + + +def test_plot_wcs_params_frame(): + order = 3 + length = 10 + ipix = np.arange(length) + pix_map = np.arange(length) + depth = np.full(length, fill_value=order) + fig, ax = plot_healpix_map( + pix_map, + ipix=ipix, + depth=depth, + fov=(100 * u.deg, 50 * u.deg), + center=SkyCoord(10, 10, unit="deg", frame="icrs"), + projection="AIT", + frame_class=EllipticalFrame, + ) + assert len(ax.collections) == 1 + col = ax.collections[0] + paths = col.get_paths() + assert len(paths) == length + wcs = WCS( + fig, + fov=(100 * u.deg, 50 * u.deg), + center=SkyCoord(10, 10, unit="deg", frame="icrs"), + coordsys=DEFAULT_COORDSYS, + rotation=DEFAULT_ROTATION, + projection="AIT", + ).w + for path, ipix in zip(paths, np.arange(len(pix_map))): + verts, codes = compute_healpix_vertices(order, np.array([ipix]), wcs) + np.testing.assert_array_equal(path.vertices, verts) + np.testing.assert_array_equal(path.codes, codes) + np.testing.assert_array_equal(col.get_array(), pix_map) + assert ax.get_transform("icrs") is not None + assert ax.frame_class == EllipticalFrame + + +def test_plot_fov_culling(): + order = 3 + length = 10 + ipix = np.arange(length) + pix_map = np.arange(length) + depth = np.full(length, fill_value=order) + fig, ax = plot_healpix_map( + pix_map, + ipix=ipix, + depth=depth, + fov=(30 * u.deg, 20 * u.deg), + center=SkyCoord(10, 10, unit="deg", frame="icrs"), + projection="AIT", + ) + assert len(ax.collections) == 1 + col = ax.collections[0] + paths = col.get_paths() + wcs = WCS( + fig, + fov=(30 * u.deg, 20 * u.deg), + center=SkyCoord(10, 10, unit="deg", frame="icrs"), + coordsys=DEFAULT_COORDSYS, + rotation=DEFAULT_ROTATION, + projection="AIT", + ).w + map_dict = {order: (ipix, pix_map)} + culled_dict = cull_to_fov(map_dict, wcs) + all_vals = [] + start_i = 0 + for iter_ord, (pixels, pix_map) in culled_dict.items(): + all_verts, all_codes = compute_healpix_vertices(iter_ord, pixels, wcs) + for i, _ in enumerate(pixels): + verts, codes = all_verts[i * 5 : (i + 1) * 5], all_codes[i * 5 : (i + 1) * 5] + path = paths[start_i + i] + np.testing.assert_array_equal(path.vertices, verts) + np.testing.assert_array_equal(path.codes, codes) + all_vals.append(pix_map) + start_i += len(pixels) + assert start_i == len(paths) + np.testing.assert_array_equal(np.concatenate(all_vals), col.get_array()) + assert ax.get_transform("icrs") is not None + assert ax.frame_class == RectangularFrame + + +def test_plot_wcs(): + order = 3 + length = 10 + ipix = np.arange(length) + pix_map = np.arange(length) + depth = np.full(length, fill_value=order) + fig = plt.figure(figsize=(10, 5)) + wcs = WCS( + fig, + fov=(100 * u.deg, 50 * u.deg), + center=SkyCoord(10, 10, unit="deg", frame="icrs"), + coordsys=DEFAULT_COORDSYS, + rotation=DEFAULT_ROTATION, + projection="AIT", + ).w + fig2, ax = plot_healpix_map(pix_map, ipix=ipix, depth=depth, fig=fig, wcs=wcs) + assert fig2 is fig + assert len(ax.collections) == 1 + col = ax.collections[0] + paths = col.get_paths() + assert len(paths) == length + for path, ipix in zip(paths, np.arange(len(pix_map))): + verts, codes = compute_healpix_vertices(order, np.array([ipix]), wcs) + np.testing.assert_array_equal(path.vertices, verts) + np.testing.assert_array_equal(path.codes, codes) + np.testing.assert_array_equal(col.get_array(), pix_map) + assert ax.get_transform("icrs") is not None + + +def test_plot_wcs_and_ax(): + order = 3 + length = 10 + ipix = np.arange(length) + pix_map = np.arange(length) + depth = np.full(length, fill_value=order) + fig = plt.figure(figsize=(10, 5)) + wcs = WCS( + fig, + fov=(100 * u.deg, 50 * u.deg), + center=SkyCoord(10, 10, unit="deg", frame="icrs"), + coordsys=DEFAULT_COORDSYS, + rotation=DEFAULT_ROTATION, + projection="AIT", + ).w + ax = fig.add_subplot(1, 1, 1, projection=wcs, frame_class=EllipticalFrame) + assert len(ax.collections) == 0 + fig2, ax2 = plot_healpix_map(pix_map, ipix=ipix, depth=depth, fig=fig, wcs=wcs, ax=ax) + assert fig2 is fig + assert ax2 is ax + assert len(ax.collections) == 1 + col = ax.collections[0] + paths = col.get_paths() + assert len(paths) == length + for path, ipix in zip(paths, np.arange(len(pix_map))): + verts, codes = compute_healpix_vertices(order, np.array([ipix]), wcs) + np.testing.assert_array_equal(path.vertices, verts) + np.testing.assert_array_equal(path.codes, codes) + np.testing.assert_array_equal(col.get_array(), pix_map) + assert ax.get_transform("icrs") is not None + + +def test_plot_ax_no_wcs(): + order = 3 + length = 10 + ipix = np.arange(length) + pix_map = np.arange(length) + depth = np.full(length, fill_value=order) + fig = plt.figure(figsize=(10, 5)) + wcs = WCS( + fig, + fov=(100 * u.deg, 50 * u.deg), + center=SkyCoord(10, 10, unit="deg", frame="icrs"), + coordsys=DEFAULT_COORDSYS, + rotation=DEFAULT_ROTATION, + projection="AIT", + ).w + ax = fig.add_subplot(1, 1, 1, projection=wcs, frame_class=EllipticalFrame) + with pytest.raises(ValueError): + plot_healpix_map(pix_map, ipix=ipix, depth=depth, fig=fig, ax=ax) + + +def test_plot_cmaps(): + order = 3 + length = 10 + cmap_name = "plasma" + ipix = np.arange(length) + pix_map = np.arange(length) + depth = np.full(length, fill_value=order) + fig, ax = plot_healpix_map(pix_map, ipix=ipix, depth=depth, cmap=cmap_name) + assert len(ax.collections) == 1 + col = ax.collections[0] + paths = col.get_paths() + assert col.get_cmap() == get_cmap(cmap_name) + assert col.colorbar is not None + assert col.colorbar.cmap == get_cmap(cmap_name) + assert len(paths) == length + wcs = WCS( + fig, + fov=DEFAULT_FOV, + center=DEFAULT_CENTER, + coordsys=DEFAULT_COORDSYS, + rotation=DEFAULT_ROTATION, + projection=DEFAULT_PROJECTION, + ).w + for path, ipix in zip(paths, np.arange(len(pix_map))): + verts, codes = compute_healpix_vertices(order, np.array([ipix]), wcs) + np.testing.assert_array_equal(path.vertices, verts) + np.testing.assert_array_equal(path.codes, codes) + np.testing.assert_array_equal(col.get_array(), pix_map) + + ipix = np.arange(length) + pix_map = np.arange(length) + depth = np.full(length, fill_value=order) + fig, ax = plot_healpix_map(pix_map, ipix=ipix, depth=depth, cmap=get_cmap(cmap_name)) + assert len(ax.collections) == 1 + col = ax.collections[0] + paths = col.get_paths() + assert col.get_cmap() == get_cmap(cmap_name) + assert col.colorbar is not None + assert col.colorbar.cmap == get_cmap(cmap_name) + assert len(paths) == length + for path, ipix in zip(paths, np.arange(len(pix_map))): + verts, codes = compute_healpix_vertices(order, np.array([ipix]), wcs) + np.testing.assert_array_equal(path.vertices, verts) + np.testing.assert_array_equal(path.codes, codes) + np.testing.assert_array_equal(col.get_array(), pix_map) + + +def test_plot_norm(): + order = 3 + length = 10 + norm = LogNorm() + ipix = np.arange(length) + pix_map = np.arange(length) + depth = np.full(length, fill_value=order) + fig, ax = plot_healpix_map(pix_map, ipix=ipix, depth=depth, norm=norm) + assert len(ax.collections) == 1 + col = ax.collections[0] + paths = col.get_paths() + assert col.norm == norm + assert col.colorbar is not None + assert col.colorbar.norm == norm + assert len(paths) == length + wcs = WCS( + fig, + fov=DEFAULT_FOV, + center=DEFAULT_CENTER, + coordsys=DEFAULT_COORDSYS, + rotation=DEFAULT_ROTATION, + projection=DEFAULT_PROJECTION, + ).w + for path, ipix in zip(paths, np.arange(len(pix_map))): + verts, codes = compute_healpix_vertices(order, np.array([ipix]), wcs) + np.testing.assert_array_equal(path.vertices, verts) + np.testing.assert_array_equal(path.codes, codes) + np.testing.assert_array_equal(col.get_array(), pix_map) + + +def test_plot_no_cbar(): + order = 3 + length = 10 + ipix = np.arange(length) + pix_map = np.arange(length) + depth = np.full(length, fill_value=order) + fig, ax = plot_healpix_map(pix_map, ipix=ipix, depth=depth, cbar=False) + assert len(ax.collections) == 1 + col = ax.collections[0] + assert col.colorbar is None + paths = col.get_paths() + assert len(paths) == length + wcs = WCS( + fig, + fov=DEFAULT_FOV, + center=DEFAULT_CENTER, + coordsys=DEFAULT_COORDSYS, + rotation=DEFAULT_ROTATION, + projection=DEFAULT_PROJECTION, + ).w + for path, ipix in zip(paths, np.arange(len(pix_map))): + verts, codes = compute_healpix_vertices(order, np.array([ipix]), wcs) + np.testing.assert_array_equal(path.vertices, verts) + np.testing.assert_array_equal(path.codes, codes) + np.testing.assert_array_equal(col.get_array(), pix_map) + + +def test_plot_kwargs(): + order = 3 + length = 10 + label = "test" + ipix = np.arange(length) + pix_map = np.arange(length) + depth = np.full(length, fill_value=order) + fig, ax = plot_healpix_map(pix_map, ipix=ipix, depth=depth, label=label) + assert len(ax.collections) == 1 + col = ax.collections[0] + assert col.get_label() == label + paths = col.get_paths() + assert len(paths) == length + wcs = WCS( + fig, + fov=DEFAULT_FOV, + center=DEFAULT_CENTER, + coordsys=DEFAULT_COORDSYS, + rotation=DEFAULT_ROTATION, + projection=DEFAULT_PROJECTION, + ).w + for path, ipix in zip(paths, np.arange(len(pix_map))): + verts, codes = compute_healpix_vertices(order, np.array([ipix]), wcs) + np.testing.assert_array_equal(path.vertices, verts) + np.testing.assert_array_equal(path.codes, codes) + np.testing.assert_array_equal(col.get_array(), pix_map) + + +def test_catalog_plot(small_sky_order1_catalog): + fig, ax = plot_pixels(small_sky_order1_catalog) + pixels = sorted(small_sky_order1_catalog.get_healpix_pixels()) + order_3_pixels = [p for pix in pixels for p in pix.convert_to_higher_order(3 - pix.order)] + order_3_orders = [pix.order for pix in pixels for _ in pix.convert_to_higher_order(3 - pix.order)] + col = ax.collections[0] + paths = col.get_paths() + assert len(paths) == len(order_3_pixels) + wcs = WCS( + fig, + fov=DEFAULT_FOV, + center=DEFAULT_CENTER, + coordsys=DEFAULT_COORDSYS, + rotation=DEFAULT_ROTATION, + projection=DEFAULT_PROJECTION, + ).w + for p, path in zip(order_3_pixels, paths): + verts, codes = compute_healpix_vertices(p.order, np.array([p.pixel]), wcs) + np.testing.assert_array_equal(path.vertices, verts) + np.testing.assert_array_equal(path.codes, codes) + np.testing.assert_array_equal(col.get_array(), np.array(order_3_orders)) + assert ax.get_title() == f"Catalog pixel density map - {small_sky_order1_catalog.catalog_name}" diff --git a/tests/hats/pixel_tree/test_pixel_tree.py b/tests/hats/pixel_tree/test_pixel_tree.py index 0531f1cd..621b4ced 100644 --- a/tests/hats/pixel_tree/test_pixel_tree.py +++ b/tests/hats/pixel_tree/test_pixel_tree.py @@ -43,6 +43,12 @@ def test_pixel_tree_to_moc(pixel_tree_2): assert np.all(moc.flatten() == moc_order_pixels) +def test_pixel_tree_to_depth29_ranges(pixel_tree_2): + np.testing.assert_array_equal( + pixel_tree_2.to_depth29_ranges(), pixel_tree_2.tree << (2 * (29 - pixel_tree_2.tree_order)) + ) + + def test_pixel_tree_contains(): tree = PixelTree.from_healpix([(0, 0)]) assert tree.contains((0, 0))