22from collections
import defaultdict
27import healsparse
as hsp
33from lsst.daf.butler
import Formatter
35from lsst.utils.timer
import timeMethod
36from .healSparseMappingProperties
import (BasePropertyMap, BasePropertyMapConfig,
37 PropertyMapMap, compute_approx_psf_size_and_shape)
40__all__ = [
"HealSparseInputMapTask",
"HealSparseInputMapConfig",
41 "HealSparseMapFormatter",
"HealSparsePropertyMapConnections",
42 "HealSparsePropertyMapConfig",
"HealSparsePropertyMapTask",
43 "ConsolidateHealSparsePropertyMapConnections",
44 "ConsolidateHealSparsePropertyMapConfig",
45 "ConsolidateHealSparsePropertyMapTask"]
49 """Interface for reading and writing healsparse.HealSparseMap files."""
50 unsupportedParameters = frozenset()
51 supportedExtensions = frozenset({
".hsp",
".fit",
".fits"})
54 def read(self, component=None):
56 path = self.fileDescriptor.location.path
58 if component ==
'coverage':
60 data = hsp.HealSparseCoverage.read(path)
61 except (OSError, RuntimeError):
62 raise ValueError(f
"Unable to read healsparse map with URI {self.fileDescriptor.location.uri}")
66 if self.fileDescriptor.parameters
is None:
70 pixels = self.fileDescriptor.parameters.get(
'pixels',
None)
71 degrade_nside = self.fileDescriptor.parameters.get(
'degrade_nside',
None)
73 data = hsp.HealSparseMap.read(path, pixels=pixels, degrade_nside=degrade_nside)
74 except (OSError, RuntimeError):
75 raise ValueError(f
"Unable to read healsparse map with URI {self.fileDescriptor.location.uri}")
79 def write(self, inMemoryDataset):
82 self.fileDescriptor.location.updateExtension(self.
extension)
83 inMemoryDataset.write(self.fileDescriptor.location.path, clobber=
True)
86def _is_power_of_two(value):
87 """Check that value is a power of two.
96 is_power_of_two : `bool`
97 True if value
is a power of two;
False otherwise,
or
98 if value
is not an integer.
100 if not isinstance(value, numbers.Integral):
107 return (value & (value - 1) == 0)
and value != 0
111 """Configuration parameters for HealSparseInputMapTask"""
112 nside = pexConfig.Field(
113 doc=
"Mapping healpix nside. Must be power of 2.",
116 check=_is_power_of_two,
118 nside_coverage = pexConfig.Field(
119 doc=
"HealSparse coverage map nside. Must be power of 2.",
122 check=_is_power_of_two,
124 bad_mask_min_coverage = pexConfig.Field(
125 doc=(
"Minimum area fraction of a map healpixel pixel that must be "
126 "covered by bad pixels to be removed from the input map. "
127 "This is approximate."),
134 """Task for making a HealSparse input map."""
136 ConfigClass = HealSparseInputMapConfig
137 _DefaultName =
"healSparseInputMap"
140 pipeBase.Task.__init__(self, **kwargs)
145 """Build a map from ccd valid polygons or bounding boxes.
150 Bounding box for region to build input map.
152 WCS object
for region to build input map.
154 Exposure catalog
with ccd data
from coadd inputs.
156 with warnings.catch_warnings():
161 warnings.simplefilter(
"ignore")
162 self.
ccd_input_map = hsp.HealSparseMap.make_empty(nside_coverage=self.config.nside_coverage,
163 nside_sparse=self.config.nside,
165 wide_mask_maxbits=len(ccds))
170 pixel_scale = wcs.getPixelScale().asArcseconds()
171 hpix_area_arcsec2 = hpg.nside_to_pixel_area(self.config.nside, degrees=
True)*(3600.**2.)
172 self.
_min_bad = self.config.bad_mask_min_coverage*hpix_area_arcsec2/(pixel_scale**2.)
177 for bit, ccd_row
in enumerate(ccds):
178 metadata[f
"B{bit:04d}CCD"] = ccd_row[
"ccd"]
179 metadata[f
"B{bit:04d}VIS"] = ccd_row[
"visit"]
180 metadata[f
"B{bit:04d}WT"] = ccd_row[
"weight"]
185 ccd_poly = ccd_row.getValidPolygon()
189 ccd_poly_radec = self.
_pixels_to_radec(ccd_row.getWcs(), ccd_poly.convexHull().getVertices())
192 poly = hsp.Polygon(ra=ccd_poly_radec[: -1, 0],
193 dec=ccd_poly_radec[: -1, 1],
201 bbox_afw_poly.convexHull().getVertices())
202 bbox_poly = hsp.Polygon(ra=bbox_poly_radec[: -1, 0], dec=bbox_poly_radec[: -1, 1],
211 dtype = [(f
"v{visit}", np.int64)
for visit
in self.
_bits_per_visit.keys()]
213 with warnings.catch_warnings():
218 warnings.simplefilter(
"ignore")
220 nside_coverage=self.config.nside_coverage,
221 nside_sparse=self.config.nside,
230 """Mask a subregion from a visit.
231 This must be run after build_ccd_input_map initializes
237 Bounding box from region to mask.
239 Visit number corresponding to warp
with mask.
240 mask : `lsst.afw.image.MaskX`
241 Mask plane
from warp exposure.
242 bit_mask_value : `int`
243 Bit mask to check
for bad pixels.
247 RuntimeError : Raised
if build_ccd_input_map was
not run first.
250 raise RuntimeError(
"Must run build_ccd_input_map before mask_warp_bbox")
253 bad_pixels = np.where(mask.array & bit_mask_value)
254 if len(bad_pixels[0]) == 0:
259 bad_ra, bad_dec = self.
_wcs.pixelToSkyArray(bad_pixels[1].astype(np.float64),
260 bad_pixels[0].astype(np.float64),
262 bad_hpix = hpg.angle_to_pixel(self.config.nside, bad_ra, bad_dec)
265 min_bad_hpix = bad_hpix.min()
266 bad_hpix_count = np.zeros(bad_hpix.max() - min_bad_hpix + 1, dtype=np.int32)
267 np.add.at(bad_hpix_count, bad_hpix - min_bad_hpix, 1)
272 pix_to_add, = np.where(bad_hpix_count > 0)
275 count_map_arr[primary] = np.clip(count_map_arr[primary], 0,
None)
277 count_map_arr[f
"v{visit}"] = np.clip(count_map_arr[f
"v{visit}"], 0,
None)
278 count_map_arr[f
"v{visit}"] += bad_hpix_count[pix_to_add]
283 """Use accumulated mask information to finalize the masking of
288 RuntimeError : Raised if build_ccd_input_map was
not run first.
291 raise RuntimeError(
"Must run build_ccd_input_map before finalize_ccd_input_map_mask.")
295 to_mask, = np.where(count_map_arr[f
"v{visit}"] > self.
_min_bad)
296 if to_mask.size == 0:
304 def _pixels_to_radec(self, wcs, pixels):
305 """Convert pixels to ra/dec positions using a wcs.
312 List of pixels to convert.
316 radec : `numpy.ndarray`
317 Nx2 array of ra/dec positions associated with pixels.
319 sph_pts = wcs.pixelToSky(pixels)
320 return np.array([(sph.getRa().asDegrees(), sph.getDec().asDegrees())
325 dimensions=(
"tract",
"band",
"skymap",),
326 defaultTemplates={
"coaddName":
"deep",
328 input_maps = pipeBase.connectionTypes.Input(
329 doc=
"Healsparse bit-wise coadd input maps",
330 name=
"{coaddName}Coadd_inputMap",
331 storageClass=
"HealSparseMap",
332 dimensions=(
"tract",
"patch",
"skymap",
"band"),
336 coadd_exposures = pipeBase.connectionTypes.Input(
337 doc=
"Coadded exposures associated with input_maps",
338 name=
"{coaddName}Coadd",
339 storageClass=
"ExposureF",
340 dimensions=(
"tract",
"patch",
"skymap",
"band"),
344 visit_summaries = pipeBase.connectionTypes.Input(
345 doc=
"Visit summary tables with aggregated statistics",
346 name=
"{calexpType}visitSummary",
347 storageClass=
"ExposureCatalog",
348 dimensions=(
"instrument",
"visit"),
352 sky_map = pipeBase.connectionTypes.Input(
353 doc=
"Input definition of geometry/bbox and projection/wcs for coadded exposures",
354 name=BaseSkyMap.SKYMAP_DATASET_TYPE_NAME,
355 storageClass=
"SkyMap",
356 dimensions=(
"skymap",),
364 for name
in BasePropertyMap.registry:
365 vars()[f
"{name}_map_min"] = pipeBase.connectionTypes.Output(
366 doc=f
"Minimum-value map of {name}",
367 name=f
"{{coaddName}}Coadd_{name}_map_min",
368 storageClass=
"HealSparseMap",
369 dimensions=(
"tract",
"skymap",
"band"),
371 vars()[f
"{name}_map_max"] = pipeBase.connectionTypes.Output(
372 doc=f
"Maximum-value map of {name}",
373 name=f
"{{coaddName}}Coadd_{name}_map_max",
374 storageClass=
"HealSparseMap",
375 dimensions=(
"tract",
"skymap",
"band"),
377 vars()[f
"{name}_map_mean"] = pipeBase.connectionTypes.Output(
378 doc=f
"Mean-value map of {name}",
379 name=f
"{{coaddName}}Coadd_{name}_map_mean",
380 storageClass=
"HealSparseMap",
381 dimensions=(
"tract",
"skymap",
"band"),
383 vars()[f
"{name}_map_weighted_mean"] = pipeBase.connectionTypes.Output(
384 doc=f
"Weighted mean-value map of {name}",
385 name=f
"{{coaddName}}Coadd_{name}_map_weighted_mean",
386 storageClass=
"HealSparseMap",
387 dimensions=(
"tract",
"skymap",
"band"),
389 vars()[f
"{name}_map_sum"] = pipeBase.connectionTypes.Output(
390 doc=f
"Sum-value map of {name}",
391 name=f
"{{coaddName}}Coadd_{name}_map_sum",
392 storageClass=
"HealSparseMap",
393 dimensions=(
"tract",
"skymap",
"band"),
396 def __init__(self, *, config=None):
397 super().__init__(config=config)
401 for name
in BasePropertyMap.registry:
402 if name
not in config.property_maps:
404 prop_config.do_min =
False
405 prop_config.do_max =
False
406 prop_config.do_mean =
False
407 prop_config.do_weighted_mean =
False
408 prop_config.do_sum =
False
410 prop_config = config.property_maps[name]
412 if not prop_config.do_min:
413 self.outputs.remove(f
"{name}_map_min")
414 if not prop_config.do_max:
415 self.outputs.remove(f
"{name}_map_max")
416 if not prop_config.do_mean:
417 self.outputs.remove(f
"{name}_map_mean")
418 if not prop_config.do_weighted_mean:
419 self.outputs.remove(f
"{name}_map_weighted_mean")
420 if not prop_config.do_sum:
421 self.outputs.remove(f
"{name}_map_sum")
424class HealSparsePropertyMapConfig(pipeBase.PipelineTaskConfig,
425 pipelineConnections=HealSparsePropertyMapConnections):
426 """Configuration parameters for HealSparsePropertyMapTask"""
427 property_maps = BasePropertyMap.registry.makeField(
429 default=[
"exposure_time",
440 doc=
"Property map computation objects",
443 def setDefaults(self):
444 self.property_maps[
"exposure_time"].do_sum =
True
445 self.property_maps[
"psf_size"].do_weighted_mean =
True
446 self.property_maps[
"psf_e1"].do_weighted_mean =
True
447 self.property_maps[
"psf_e2"].do_weighted_mean =
True
448 self.property_maps[
"psf_maglim"].do_weighted_mean =
True
449 self.property_maps[
"sky_noise"].do_weighted_mean =
True
450 self.property_maps[
"sky_background"].do_weighted_mean =
True
451 self.property_maps[
"dcr_dra"].do_weighted_mean =
True
452 self.property_maps[
"dcr_ddec"].do_weighted_mean =
True
453 self.property_maps[
"dcr_e1"].do_weighted_mean =
True
454 self.property_maps[
"dcr_e2"].do_weighted_mean =
True
457class HealSparsePropertyMapTask(pipeBase.PipelineTask):
458 """Task to compute Healsparse property maps.
460 This task will compute individual property maps (per tract, per
461 map type, per band). These maps cover the full coadd tract, and
462 are
not truncated to the inner tract region.
464 ConfigClass = HealSparsePropertyMapConfig
465 _DefaultName = "healSparsePropertyMapTask"
467 def __init__(self, **kwargs):
468 super().__init__(**kwargs)
470 for name, config, PropertyMapClass
in self.config.property_maps.apply():
471 self.property_maps[name] = PropertyMapClass(config, name)
474 def runQuantum(self, butlerQC, inputRefs, outputRefs):
475 inputs = butlerQC.get(inputRefs)
477 sky_map = inputs.pop(
"sky_map")
479 tract = butlerQC.quantum.dataId[
"tract"]
480 band = butlerQC.quantum.dataId[
"band"]
482 input_map_dict = {ref.dataId[
"patch"]: ref
for ref
in inputs[
"input_maps"]}
483 coadd_dict = {ref.dataId[
"patch"]: ref
for ref
in inputs[
"coadd_exposures"]}
485 visit_summary_dict = {ref.dataId[
"visit"]: ref.get()
486 for ref
in inputs[
"visit_summaries"]}
488 self.run(sky_map, tract, band, coadd_dict, input_map_dict, visit_summary_dict)
491 for name, property_map
in self.property_maps.items():
492 if property_map.config.do_min:
493 butlerQC.put(property_map.min_map,
494 getattr(outputRefs, f
"{name}_map_min"))
495 if property_map.config.do_max:
496 butlerQC.put(property_map.max_map,
497 getattr(outputRefs, f
"{name}_map_max"))
498 if property_map.config.do_mean:
499 butlerQC.put(property_map.mean_map,
500 getattr(outputRefs, f
"{name}_map_mean"))
501 if property_map.config.do_weighted_mean:
502 butlerQC.put(property_map.weighted_mean_map,
503 getattr(outputRefs, f
"{name}_map_weighted_mean"))
504 if property_map.config.do_sum:
505 butlerQC.put(property_map.sum_map,
506 getattr(outputRefs, f
"{name}_map_sum"))
508 def run(self, sky_map, tract, band, coadd_dict, input_map_dict, visit_summary_dict):
509 """Run the healsparse property task.
513 sky_map : Sky map object
517 Band name for logging.
518 coadd_dict : `dict` [`int`: `lsst.daf.butler.DeferredDatasetHandle`]
519 Dictionary of coadd exposure references. Keys are patch numbers.
520 input_map_dict : `dict` [`int`: `lsst.daf.butler.DeferredDatasetHandle`]
521 Dictionary of input map references. Keys are patch numbers.
523 Dictionary of visit summary tables. Keys are visit numbers.
527 RepeatableQuantumError
528 If visit_summary_dict
is missing any visits
or detectors found
in an
529 input map. This leads to an inconsistency between what
is in the coadd
530 (via the input map)
and the visit summary tables which contain data
533 tract_info = sky_map[tract]
535 tract_maps_initialized = False
537 for patch
in input_map_dict.keys():
538 self.log.info(
"Making maps for band %s, tract %d, patch %d.",
541 patch_info = tract_info[patch]
543 input_map = input_map_dict[patch].get()
544 coadd_photo_calib = coadd_dict[patch].get(component=
"photoCalib")
545 coadd_inputs = coadd_dict[patch].get(component=
"coaddInputs")
547 coadd_zeropoint = 2.5*np.log10(coadd_photo_calib.getInstFluxAtZeroMagnitude())
550 poly_vertices = patch_info.getInnerSkyPolygon(tract_info.getWcs()).getVertices()
551 patch_radec = self._vertices_to_radec(poly_vertices)
552 patch_poly = hsp.Polygon(ra=patch_radec[:, 0], dec=patch_radec[:, 1],
553 value=np.arange(input_map.wide_mask_maxbits))
554 patch_poly_map = patch_poly.get_map_like(input_map)
555 input_map = hsp.and_intersection([input_map, patch_poly_map])
557 if not tract_maps_initialized:
560 nside_coverage = self._compute_nside_coverage_tract(tract_info)
561 nside = input_map.nside_sparse
563 do_compute_approx_psf =
False
565 for property_map
in self.property_maps:
566 property_map.initialize_tract_maps(nside_coverage, nside)
567 if property_map.requires_psf:
568 do_compute_approx_psf =
True
570 tract_maps_initialized =
True
572 valid_pixels, vpix_ra, vpix_dec = input_map.valid_pixels_pos(return_pixels=
True)
575 if valid_pixels.size == 0:
579 for property_map
in self.property_maps:
580 property_map.initialize_values(valid_pixels.size)
581 property_map.zeropoint = coadd_zeropoint
584 total_weights = np.zeros(valid_pixels.size)
585 total_inputs = np.zeros(valid_pixels.size, dtype=np.int32)
587 for bit, ccd_row
in enumerate(coadd_inputs.ccds):
589 inmap, = np.where(input_map.check_bits_pix(valid_pixels, [bit]))
596 visit = ccd_row[
"visit"]
597 detector_id = ccd_row[
"ccd"]
598 weight = ccd_row[
"weight"]
600 x, y = ccd_row.getWcs().skyToPixelArray(vpix_ra[inmap], vpix_dec[inmap], degrees=
True)
601 scalings = self._compute_calib_scale(ccd_row, x, y)
603 if do_compute_approx_psf:
604 psf_array = compute_approx_psf_size_and_shape(ccd_row, vpix_ra[inmap], vpix_dec[inmap])
608 total_weights[inmap] += weight
609 total_inputs[inmap] += 1
612 if visit
not in visit_summary_dict:
613 msg = f
"Visit {visit} not found in visit_summaries."
614 raise pipeBase.RepeatableQuantumError(msg)
615 row = visit_summary_dict[visit].find(detector_id)
617 msg = f
"Visit {visit} / detector_id {detector_id} not found in visit_summaries."
618 raise pipeBase.RepeatableQuantumError(msg)
621 for property_map
in self.property_maps:
622 property_map.accumulate_values(inmap,
631 for property_map
in self.property_maps:
632 property_map.finalize_mean_values(total_weights, total_inputs)
633 property_map.set_map_values(valid_pixels)
635 def _compute_calib_scale(self, ccd_row, x, y):
636 """Compute calibration scaling values.
641 Exposure metadata for a given detector exposure.
643 Array of x positions.
645 Array of y positions.
649 calib_scale : `np.ndarray`
650 Array of calibration scale values.
652 photo_calib = ccd_row.getPhotoCalib()
653 bf = photo_calib.computeScaledCalibration()
654 if bf.getBBox() == ccd_row.getBBox():
656 calib_scale = photo_calib.getCalibrationMean()*bf.evaluate(x, y)
659 calib_scale = photo_calib.getCalibrationMean()
663 def _vertices_to_radec(self, vertices):
664 """Convert polygon vertices to ra/dec.
669 Vertices for bounding polygon.
673 radec : `numpy.ndarray`
674 Nx2 array of ra/dec positions (
in degrees) associated
with vertices.
677 radec = np.array([(x.getLon().asDegrees(), x.getLat().asDegrees())
for
681 def _compute_nside_coverage_tract(self, tract_info):
682 """Compute the optimal coverage nside for a tract.
687 Tract information object.
691 nside_coverage : `int`
692 Optimal coverage nside for a tract map.
694 num_patches = tract_info.getNumPatches()
697 patch_info = tract_info.getPatchInfo(0)
698 vertices = patch_info.getInnerSkyPolygon(tract_info.getWcs()).getVertices()
699 radec = self._vertices_to_radec(vertices)
700 delta_ra = np.max(radec[:, 0]) - np.min(radec[:, 0])
701 delta_dec = np.max(radec[:, 1]) - np.min(radec[:, 1])
702 patch_area = delta_ra*delta_dec*np.cos(np.deg2rad(np.mean(radec[:, 1])))
704 tract_area = num_patches[0]*num_patches[1]*patch_area
706 nside_coverage_tract = 32
707 while hpg.nside_to_pixel_area(nside_coverage_tract, degrees=
True) > tract_area:
708 nside_coverage_tract = 2*nside_coverage_tract
711 nside_coverage_tract = int(np.clip(nside_coverage_tract/2, 32, 128))
713 return nside_coverage_tract
716class ConsolidateHealSparsePropertyMapConnections(pipeBase.PipelineTaskConnections,
717 dimensions=(
"band",
"skymap",),
718 defaultTemplates={
"coaddName":
"deep"}):
719 sky_map = pipeBase.connectionTypes.Input(
720 doc=
"Input definition of geometry/bbox and projection/wcs for coadded exposures",
721 name=BaseSkyMap.SKYMAP_DATASET_TYPE_NAME,
722 storageClass=
"SkyMap",
723 dimensions=(
"skymap",),
731 for name
in BasePropertyMap.registry:
732 vars()[f
"{name}_map_min"] = pipeBase.connectionTypes.Input(
733 doc=f
"Minimum-value map of {name}",
734 name=f
"{{coaddName}}Coadd_{name}_map_min",
735 storageClass=
"HealSparseMap",
736 dimensions=(
"tract",
"skymap",
"band"),
740 vars()[f
"{name}_consolidated_map_min"] = pipeBase.connectionTypes.Output(
741 doc=f
"Minumum-value map of {name}",
742 name=f
"{{coaddName}}Coadd_{name}_consolidated_map_min",
743 storageClass=
"HealSparseMap",
744 dimensions=(
"skymap",
"band"),
746 vars()[f
"{name}_map_max"] = pipeBase.connectionTypes.Input(
747 doc=f
"Maximum-value map of {name}",
748 name=f
"{{coaddName}}Coadd_{name}_map_max",
749 storageClass=
"HealSparseMap",
750 dimensions=(
"tract",
"skymap",
"band"),
754 vars()[f
"{name}_consolidated_map_max"] = pipeBase.connectionTypes.Output(
755 doc=f
"Minumum-value map of {name}",
756 name=f
"{{coaddName}}Coadd_{name}_consolidated_map_max",
757 storageClass=
"HealSparseMap",
758 dimensions=(
"skymap",
"band"),
760 vars()[f
"{name}_map_mean"] = pipeBase.connectionTypes.Input(
761 doc=f
"Mean-value map of {name}",
762 name=f
"{{coaddName}}Coadd_{name}_map_mean",
763 storageClass=
"HealSparseMap",
764 dimensions=(
"tract",
"skymap",
"band"),
768 vars()[f
"{name}_consolidated_map_mean"] = pipeBase.connectionTypes.Output(
769 doc=f
"Minumum-value map of {name}",
770 name=f
"{{coaddName}}Coadd_{name}_consolidated_map_mean",
771 storageClass=
"HealSparseMap",
772 dimensions=(
"skymap",
"band"),
774 vars()[f
"{name}_map_weighted_mean"] = pipeBase.connectionTypes.Input(
775 doc=f
"Weighted mean-value map of {name}",
776 name=f
"{{coaddName}}Coadd_{name}_map_weighted_mean",
777 storageClass=
"HealSparseMap",
778 dimensions=(
"tract",
"skymap",
"band"),
782 vars()[f
"{name}_consolidated_map_weighted_mean"] = pipeBase.connectionTypes.Output(
783 doc=f
"Minumum-value map of {name}",
784 name=f
"{{coaddName}}Coadd_{name}_consolidated_map_weighted_mean",
785 storageClass=
"HealSparseMap",
786 dimensions=(
"skymap",
"band"),
788 vars()[f
"{name}_map_sum"] = pipeBase.connectionTypes.Input(
789 doc=f
"Sum-value map of {name}",
790 name=f
"{{coaddName}}Coadd_{name}_map_sum",
791 storageClass=
"HealSparseMap",
792 dimensions=(
"tract",
"skymap",
"band"),
796 vars()[f
"{name}_consolidated_map_sum"] = pipeBase.connectionTypes.Output(
797 doc=f
"Minumum-value map of {name}",
798 name=f
"{{coaddName}}Coadd_{name}_consolidated_map_sum",
799 storageClass=
"HealSparseMap",
800 dimensions=(
"skymap",
"band"),
803 def __init__(self, *, config=None):
804 super().__init__(config=config)
808 for name
in BasePropertyMap.registry:
809 if name
not in config.property_maps:
811 prop_config.do_min =
False
812 prop_config.do_max =
False
813 prop_config.do_mean =
False
814 prop_config.do_weighted_mean =
False
815 prop_config.do_sum =
False
817 prop_config = config.property_maps[name]
819 if not prop_config.do_min:
820 self.inputs.remove(f
"{name}_map_min")
821 self.outputs.remove(f
"{name}_consolidated_map_min")
822 if not prop_config.do_max:
823 self.inputs.remove(f
"{name}_map_max")
824 self.outputs.remove(f
"{name}_consolidated_map_max")
825 if not prop_config.do_mean:
826 self.inputs.remove(f
"{name}_map_mean")
827 self.outputs.remove(f
"{name}_consolidated_map_mean")
828 if not prop_config.do_weighted_mean:
829 self.inputs.remove(f
"{name}_map_weighted_mean")
830 self.outputs.remove(f
"{name}_consolidated_map_weighted_mean")
831 if not prop_config.do_sum:
832 self.inputs.remove(f
"{name}_map_sum")
833 self.outputs.remove(f
"{name}_consolidated_map_sum")
836class ConsolidateHealSparsePropertyMapConfig(pipeBase.PipelineTaskConfig,
837 pipelineConnections=ConsolidateHealSparsePropertyMapConnections):
838 """Configuration parameters for ConsolidateHealSparsePropertyMapTask"""
839 property_maps = BasePropertyMap.registry.makeField(
841 default=[
"exposure_time",
852 doc=
"Property map computation objects",
854 nside_coverage = pexConfig.Field(
855 doc=
"Consolidated HealSparse coverage map nside. Must be power of 2.",
858 check=_is_power_of_two,
861 def setDefaults(self):
862 self.property_maps[
"exposure_time"].do_sum =
True
863 self.property_maps[
"psf_size"].do_weighted_mean =
True
864 self.property_maps[
"psf_e1"].do_weighted_mean =
True
865 self.property_maps[
"psf_e2"].do_weighted_mean =
True
866 self.property_maps[
"psf_maglim"].do_weighted_mean =
True
867 self.property_maps[
"sky_noise"].do_weighted_mean =
True
868 self.property_maps[
"sky_background"].do_weighted_mean =
True
869 self.property_maps[
"dcr_dra"].do_weighted_mean =
True
870 self.property_maps[
"dcr_ddec"].do_weighted_mean =
True
871 self.property_maps[
"dcr_e1"].do_weighted_mean =
True
872 self.property_maps[
"dcr_e2"].do_weighted_mean =
True
875class ConsolidateHealSparsePropertyMapTask(pipeBase.PipelineTask):
876 """Task to consolidate HealSparse property maps.
878 This task will take all the individual tract-based maps (per map type,
879 per band) and consolidate them into one survey-wide map (per map type,
880 per band). Each tract map
is truncated to its inner region before
883 ConfigClass = ConsolidateHealSparsePropertyMapConfig
884 _DefaultName = "consolidateHealSparsePropertyMapTask"
886 def __init__(self, **kwargs):
887 super().__init__(**kwargs)
889 for name, config, PropertyMapClass
in self.config.property_maps.apply():
890 self.property_maps[name] = PropertyMapClass(config, name)
893 def runQuantum(self, butlerQC, inputRefs, outputRefs):
894 inputs = butlerQC.get(inputRefs)
896 sky_map = inputs.pop(
"sky_map")
899 for name
in self.config.property_maps.names:
900 for type_
in [
'min',
'max',
'mean',
'weighted_mean',
'sum']:
901 map_type = f
"{name}_map_{type_}"
902 if map_type
in inputs:
903 input_refs = {ref.dataId[
'tract']: ref
904 for ref
in inputs[map_type]}
905 consolidated_map = self.consolidate_map(sky_map, input_refs)
906 butlerQC.put(consolidated_map,
907 getattr(outputRefs, f
"{name}_consolidated_map_{type_}"))
909 def consolidate_map(self, sky_map, input_refs):
910 """Consolidate the healsparse property maps.
914 sky_map : Sky map object
915 input_refs : `dict` [`int`: `lsst.daf.butler.DeferredDatasetHandle`]
916 Dictionary of tract_id mapping to dataref.
920 consolidated_map : `healsparse.HealSparseMap`
921 Consolidated HealSparse map.
926 nside_coverage_inputs =
None
927 for tract_id
in input_refs:
928 cov = input_refs[tract_id].get(component=
'coverage')
930 cov_mask = cov.coverage_mask
931 nside_coverage_inputs = cov.nside_coverage
933 cov_mask |= cov.coverage_mask
935 cov_pix_inputs, = np.where(cov_mask)
938 if nside_coverage_inputs == self.config.nside_coverage:
939 cov_pix = cov_pix_inputs
940 elif nside_coverage_inputs > self.config.nside_coverage:
943 bit_shift = hsp.utils._compute_bitshift(self.config.nside_coverage,
944 nside_coverage_inputs)
945 cov_pix = np.right_shift(cov_pix_inputs, bit_shift)
949 bit_shift = hsp.utils._compute_bitshift(nside_coverage_inputs,
950 self.config.nside_coverage)
951 cov_pix = np.left_shift(cov_pix_inputs, bit_shift)
954 consolidated_map =
None
955 for tract_id
in input_refs:
956 input_map = input_refs[tract_id].get()
957 if consolidated_map
is None:
958 consolidated_map = hsp.HealSparseMap.make_empty(
959 self.config.nside_coverage,
960 input_map.nside_sparse,
962 sentinel=input_map._sentinel,
966 vpix, ra, dec = input_map.valid_pixels_pos(return_pixels=
True)
967 vpix_tract_ids = sky_map.findTractIdArray(ra, dec, degrees=
True)
969 in_tract = (vpix_tract_ids == tract_id)
971 consolidated_map[vpix[in_tract]] = input_map[vpix[in_tract]]
973 return consolidated_map