22"""Tasks for making and manipulating HIPS images."""
24__all__ = [
"HighResolutionHipsTask",
"HighResolutionHipsConfig",
"HighResolutionHipsConnections",
25 "HighResolutionHipsQuantumGraphBuilder",
26 "GenerateHipsTask",
"GenerateHipsConfig",
"GenerateColorHipsTask",
"GenerateColorHipsConfig"]
28from collections
import defaultdict
36from datetime
import datetime
38import healsparse
as hsp
39from astropy.io
import fits
40from astropy.visualization.lupton_rgb
import AsinhMapping
44from lsst.utils.timer
import timeMethod
45from lsst.daf.butler
import Butler
48from lsst.pipe.base.quantum_graph_builder
import QuantumGraphBuilder
49from lsst.pipe.base.quantum_graph_skeleton
import QuantumGraphSkeleton, DatasetKey
55from lsst.resources
import ResourcePath
57from .healSparseMapping
import _is_power_of_two
61 dimensions=(
"healpix9",
"band"),
62 defaultTemplates={
"coaddName":
"deep"}):
63 coadd_exposure_handles = pipeBase.connectionTypes.Input(
64 doc=
"Coadded exposures to convert to HIPS format.",
65 name=
"{coaddName}Coadd_calexp",
66 storageClass=
"ExposureF",
67 dimensions=(
"tract",
"patch",
"skymap",
"band"),
71 hips_exposures = pipeBase.connectionTypes.Output(
72 doc=
"HiPS-compatible HPX image.",
73 name=
"{coaddName}Coadd_hpx",
74 storageClass=
"ExposureF",
75 dimensions=(
"healpix11",
"band"),
79 def __init__(self, *, config=None):
80 super().__init__(config=config)
83 for dim
in self.dimensions:
85 if quantum_order
is not None:
86 raise ValueError(
"Must not specify more than one quantum healpix dimension.")
87 quantum_order = int(dim.split(
"healpix")[1])
88 if quantum_order
is None:
89 raise ValueError(
"Must specify a healpix dimension in quantum dimensions.")
91 if quantum_order > config.hips_order:
92 raise ValueError(
"Quantum healpix dimension order must not be greater than hips_order")
95 for dim
in self.hips_exposures.dimensions:
98 raise ValueError(
"Must not specify more than one healpix dimension.")
99 order = int(dim.split(
"healpix")[1])
101 raise ValueError(
"Must specify a healpix dimension in hips_exposure dimensions.")
103 if order != config.hips_order:
104 raise ValueError(
"healpix dimension order must match config.hips_order.")
107class HighResolutionHipsConfig(pipeBase.PipelineTaskConfig,
108 pipelineConnections=HighResolutionHipsConnections):
109 """Configuration parameters for HighResolutionHipsTask.
113 A HiPS image covers one HEALPix cell, with the HEALPix nside equal to
114 2**hips_order. Each cell is 'shift_order' orders deeper than the HEALPix
115 cell, with 2**shift_order x 2**shift_order sub-pixels on a side, which
116 defines the target resolution of the HiPS image. The IVOA recommends
117 shift_order=9, for 2**9=512 pixels on a side.
120 https://www.ivoa.net/documents/HiPS/20170519/REC-HIPS-1.0-20170519.pdf
121 shows the relationship between hips_order, number of tiles (full
122 sky coverage), cell size, and sub-pixel size/image resolution (with
123 the default shift_order=9):
124 +------------+-----------------+--------------+------------------+
125 | hips_order | Number of Tiles | Cell Size | Image Resolution |
126 +============+=================+==============+==================+
127 | 0 | 12 | 58.63 deg | 6.871 arcmin |
128 | 1 | 48 | 29.32 deg | 3.435 arcmin |
129 | 2 | 192 | 14.66 deg | 1.718 arcmin |
130 | 3 | 768 | 7.329 deg | 51.53 arcsec |
131 | 4 | 3072 | 3.665 deg | 25.77 arcsec |
132 | 5 | 12288 | 1.832 deg | 12.88 arcsec |
133 | 6 | 49152 | 54.97 arcmin | 6.442 arcsec |
134 | 7 | 196608 | 27.48 arcmin | 3.221 arcsec |
135 | 8 | 786432 | 13.74 arcmin | 1.61 arcsec |
136 | 9 | 3145728 | 6.871 arcmin | 805.2mas |
137 | 10 | 12582912 | 3.435 arcmin | 402.6mas |
138 | 11 | 50331648 | 1.718 arcmin | 201.3mas |
139 | 12 | 201326592 | 51.53 arcsec | 100.6mas |
140 | 13 | 805306368 | 25.77 arcsec | 50.32mas |
141 +------------+-----------------+--------------+------------------+
143 hips_order = pexConfig.Field(
144 doc=
"HIPS image order.",
148 shift_order = pexConfig.Field(
149 doc=
"HIPS shift order (such that each tile is 2**shift_order pixels on a side)",
153 warp = pexConfig.ConfigField(
154 dtype=afwMath.Warper.ConfigClass,
155 doc=
"Warper configuration",
158 def setDefaults(self):
159 self.warp.warpingKernelName =
"lanczos5"
162class HipsTaskNameDescriptor:
163 """Descriptor used create a DefaultName that matches the order of
164 the defined dimensions in the connections class.
169 The prefix of the Default name, to which the order will be
172 def __init__(self, prefix):
174 self._defaultName = f
"{prefix}{{}}"
177 def __get__(self, obj, klass=None):
180 "HipsTaskDescriptor was used in an unexpected context"
182 if self._order
is None:
183 klassDimensions = klass.ConfigClass.ConnectionsClass.dimensions
184 for dim
in klassDimensions:
185 if (match := re.match(
r"^healpix(\d*)$", dim))
is not None:
186 self._order = int(match.group(1))
190 "Could not find healpix dimension in connections class"
192 return self._defaultName.format(self._order)
195class HighResolutionHipsTask(pipeBase.PipelineTask):
196 """Task for making high resolution HiPS images."""
197 ConfigClass = HighResolutionHipsConfig
198 _DefaultName = HipsTaskNameDescriptor(
"highResolutionHips")
200 def __init__(self, **kwargs):
201 super().__init__(**kwargs)
202 self.warper = afwMath.Warper.fromConfig(self.config.warp)
205 def runQuantum(self, butlerQC, inputRefs, outputRefs):
206 inputs = butlerQC.get(inputRefs)
208 healpix_dim = f
"healpix{self.config.hips_order}"
210 pixels = [hips_exposure.dataId[healpix_dim]
211 for hips_exposure
in outputRefs.hips_exposures]
213 outputs = self.run(pixels=pixels, coadd_exposure_handles=inputs[
"coadd_exposure_handles"])
215 hips_exposure_ref_dict = {hips_exposure_ref.dataId[healpix_dim]:
216 hips_exposure_ref
for hips_exposure_ref
in outputRefs.hips_exposures}
217 for pixel, hips_exposure
in outputs.hips_exposures.items():
218 butlerQC.put(hips_exposure, hips_exposure_ref_dict[pixel])
220 def run(self, pixels, coadd_exposure_handles):
221 """Run the HighResolutionHipsTask.
225 pixels : `Iterable` [ `int` ]
226 Iterable of healpix pixels (nest ordering) to warp to.
227 coadd_exposure_handles : `list` [`lsst.daf.butler.DeferredDatasetHandle`]
228 Handles for the coadd exposures.
232 outputs : `lsst.pipe.base.Struct`
233 ``hips_exposures`` is a dict with pixel (key) and hips_exposure (value)
235 self.log.info(
"Generating HPX images for %d pixels at order %d", len(pixels), self.config.hips_order)
237 npix = 2**self.config.shift_order
247 wcs_hpx = afwGeom.makeHpxWcs(self.config.hips_order, pixel, shift_order=self.config.shift_order)
248 exp_hpx = afwImage.ExposureF(bbox_hpx, wcs_hpx)
249 exp_hpx_dict[pixel] = exp_hpx
250 warp_dict[pixel] = []
255 for handle
in coadd_exposure_handles:
256 coadd_exp = handle.get()
260 warped = self.warper.warpExposure(exp_hpx_dict[pixel].getWcs(), coadd_exp, maxBBox=bbox_hpx)
262 exp = afwImage.ExposureF(exp_hpx_dict[pixel].getBBox(), exp_hpx_dict[pixel].getWcs())
263 exp.maskedImage.set(np.nan, afwImage.Mask.getPlaneBitMask(
"NO_DATA"), np.nan)
268 exp_hpx_dict[pixel].mask.conformMaskPlanes(coadd_exp.mask.getMaskPlaneDict())
269 exp_hpx_dict[pixel].setFilter(coadd_exp.getFilter())
270 exp_hpx_dict[pixel].setPhotoCalib(coadd_exp.getPhotoCalib())
272 if warped.getBBox().getArea() == 0
or not np.any(np.isfinite(warped.image.array)):
275 "No overlap between output HPX %d and input exposure %s",
281 exp.maskedImage.assign(warped.maskedImage, warped.getBBox())
282 warp_dict[pixel].append(exp.maskedImage)
288 stats_ctrl.setNanSafe(
True)
289 stats_ctrl.setWeighted(
True)
290 stats_ctrl.setCalcErrorFromInputVariance(
True)
296 exp_hpx_dict[pixel].maskedImage.set(np.nan, afwImage.Mask.getPlaneBitMask(
"NO_DATA"), np.nan)
298 if not warp_dict[pixel]:
300 self.log.debug(
"No data in HPX pixel %d", pixel)
303 exp_hpx_dict.pop(pixel)
310 [1.0]*len(warp_dict[pixel]),
315 return pipeBase.Struct(hips_exposures=exp_hpx_dict)
318 def build_quantum_graph_cli(cls, argv):
319 """A command-line interface entry point to `build_quantum_graph`.
320 This method provides the implementation for the
321 ``build-high-resolution-hips-qg`` script.
325 argv : `Sequence` [ `str` ]
326 Command-line arguments (e.g. ``sys.argv[1:]``).
328 parser = cls._make_cli_parser()
330 args = parser.parse_args(argv)
332 if args.subparser_name
is None:
336 pipeline = pipeBase.Pipeline.from_uri(args.pipeline)
337 pipeline_graph = pipeline.to_graph()
339 if len(pipeline_graph.tasks) != 1:
340 raise RuntimeError(f
"Pipeline file {args.pipeline} may only contain one task.")
342 (task_node,) = pipeline_graph.tasks.values()
344 butler = Butler(args.butler_config, collections=args.input)
346 if args.subparser_name ==
"segment":
349 dataset = task_node.inputs[
"coadd_exposure_handles"].dataset_type_name
350 data_ids =
set(butler.registry.queryDataIds(
"tract", datasets=dataset).expanded())
352 for data_id
in data_ids:
353 region = data_id.region
354 pixel_range = hpix_pixelization.envelope(region)
355 for r
in pixel_range.ranges():
356 region_pixels.extend(range(r[0], r[1]))
357 indices = np.unique(region_pixels)
359 print(f
"Pixels to run at HEALPix order --hpix_build_order {args.hpix_build_order}:")
360 for pixel
in indices:
363 elif args.subparser_name ==
"build":
367 if args.output_run
is None:
368 if args.output
is None:
369 raise ValueError(
"At least one of --output or --output-run options is required.")
370 args.output_run =
"{}/{}".format(args.output, pipeBase.Instrument.makeCollectionTimestamp())
372 build_ranges =
RangeSet(sorted(args.pixels))
377 "butler_argument": args.butler_config,
378 "output": args.output,
379 "output_run": args.output_run,
380 "data_query": args.where,
381 "time": f
"{datetime.now()}",
384 builder = HighResolutionHipsQuantumGraphBuilder(
387 input_collections=args.input,
388 output_run=args.output_run,
389 constraint_order=args.hpix_build_order,
390 constraint_ranges=build_ranges,
393 qg = builder.build(metadata, attach_datastore_records=
True)
394 qg.saveUri(args.save_qgraph)
397 def _make_cli_parser(cls):
398 """Make the command-line parser.
402 parser : `argparse.ArgumentParser`
404 parser = argparse.ArgumentParser(
406 "Build a QuantumGraph that runs HighResolutionHipsTask on existing coadd datasets."
409 subparsers = parser.add_subparsers(help=
"sub-command help", dest=
"subparser_name")
411 parser_segment = subparsers.add_parser(
"segment",
412 help=
"Determine survey segments for workflow.")
413 parser_build = subparsers.add_parser(
"build",
414 help=
"Build quantum graph for HighResolutionHipsTask")
416 for sub
in [parser_segment, parser_build]:
422 help=
"Path to data repository or butler configuration.",
429 help=
"Pipeline file, limited to one task.",
437 help=
"Input collection(s) to search for coadd exposures.",
442 "--hpix_build_order",
445 help=
"HEALPix order to segment sky for building quantum graph files.",
452 help=
"Data ID expression used when querying for input coadd datasets.",
455 parser_build.add_argument(
459 "Name of the output CHAINED collection. If this options is specified and "
460 "--output-run is not, then a new RUN collection will be created by appending "
461 "a timestamp to the value of this option."
466 parser_build.add_argument(
470 "Output RUN collection to write resulting images. If not provided "
471 "then --output must be provided and a new RUN collection will be created "
472 "by appending a timestamp to the value passed with --output."
477 parser_build.add_argument(
481 help=
"Output filename for QuantumGraph.",
484 parser_build.add_argument(
489 help=
"Pixels at --hpix_build_order to generate quantum graph.",
496class HighResolutionHipsQuantumGraphBuilder(QuantumGraphBuilder):
497 """A custom a `lsst.pipe.base.QuantumGraphBuilder` for running
498 `HighResolutionHipsTask` only.
500 This is a workaround for incomplete butler query support for HEALPix
505 pipeline_graph : `lsst.pipe.base.PipelineGraph`
506 Pipeline graph with exactly one task, which must be a configuration
507 of `HighResolutionHipsTask`.
508 butler : `lsst.daf.butler.Butler`
509 Client for the butler data repository. May be read-only.
510 input_collections : `str` or `Iterable` [ `str` ], optional
511 Collection or collections to search for input datasets, in order.
512 If not provided, ``butler.collections`` will be searched.
513 output_run : `str`, optional
514 Name of the output collection. If not provided, ``butler.run`` will
516 constraint_order : `int`
517 HEALPix order used to constrain which quanta are generated, via
518 ``constraint_indices``. This should be a coarser grid (smaller
519 order) than the order used for the task's quantum and output data
520 IDs, and ideally something between the spatial scale of a patch or
521 the data repository's "common skypix" system (usually ``htm7``).
522 constraint_ranges : `lsst.sphgeom.RangeSet`
523 RangeSet that describes constraint pixels (HEALPix NEST, with order
524 ``constraint_order``) to constrain generated quanta.
525 where : `str`, optional
526 A boolean `str` expression of the form accepted by
527 `Registry.queryDatasets` to constrain input datasets. This may
528 contain a constraint on tracts, patches, or bands, but not HEALPix
529 indices. Constraints on tracts and patches should usually be
530 unnecessary, however - existing coadds that overlap the given
531 HEALpix indices will be selected without such a constraint, and
532 providing one may reject some that should normally be included.
540 input_collections=None,
546 super().__init__(pipeline_graph, butler, input_collections=input_collections, output_run=output_run)
547 self.constraint_order = constraint_order
548 self.constraint_ranges = constraint_ranges
551 def process_subgraph(self, subgraph):
553 (task_node,) = subgraph.tasks.values()
557 (input_dataset_type_node,) = subgraph.inputs_of(task_node.label).values()
558 assert input_dataset_type_node
is not None,
"PipelineGraph should be resolved by base class."
559 (output_edge,) = task_node.outputs.values()
560 output_dataset_type_node = subgraph.dataset_types[output_edge.parent_dataset_type_name]
561 (hpx_output_dimension,) = (
562 self.butler.dimensions.skypix_dimensions[d]
563 for d
in output_dataset_type_node.dimensions.skypix.names
565 constraint_hpx_pixelization = (
566 self.butler.dimensions.skypix_dimensions[f
"healpix{self.constraint_order}"].pixelization
568 common_skypix_name = self.butler.dimensions.commonSkyPix.name
569 common_skypix_pixelization = self.butler.dimensions.commonSkyPix.pixelization
575 self.butler.dimensions.skypix_dimensions[d]
for d
in task_node.dimensions.names
if d !=
"band"
577 hpx_pixelization = hpx_dimension.pixelization
578 if hpx_pixelization.level < self.constraint_order:
579 raise ValueError(f
"Quantum order {hpx_pixelization.level} must be < {self.constraint_order}")
580 hpx_ranges = self.constraint_ranges.scaled(4**(hpx_pixelization.level - self.constraint_order))
585 for begin, end
in self.constraint_ranges:
586 for hpx_index
in range(begin, end):
587 constraint_hpx_region = constraint_hpx_pixelization.pixel(hpx_index)
588 common_skypix_ranges |= common_skypix_pixelization.envelope(constraint_hpx_region)
592 for simp
in range(1, 10):
593 if len(common_skypix_ranges) < 100:
595 common_skypix_ranges.simplify(simp)
602 for n, (begin, end)
in enumerate(common_skypix_ranges):
605 where_terms.append(f
"{common_skypix_name} = cpx{n}")
606 bind[f
"cpx{n}"] = begin
608 where_terms.append(f
"({common_skypix_name} >= cpx{n}a AND {common_skypix_name} <= cpx{n}b)")
609 bind[f
"cpx{n}a"] = begin
610 bind[f
"cpx{n}b"] = stop
611 where =
" OR ".join(where_terms)
613 where = f
"({self.where}) AND ({where})"
617 input_refs = self.butler.registry.queryDatasets(
618 input_dataset_type_node.dataset_type,
621 collections=self.input_collections,
624 inputs_by_patch = defaultdict(set)
625 patch_dimensions = self.butler.dimensions.conform([
"patch"])
626 for input_ref
in input_refs:
627 dataset_key = DatasetKey(input_ref.datasetType.name, input_ref.dataId.required_values)
628 self.existing_datasets.inputs[dataset_key] = input_ref
629 inputs_by_patch[input_ref.dataId.subset(patch_dimensions)].add(dataset_key)
630 if not inputs_by_patch:
631 message_body =
"\n".join(input_refs.explain_no_results())
632 raise RuntimeError(f
"No inputs found:\n{message_body}")
637 inputs_by_hpx = defaultdict(set)
638 for patch_data_id, input_keys_for_patch
in inputs_by_patch.items():
639 patch_hpx_ranges = hpx_pixelization.envelope(patch_data_id.region)
640 for begin, end
in patch_hpx_ranges & hpx_ranges:
641 for hpx_index
in range(begin, end):
642 inputs_by_hpx[hpx_index].update(input_keys_for_patch)
645 skeleton = QuantumGraphSkeleton([task_node.label])
646 for hpx_index, input_keys_for_hpx_index
in inputs_by_hpx.items():
648 input_keys_by_band = defaultdict(list)
649 for input_key
in input_keys_for_hpx_index:
650 input_ref = self.existing_datasets.inputs[input_key]
651 input_keys_by_band[input_ref.dataId[
"band"]].append(input_key)
653 for band, input_keys_for_band
in input_keys_by_band.items():
654 data_id = self.butler.registry.expandDataId({hpx_dimension.name: hpx_index,
"band": band})
655 quantum_key = skeleton.add_quantum_node(task_node.label, data_id)
657 skeleton.add_input_edges(quantum_key, input_keys_for_band)
659 hpx_pixel_ranges =
RangeSet(hpx_index)
660 hpx_output_ranges = hpx_pixel_ranges.scaled(
661 4**(task_node.config.hips_order - hpx_pixelization.level)
663 for begin, end
in hpx_output_ranges:
664 for hpx_output_index
in range(begin, end):
665 dataset_key = skeleton.add_dataset_node(
666 output_dataset_type_node.name,
667 self.butler.registry.expandDataId(
668 {hpx_output_dimension: hpx_output_index,
"band": band}
671 skeleton.add_output_edge(quantum_key, dataset_key)
673 for write_edge
in task_node.iter_all_outputs():
674 if write_edge.connection_name == output_edge.connection_name:
676 dataset_key = skeleton.add_dataset_node(write_edge.parent_dataset_type_name, data_id)
677 skeleton.add_output_edge(quantum_key, dataset_key)
681class HipsPropertiesSpectralTerm(pexConfig.Config):
682 lambda_min = pexConfig.Field(
683 doc=
"Minimum wavelength (nm)",
686 lambda_max = pexConfig.Field(
687 doc=
"Maximum wavelength (nm)",
692class HipsPropertiesConfig(pexConfig.Config):
693 """Configuration parameters for writing a HiPS properties file."""
694 creator_did_template = pexConfig.Field(
695 doc=(
"Unique identifier of the HiPS - Format: IVOID. "
696 "Use ``{band}`` to substitute the band name."),
700 obs_collection = pexConfig.Field(
701 doc=
"Short name of original data set - Format: one word",
705 obs_description_template = pexConfig.Field(
706 doc=(
"Data set description - Format: free text, longer free text "
707 "description of the dataset. Use ``{band}`` to substitute "
711 prov_progenitor = pexConfig.ListField(
712 doc=
"Provenance of the original data - Format: free text",
716 obs_title_template = pexConfig.Field(
717 doc=(
"Data set title format: free text, but should be short. "
718 "Use ``{band}`` to substitute the band name."),
722 spectral_ranges = pexConfig.ConfigDictField(
723 doc=(
"Mapping from band to lambda_min, lamba_max (nm). May be approximate."),
725 itemtype=HipsPropertiesSpectralTerm,
728 initial_ra = pexConfig.Field(
729 doc=
"Initial RA (deg) (default for HiPS viewer). If not set will use a point in MOC.",
733 initial_dec = pexConfig.Field(
734 doc=
"Initial Declination (deg) (default for HiPS viewer). If not set will use a point in MOC.",
738 initial_fov = pexConfig.Field(
739 doc=
"Initial field-of-view (deg). If not set will use ~1 healpix tile.",
743 obs_ack = pexConfig.Field(
744 doc=
"Observation acknowledgements (free text).",
748 t_min = pexConfig.Field(
749 doc=
"Time (MJD) of earliest observation included in HiPS",
753 t_max = pexConfig.Field(
754 doc=
"Time (MJD) of latest observation included in HiPS",
762 if self.obs_collection
is not None:
763 if re.search(
r"\s", self.obs_collection):
764 raise ValueError(
"obs_collection cannot contain any space characters.")
766 def setDefaults(self):
769 u_term = HipsPropertiesSpectralTerm()
770 u_term.lambda_min = 330.
771 u_term.lambda_max = 400.
772 self.spectral_ranges[
"u"] = u_term
773 g_term = HipsPropertiesSpectralTerm()
774 g_term.lambda_min = 402.
775 g_term.lambda_max = 552.
776 self.spectral_ranges[
"g"] = g_term
777 r_term = HipsPropertiesSpectralTerm()
778 r_term.lambda_min = 552.
779 r_term.lambda_max = 691.
780 self.spectral_ranges[
"r"] = r_term
781 i_term = HipsPropertiesSpectralTerm()
782 i_term.lambda_min = 691.
783 i_term.lambda_max = 818.
784 self.spectral_ranges[
"i"] = i_term
785 z_term = HipsPropertiesSpectralTerm()
786 z_term.lambda_min = 818.
787 z_term.lambda_max = 922.
788 self.spectral_ranges[
"z"] = z_term
789 y_term = HipsPropertiesSpectralTerm()
790 y_term.lambda_min = 970.
791 y_term.lambda_max = 1060.
792 self.spectral_ranges[
"y"] = y_term
795class GenerateHipsConnections(pipeBase.PipelineTaskConnections,
796 dimensions=(
"instrument",
"band"),
797 defaultTemplates={
"coaddName":
"deep"}):
798 hips_exposure_handles = pipeBase.connectionTypes.Input(
799 doc=
"HiPS-compatible HPX images.",
800 name=
"{coaddName}Coadd_hpx",
801 storageClass=
"ExposureF",
802 dimensions=(
"healpix11",
"band"),
808class GenerateHipsConfig(pipeBase.PipelineTaskConfig,
809 pipelineConnections=GenerateHipsConnections):
810 """Configuration parameters for GenerateHipsTask."""
815 hips_base_uri = pexConfig.Field(
816 doc=
"URI to HiPS base for output.",
820 min_order = pexConfig.Field(
821 doc=
"Minimum healpix order for HiPS tree.",
825 properties = pexConfig.ConfigField(
826 dtype=HipsPropertiesConfig,
827 doc=
"Configuration for properties file.",
829 allsky_tilesize = pexConfig.Field(
831 doc=
"Allsky tile size; must be power of 2. HiPS standard recommends 64x64 tiles.",
833 check=_is_power_of_two,
835 png_gray_asinh_minimum = pexConfig.Field(
836 doc=
"AsinhMapping intensity to be mapped to black for grayscale png scaling.",
840 png_gray_asinh_stretch = pexConfig.Field(
841 doc=
"AsinhMapping linear stretch for grayscale png scaling.",
845 png_gray_asinh_softening = pexConfig.Field(
846 doc=
"AsinhMapping softening parameter (Q) for grayscale png scaling.",
852class GenerateHipsTask(pipeBase.PipelineTask):
853 """Task for making a HiPS tree with FITS and grayscale PNGs."""
854 ConfigClass = GenerateHipsConfig
855 _DefaultName =
"generateHips"
859 def runQuantum(self, butlerQC, inputRefs, outputRefs):
860 inputs = butlerQC.get(inputRefs)
862 dims = inputRefs.hips_exposure_handles[0].dataId.dimensions.names
866 order = int(dim.split(
"healpix")[1])
870 raise RuntimeError(
"Could not determine healpix order for input exposures.")
872 hips_exposure_handle_dict = {
873 (hips_exposure_handle.dataId[healpix_dim],
874 hips_exposure_handle.dataId[
"band"]): hips_exposure_handle
875 for hips_exposure_handle
in inputs[
"hips_exposure_handles"]
878 data_bands = {hips_exposure_handle.dataId[
"band"]
879 for hips_exposure_handle
in inputs[
"hips_exposure_handles"]}
880 bands = self._check_data_bands(data_bands)
885 hips_exposure_handle_dict=hips_exposure_handle_dict,
886 do_color=self.color_task,
889 def _check_data_bands(self, data_bands):
890 """Check that the data has only a single band.
894 data_bands : `set` [`str`]
895 Bands from the input data.
899 bands : `list` [`str`]
900 List of single band to process.
904 RuntimeError if there is not exactly one band.
906 if len(data_bands) != 1:
907 raise RuntimeError(
"GenerateHipsTask can only use data from a single band.")
909 return list(data_bands)
912 def run(self, bands, max_order, hips_exposure_handle_dict, do_color=False):
913 """Run the GenerateHipsTask.
917 bands : `list [ `str` ]
918 List of bands to be processed (or single band).
920 HEALPix order of the maximum (native) HPX exposures.
921 hips_exposure_handle_dict : `dict` {`int`: `lsst.daf.butler.DeferredDatasetHandle`}
922 Dict of handles for the HiPS high-resolution exposures.
923 Key is (pixel number, ``band``).
924 do_color : `bool`, optional
925 Do color pngs instead of per-band grayscale.
927 min_order = self.config.min_order
930 png_grayscale_mapping = AsinhMapping(
931 self.config.png_gray_asinh_minimum,
932 self.config.png_gray_asinh_stretch,
933 Q=self.config.png_gray_asinh_softening,
936 png_color_mapping = AsinhMapping(
937 self.config.png_color_asinh_minimum,
938 self.config.png_color_asinh_stretch,
939 Q=self.config.png_color_asinh_softening,
942 bcb = self.config.blue_channel_band
943 gcb = self.config.green_channel_band
944 rcb = self.config.red_channel_band
945 colorstr = f
"{bcb}{gcb}{rcb}"
948 hips_base_path = ResourcePath(self.config.hips_base_uri, forceDirectory=
True)
952 pixels = np.unique(np.array([pixel
for pixel, _
in hips_exposure_handle_dict.keys()]))
955 pixels = np.append(pixels, [0])
959 pixels_shifted[max_order] = pixels
960 for order
in range(max_order - 1, min_order - 1, -1):
961 pixels_shifted[order] = np.right_shift(pixels_shifted[order + 1], 2)
964 for order
in range(min_order, max_order + 1):
965 pixels_shifted[order][-1] = -1
968 exp0 = list(hips_exposure_handle_dict.values())[0].get()
969 bbox = exp0.getBBox()
970 npix = bbox.getWidth()
971 shift_order = int(np.round(np.log2(npix)))
978 for order
in range(min_order, max_order + 1):
979 exp = exp0.Factory(bbox=bbox)
980 exp.image.array[:, :] = np.nan
981 exposures[(band, order)] = exp
984 for pixel_counter, pixel
in enumerate(pixels[:-1]):
985 self.log.debug(
"Working on high resolution pixel %d", pixel)
991 if (pixel, band)
in hips_exposure_handle_dict:
992 exposures[(band, max_order)] = hips_exposure_handle_dict[(pixel, band)].get()
998 for order
in range(max_order, min_order - 1, -1):
999 if pixels_shifted[order][pixel_counter + 1] == pixels_shifted[order][pixel_counter]:
1007 self._write_hips_image(
1008 hips_base_path.join(f
"band_{band}", forceDirectory=
True),
1010 pixels_shifted[order][pixel_counter],
1011 exposures[(band, order)].image,
1012 png_grayscale_mapping,
1013 shift_order=shift_order,
1017 self._write_hips_color_png(
1018 hips_base_path.join(f
"color_{colorstr}", forceDirectory=
True),
1020 pixels_shifted[order][pixel_counter],
1021 exposures[(self.config.red_channel_band, order)].image,
1022 exposures[(self.config.green_channel_band, order)].image,
1023 exposures[(self.config.blue_channel_band, order)].image,
1027 log_level = self.log.INFO
if order == (max_order - 3)
else self.log.DEBUG
1030 "Completed HiPS generation for %s, order %d, pixel %d (%d/%d)",
1033 pixels_shifted[order][pixel_counter],
1039 if order == min_order:
1041 exposures[(band, order)].image.array[:, :] = np.nan
1046 arr = exposures[(band, order)].image.array.reshape(npix//2, 2, npix//2, 2)
1047 with warnings.catch_warnings():
1048 warnings.simplefilter(
"ignore")
1049 binned_image_arr = np.nanmean(arr, axis=(1, 3))
1053 sub_index = (pixels_shifted[order][pixel_counter]
1054 - np.left_shift(pixels_shifted[order - 1][pixel_counter], 2))
1057 exp = exposures[(band, order - 1)]
1061 exp.image.array[npix//2:, 0: npix//2] = binned_image_arr
1062 elif sub_index == 1:
1063 exp.image.array[0: npix//2, 0: npix//2] = binned_image_arr
1064 elif sub_index == 2:
1065 exp.image.array[npix//2:, npix//2:] = binned_image_arr
1066 elif sub_index == 3:
1067 exp.image.array[0: npix//2, npix//2:] = binned_image_arr
1070 raise ValueError(
"Illegal pixel sub index")
1073 if order < max_order:
1074 exposures[(band, order)].image.array[:, :] = np.nan
1079 band_pixels = np.array([pixel
1080 for pixel, band_
in hips_exposure_handle_dict.keys()
1082 band_pixels = np.sort(band_pixels)
1084 self._write_properties_and_moc(
1085 hips_base_path.join(f
"band_{band}", forceDirectory=
True),
1093 self._write_allsky_file(
1094 hips_base_path.join(f
"band_{band}", forceDirectory=
True),
1098 self._write_properties_and_moc(
1099 hips_base_path.join(f
"color_{colorstr}", forceDirectory=
True),
1107 self._write_allsky_file(
1108 hips_base_path.join(f
"color_{colorstr}", forceDirectory=
True),
1112 def _write_hips_image(self, hips_base_path, order, pixel, image, png_mapping, shift_order=9):
1113 """Write a HiPS image.
1117 hips_base_path : `lsst.resources.ResourcePath`
1118 Resource path to the base of the HiPS directory tree.
1120 HEALPix order of the HiPS image to write.
1122 HEALPix pixel of the HiPS image.
1123 image : `lsst.afw.image.Image`
1125 png_mapping : `astropy.visualization.lupton_rgb.AsinhMapping`
1126 Mapping to convert image to scaled png.
1127 shift_order : `int`, optional
1135 dir_number = self._get_dir_number(pixel)
1136 hips_dir = hips_base_path.join(
1144 wcs = makeHpxWcs(order, pixel, shift_order=shift_order)
1146 uri = hips_dir.join(f
"Npix{pixel}.fits")
1148 with ResourcePath.temporary_uri(suffix=uri.getExtension())
as temporary_uri:
1149 image.writeFits(temporary_uri.ospath, metadata=wcs.getFitsMetadata())
1151 uri.transfer_from(temporary_uri, transfer=
"copy", overwrite=
True)
1155 with np.errstate(invalid=
"ignore"):
1156 vals = 255 - png_mapping.map_intensity_to_uint8(image.array).astype(np.uint8)
1158 vals[~np.isfinite(image.array) | (image.array < 0)] = 0
1159 im = Image.fromarray(vals[::-1, :],
"L")
1161 uri = hips_dir.join(f
"Npix{pixel}.png")
1163 with ResourcePath.temporary_uri(suffix=uri.getExtension())
as temporary_uri:
1164 im.save(temporary_uri.ospath)
1166 uri.transfer_from(temporary_uri, transfer=
"copy", overwrite=
True)
1168 def _write_hips_color_png(
1178 """Write a color png HiPS image.
1182 hips_base_path : `lsst.resources.ResourcePath`
1183 Resource path to the base of the HiPS directory tree.
1185 HEALPix order of the HiPS image to write.
1187 HEALPix pixel of the HiPS image.
1188 image_red : `lsst.afw.image.Image`
1189 Input for red channel of output png.
1190 image_green : `lsst.afw.image.Image`
1191 Input for green channel of output png.
1192 image_blue : `lsst.afw.image.Image`
1193 Input for blue channel of output png.
1194 png_mapping : `astropy.visualization.lupton_rgb.AsinhMapping`
1195 Mapping to convert image to scaled png.
1202 dir_number = self._get_dir_number(pixel)
1203 hips_dir = hips_base_path.join(
1212 arr_red = image_red.array.copy()
1213 arr_red[np.isnan(arr_red)] = png_mapping.minimum[0]
1214 arr_green = image_green.array.copy()
1215 arr_green[np.isnan(arr_green)] = png_mapping.minimum[1]
1216 arr_blue = image_blue.array.copy()
1217 arr_blue[np.isnan(arr_blue)] = png_mapping.minimum[2]
1219 image_array = png_mapping.make_rgb_image(arr_red, arr_green, arr_blue)
1221 im = Image.fromarray(image_array[::-1, :, :], mode=
"RGB")
1223 uri = hips_dir.join(f
"Npix{pixel}.png")
1225 with ResourcePath.temporary_uri(suffix=uri.getExtension())
as temporary_uri:
1226 im.save(temporary_uri.ospath)
1228 uri.transfer_from(temporary_uri, transfer=
"copy", overwrite=
True)
1230 def _write_properties_and_moc(
1240 """Write HiPS properties file and MOC.
1244 hips_base_path : : `lsst.resources.ResourcePath`
1245 Resource path to the base of the HiPS directory tree.
1247 Maximum HEALPix order.
1248 pixels : `np.ndarray` (N,)
1249 Array of pixels used.
1250 exposure : `lsst.afw.image.Exposure`
1251 Sample HPX exposure used for generating HiPS tiles.
1257 Is band multiband / color?
1259 area = hpg.nside_to_pixel_area(2**max_order, degrees=
True)*len(pixels)
1261 initial_ra = self.config.properties.initial_ra
1262 initial_dec = self.config.properties.initial_dec
1263 initial_fov = self.config.properties.initial_fov
1265 if initial_ra
is None or initial_dec
is None or initial_fov
is None:
1268 temp_pixels = pixels.copy()
1269 if temp_pixels.size % 2 == 0:
1270 temp_pixels = np.append(temp_pixels, [temp_pixels[0]])
1271 medpix = int(np.median(temp_pixels))
1272 _initial_ra, _initial_dec = hpg.pixel_to_angle(2**max_order, medpix)
1273 _initial_fov = hpg.nside_to_resolution(2**max_order, units=
'arcminutes')/60.
1275 if initial_ra
is None or initial_dec
is None:
1276 initial_ra = _initial_ra
1277 initial_dec = _initial_dec
1278 if initial_fov
is None:
1279 initial_fov = _initial_fov
1281 self._write_hips_properties_file(
1283 self.config.properties,
1296 self._write_hips_moc_file(
1302 def _write_hips_properties_file(
1316 """Write HiPS properties file.
1320 hips_base_path : `lsst.resources.ResourcePath`
1321 ResourcePath at top of HiPS tree. File will be written
1322 to this path as ``properties``.
1323 properties_config : `lsst.pipe.tasks.hips.HipsPropertiesConfig`
1324 Configuration for properties values.
1326 Name of band(s) for HiPS tree.
1328 Is multiband / color?
1329 exposure : `lsst.afw.image.Exposure`
1330 Sample HPX exposure used for generating HiPS tiles.
1332 Maximum HEALPix order.
1336 Coverage area in square degrees.
1337 initial_ra : `float`
1338 Initial HiPS RA position (degrees).
1339 initial_dec : `float`
1340 Initial HiPS Dec position (degrees).
1341 initial_fov : `float`
1342 Initial HiPS display size (degrees).
1348 def _write_property(fh, name, value):
1349 """Write a property name/value to a file handle.
1353 fh : file handle (blah)
1362 if re.search(
r"\s", name):
1363 raise ValueError(f
"``{name}`` cannot contain any space characters.")
1365 raise ValueError(f
"``{name}`` cannot contain an ``=``")
1367 fh.write(f
"{name:25}= {value}\n")
1369 if exposure.image.array.dtype == np.dtype(
"float32"):
1371 elif exposure.image.array.dtype == np.dtype(
"float64"):
1373 elif exposure.image.array.dtype == np.dtype(
"int32"):
1376 date_iso8601 = datetime.utcnow().isoformat(timespec=
"seconds") +
"Z"
1377 pixel_scale = hpg.nside_to_resolution(2**(max_order + shift_order), units=
'degrees')
1379 uri = hips_base_path.join(
"properties")
1380 with ResourcePath.temporary_uri(suffix=uri.getExtension())
as temporary_uri:
1381 with open(temporary_uri.ospath,
"w")
as fh:
1385 properties_config.creator_did_template.format(band=band),
1387 if properties_config.obs_collection
is not None:
1388 _write_property(fh,
"obs_collection", properties_config.obs_collection)
1392 properties_config.obs_title_template.format(band=band),
1394 if properties_config.obs_description_template
is not None:
1398 properties_config.obs_description_template.format(band=band),
1400 if len(properties_config.prov_progenitor) > 0:
1401 for prov_progenitor
in properties_config.prov_progenitor:
1402 _write_property(fh,
"prov_progenitor", prov_progenitor)
1403 if properties_config.obs_ack
is not None:
1404 _write_property(fh,
"obs_ack", properties_config.obs_ack)
1405 _write_property(fh,
"obs_regime",
"Optical")
1406 _write_property(fh,
"data_pixel_bitpix", str(bitpix))
1407 _write_property(fh,
"dataproduct_type",
"image")
1408 _write_property(fh,
"moc_sky_fraction", str(area/41253.))
1409 _write_property(fh,
"data_ucd",
"phot.flux")
1410 _write_property(fh,
"hips_creation_date", date_iso8601)
1411 _write_property(fh,
"hips_builder",
"lsst.pipe.tasks.hips.GenerateHipsTask")
1412 _write_property(fh,
"hips_creator",
"Vera C. Rubin Observatory")
1413 _write_property(fh,
"hips_version",
"1.4")
1414 _write_property(fh,
"hips_release_date", date_iso8601)
1415 _write_property(fh,
"hips_frame",
"equatorial")
1416 _write_property(fh,
"hips_order", str(max_order))
1417 _write_property(fh,
"hips_tile_width", str(exposure.getBBox().getWidth()))
1418 _write_property(fh,
"hips_status",
"private master clonableOnce")
1420 _write_property(fh,
"hips_tile_format",
"png")
1421 _write_property(fh,
"dataproduct_subtype",
"color")
1423 _write_property(fh,
"hips_tile_format",
"png fits")
1424 _write_property(fh,
"hips_pixel_bitpix", str(bitpix))
1425 _write_property(fh,
"hips_pixel_scale", str(pixel_scale))
1426 _write_property(fh,
"hips_initial_ra", str(initial_ra))
1427 _write_property(fh,
"hips_initial_dec", str(initial_dec))
1428 _write_property(fh,
"hips_initial_fov", str(initial_fov))
1430 if self.config.blue_channel_band
in properties_config.spectral_ranges:
1431 em_min = properties_config.spectral_ranges[
1432 self.config.blue_channel_band
1435 self.log.warning(
"blue band %s not in self.config.spectral_ranges.", band)
1437 if self.config.red_channel_band
in properties_config.spectral_ranges:
1438 em_max = properties_config.spectral_ranges[
1439 self.config.red_channel_band
1442 self.log.warning(
"red band %s not in self.config.spectral_ranges.", band)
1445 if band
in properties_config.spectral_ranges:
1446 em_min = properties_config.spectral_ranges[band].lambda_min/1e9
1447 em_max = properties_config.spectral_ranges[band].lambda_max/1e9
1449 self.log.warning(
"band %s not in self.config.spectral_ranges.", band)
1452 _write_property(fh,
"em_min", str(em_min))
1453 _write_property(fh,
"em_max", str(em_max))
1454 if properties_config.t_min
is not None:
1455 _write_property(fh,
"t_min", properties_config.t_min)
1456 if properties_config.t_max
is not None:
1457 _write_property(fh,
"t_max", properties_config.t_max)
1459 uri.transfer_from(temporary_uri, transfer=
"copy", overwrite=
True)
1461 def _write_hips_moc_file(self, hips_base_path, max_order, pixels, min_uniq_order=1):
1462 """Write HiPS MOC file.
1466 hips_base_path : `lsst.resources.ResourcePath`
1467 ResourcePath to top of HiPS tree. File will be written as
1468 to this path as ``Moc.fits``.
1470 Maximum HEALPix order.
1471 pixels : `np.ndarray`
1472 Array of pixels covered.
1473 min_uniq_order : `int`, optional
1474 Minimum HEALPix order for looking for fully covered pixels.
1482 uniq = 4*(4**max_order) + pixels
1485 hspmap = hsp.HealSparseMap.make_empty(2**min_uniq_order, 2**max_order, dtype=np.float32)
1486 hspmap[pixels] = 1.0
1489 for uniq_order
in range(max_order - 1, min_uniq_order - 1, -1):
1490 hspmap = hspmap.degrade(2**uniq_order, reduction=
"sum")
1491 pix_shift = np.right_shift(pixels, 2*(max_order - uniq_order))
1493 covered, = np.isclose(hspmap[pix_shift], 4**(max_order - uniq_order)).nonzero()
1494 if covered.size == 0:
1498 uniq[covered] = 4*(4**uniq_order) + pix_shift[covered]
1501 uniq = np.unique(uniq)
1504 tbl = np.zeros(uniq.size, dtype=[(
"UNIQ",
"i8")])
1507 order = np.log2(tbl[
"UNIQ"]//4).astype(np.int32)//2
1508 moc_order = np.max(order)
1510 hdu = fits.BinTableHDU(tbl)
1511 hdu.header[
"PIXTYPE"] =
"HEALPIX"
1512 hdu.header[
"ORDERING"] =
"NUNIQ"
1513 hdu.header[
"COORDSYS"] =
"C"
1514 hdu.header[
"MOCORDER"] = moc_order
1515 hdu.header[
"MOCTOOL"] =
"lsst.pipe.tasks.hips.GenerateHipsTask"
1517 uri = hips_base_path.join(
"Moc.fits")
1519 with ResourcePath.temporary_uri(suffix=uri.getExtension())
as temporary_uri:
1520 hdu.writeto(temporary_uri.ospath)
1522 uri.transfer_from(temporary_uri, transfer=
"copy", overwrite=
True)
1524 def _write_allsky_file(self, hips_base_path, allsky_order):
1525 """Write an Allsky.png file.
1529 hips_base_path : `lsst.resources.ResourcePath`
1530 Resource path to the base of the HiPS directory tree.
1531 allsky_order : `int`
1532 HEALPix order of the minimum order to make allsky file.
1534 tile_size = self.config.allsky_tilesize
1546 n_tiles = hpg.nside_to_npixel(hpg.order_to_nside(allsky_order))
1547 n_tiles_wide = int(np.floor(np.sqrt(n_tiles)))
1548 n_tiles_high = int(np.ceil(n_tiles / n_tiles_wide))
1552 allsky_order_uri = hips_base_path.join(f
"Norder{allsky_order}", forceDirectory=
True)
1553 pixel_regex = re.compile(
r"Npix([0-9]+)\.png$")
1555 ResourcePath.findFileResources(
1556 candidates=[allsky_order_uri],
1557 file_filter=pixel_regex,
1561 for png_uri
in png_uris:
1562 matches = re.match(pixel_regex, png_uri.basename())
1563 pix_num = int(matches.group(1))
1564 tile_image = Image.open(io.BytesIO(png_uri.read()))
1565 row = math.floor(pix_num//n_tiles_wide)
1566 column = pix_num % n_tiles_wide
1567 box = (column*tile_size, row*tile_size, (column + 1)*tile_size, (row + 1)*tile_size)
1568 tile_image_shrunk = tile_image.resize((tile_size, tile_size))
1570 if allsky_image
is None:
1571 allsky_image = Image.new(
1573 (n_tiles_wide*tile_size, n_tiles_high*tile_size),
1575 allsky_image.paste(tile_image_shrunk, box)
1577 uri = allsky_order_uri.join(
"Allsky.png")
1579 with ResourcePath.temporary_uri(suffix=uri.getExtension())
as temporary_uri:
1580 allsky_image.save(temporary_uri.ospath)
1582 uri.transfer_from(temporary_uri, transfer=
"copy", overwrite=
True)
1584 def _get_dir_number(self, pixel):
1585 """Compute the directory number from a pixel.
1590 HEALPix pixel number.
1595 HiPS directory number.
1597 return (pixel//10000)*10000
1600class GenerateColorHipsConnections(pipeBase.PipelineTaskConnections,
1601 dimensions=(
"instrument", ),
1602 defaultTemplates={
"coaddName":
"deep"}):
1603 hips_exposure_handles = pipeBase.connectionTypes.Input(
1604 doc=
"HiPS-compatible HPX images.",
1605 name=
"{coaddName}Coadd_hpx",
1606 storageClass=
"ExposureF",
1607 dimensions=(
"healpix11",
"band"),
1613class GenerateColorHipsConfig(GenerateHipsConfig,
1614 pipelineConnections=GenerateColorHipsConnections):
1615 """Configuration parameters for GenerateColorHipsTask."""
1616 blue_channel_band = pexConfig.Field(
1617 doc=
"Band to use for blue channel of color pngs.",
1621 green_channel_band = pexConfig.Field(
1622 doc=
"Band to use for green channel of color pngs.",
1626 red_channel_band = pexConfig.Field(
1627 doc=
"Band to use for red channel of color pngs.",
1631 png_color_asinh_minimum = pexConfig.Field(
1632 doc=
"AsinhMapping intensity to be mapped to black for color png scaling.",
1636 png_color_asinh_stretch = pexConfig.Field(
1637 doc=
"AsinhMapping linear stretch for color png scaling.",
1641 png_color_asinh_softening = pexConfig.Field(
1642 doc=
"AsinhMapping softening parameter (Q) for color png scaling.",
1648class GenerateColorHipsTask(GenerateHipsTask):
1649 """Task for making a HiPS tree with color pngs."""
1650 ConfigClass = GenerateColorHipsConfig
1651 _DefaultName =
"generateColorHips"
1654 def _check_data_bands(self, data_bands):
1655 """Check the data for configured bands.
1657 Warn if any color bands are missing data.
1661 data_bands : `set` [`str`]
1662 Bands from the input data.
1666 bands : `list` [`str`]
1667 List of bands in bgr color order.
1669 if len(data_bands) == 0:
1670 raise RuntimeError(
"GenerateColorHipsTask must have data from at least one band.")
1672 if self.config.blue_channel_band
not in data_bands:
1674 "Color png blue_channel_band %s not in dataset.",
1675 self.config.blue_channel_band
1677 if self.config.green_channel_band
not in data_bands:
1679 "Color png green_channel_band %s not in dataset.",
1680 self.config.green_channel_band
1682 if self.config.red_channel_band
not in data_bands:
1684 "Color png red_channel_band %s not in dataset.",
1685 self.config.red_channel_band
1689 self.config.blue_channel_band,
1690 self.config.green_channel_band,
1691 self.config.red_channel_band,
Pass parameters to a Statistics object.
An integer coordinate rectangle.
A RangeSet is a set of unsigned 64 bit integers.
daf::base::PropertySet * set
std::shared_ptr< lsst::afw::image::Image< PixelT > > statisticsStack(std::vector< std::shared_ptr< lsst::afw::image::Image< PixelT > > > &images, Property flags, StatisticsControl const &sctrl=StatisticsControl(), std::vector< lsst::afw::image::VariancePixel > const &wvector=std::vector< lsst::afw::image::VariancePixel >(0))
A function to compute some statistics of a stack of Images.
Property stringToStatisticsProperty(std::string const property)
Conversion function to switch a string to a Property (see Statistics.h)