lsst.pipe.tasks  21.0.0-85-g296c1d01+e289388810
healSparseMapping.py
Go to the documentation of this file.
1 #
2 # LSST Data Management System
3 # Copyright 2008-2021 AURA/LSST.
4 #
5 # This product includes software developed by the
6 # LSST Project (http://www.lsst.org/).
7 #
8 # This program is free software: you can redistribute it and/or modify
9 # it under the terms of the GNU General Public License as published by
10 # the Free Software Foundation, either version 3 of the License, or
11 # (at your option) any later version.
12 #
13 # This program is distributed in the hope that it will be useful,
14 # but WITHOUT ANY WARRANTY; without even the implied warranty of
15 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 # GNU General Public License for more details.
17 #
18 # You should have received a copy of the LSST License Statement and
19 # the GNU General Public License along with this program. If not,
20 # see <http://www.lsstcorp.org/LegalNotices/>.
21 #
22 from collections import defaultdict
23 import numbers
24 import numpy as np
25 import healpy as hp
26 import healsparse as hsp
27 
28 import lsst.pex.config as pexConfig
29 import lsst.pipe.base as pipeBase
30 import lsst.geom
31 import lsst.afw.geom as afwGeom
32 from lsst.daf.butler import Formatter
33 
34 
35 __all__ = ["HealSparseInputMapTask", "HealSparseInputMapConfig",
36  "HealSparseMapFormatter"]
37 
38 
39 class HealSparseMapFormatter(Formatter):
40  """Interface for reading and writing healsparse.HealSparseMap files
41  """
42  unsupportedParameters = frozenset()
43  supportedExtensions = frozenset({".hsp", ".fit", ".fits"})
44  extension = '.hsp'
45 
46  def read(self, component=None):
47  # Docstring inherited from Formatter.read.
48  path = self.fileDescriptor.location.path
49 
50  if component == 'coverage':
51  try:
52  data = hsp.HealSparseCoverage.read(path)
53  except (OSError, RuntimeError):
54  raise ValueError(f"Unable to read healsparse map with URI {self.fileDescriptor.location.uri}")
55 
56  return data
57 
58  if self.fileDescriptor.parameters is None:
59  pixels = None
60  degrade_nside = None
61  else:
62  pixels = self.fileDescriptor.parameters.get('pixels', None)
63  degrade_nside = self.fileDescriptor.parameters.get('degrade_nside', None)
64  try:
65  data = hsp.HealSparseMap.read(path, pixels=pixels, degrade_nside=degrade_nside)
66  except (OSError, RuntimeError):
67  raise ValueError(f"Unable to read healsparse map with URI {self.fileDescriptor.location.uri}")
68 
69  return data
70 
71  def write(self, inMemoryDataset):
72  # Docstring inherited from Formatter.write.
73  # Update the location with the formatter-preferred file extension
74  self.fileDescriptor.location.updateExtension(self.extensionextension)
75  inMemoryDataset.write(self.fileDescriptor.location.path, clobber=True)
76 
77 
78 def _is_power_of_two(value):
79  """Check that value is a power of two.
80 
81  Parameters
82  ----------
83  value : `int`
84  Value to check.
85 
86  Returns
87  -------
88  is_power_of_two : `bool`
89  True if value is a power of two; False otherwise, or
90  if value is not an integer.
91  """
92  if not isinstance(value, numbers.Integral):
93  return False
94 
95  # See https://stackoverflow.com/questions/57025836
96  # Every power of 2 has exactly 1 bit set to 1; subtracting
97  # 1 therefore flips every preceding bit. If you and that
98  # together with the original value it must be 0.
99  return (value & (value - 1) == 0) and value != 0
100 
101 
102 class HealSparseInputMapConfig(pexConfig.Config):
103  """Configuration parameters for HealSparseInputMapTask"""
104  nside = pexConfig.Field(
105  doc="Mapping healpix nside. Must be power of 2.",
106  dtype=int,
107  default=32768,
108  check=_is_power_of_two,
109  )
110  nside_coverage = pexConfig.Field(
111  doc="HealSparse coverage map nside. Must be power of 2.",
112  dtype=int,
113  default=256,
114  check=_is_power_of_two,
115  )
116  bad_mask_min_coverage = pexConfig.Field(
117  doc=("Minimum area fraction of a map healpixel pixel that must be "
118  "covered by bad pixels to be removed from the input map. "
119  "This is approximate."),
120  dtype=float,
121  default=0.5,
122  )
123 
124 
125 class HealSparseInputMapTask(pipeBase.Task):
126  """Task for making a HealSparse input map."""
127 
128  ConfigClass = HealSparseInputMapConfig
129  _DefaultName = "healSparseInputMap"
130 
131  def __init__(self, **kwargs):
132  pipeBase.Task.__init__(self, **kwargs)
133 
134  self.ccd_input_mapccd_input_map = None
135 
136  def build_ccd_input_map(self, bbox, wcs, ccds):
137  """Build a map from ccd valid polygons or bounding boxes.
138 
139  Parameters
140  ----------
141  bbox : `lsst.geom.Box2I`
142  Bounding box for region to build input map.
143  wcs : `lsst.afw.geom.SkyWcs`
144  WCS object for region to build input map.
145  ccds : `lsst.afw.table.ExposureCatalog`
146  Exposure catalog with ccd data from coadd inputs.
147  """
148  self.ccd_input_mapccd_input_map = hsp.HealSparseMap.make_empty(nside_coverage=self.config.nside_coverage,
149  nside_sparse=self.config.nside,
150  dtype=hsp.WIDE_MASK,
151  wide_mask_maxbits=len(ccds))
152  self._wcs_wcs = wcs
153  self._bbox_bbox = bbox
154  self._ccds_ccds = ccds
155 
156  pixel_scale = wcs.getPixelScale().asArcseconds()
157  hpix_area_arcsec2 = hp.nside2pixarea(self.config.nside, degrees=True)*(3600.**2.)
158  self._min_bad_min_bad = self.config.bad_mask_min_coverage*hpix_area_arcsec2/(pixel_scale**2.)
159 
160  metadata = {}
161  self._bits_per_visit_ccd_bits_per_visit_ccd = {}
162  self._bits_per_visit_bits_per_visit = defaultdict(list)
163  for bit, ccd in enumerate(ccds):
164  metadata[f'B{bit:04d}CCD'] = ccd['ccd']
165  metadata[f'B{bit:04d}VIS'] = ccd['visit']
166  metadata[f'B{bit:04d}WT'] = ccd['weight']
167 
168  self._bits_per_visit_ccd_bits_per_visit_ccd[(ccd['visit'], ccd['ccd'])] = bit
169  self._bits_per_visit_bits_per_visit[ccd['visit']].append(bit)
170 
171  ccd_poly = ccd.getValidPolygon()
172  if ccd_poly is None:
173  ccd_poly = afwGeom.Polygon(lsst.geom.Box2D(ccd.getBBox()))
174  # Detectors need to be rendered with their own wcs.
175  ccd_poly_radec = self._pixels_to_radec_pixels_to_radec(ccd.getWcs(), ccd_poly.convexHull().getVertices())
176 
177  # Create a ccd healsparse polygon
178  poly = hsp.Polygon(ra=ccd_poly_radec[: -1, 0],
179  dec=ccd_poly_radec[: -1, 1],
180  value=[bit])
181  self.ccd_input_mapccd_input_map.set_bits_pix(poly.get_pixels(nside=self.ccd_input_mapccd_input_map.nside_sparse),
182  [bit])
183 
184  # Cut down to the overall bounding box with associated wcs.
185  bbox_afw_poly = afwGeom.Polygon(lsst.geom.Box2D(bbox))
186  bbox_poly_radec = self._pixels_to_radec_pixels_to_radec(self._wcs_wcs,
187  bbox_afw_poly.convexHull().getVertices())
188  bbox_poly = hsp.Polygon(ra=bbox_poly_radec[: -1, 0], dec=bbox_poly_radec[: -1, 1],
189  value=np.arange(self.ccd_input_mapccd_input_map.wide_mask_maxbits))
190  bbox_poly_map = bbox_poly.get_map_like(self.ccd_input_mapccd_input_map)
191  self.ccd_input_mapccd_input_map = hsp.and_intersection([self.ccd_input_mapccd_input_map, bbox_poly_map])
192  self.ccd_input_mapccd_input_map.metadata = metadata
193 
194  # Create a temporary map to hold the count of bad pixels in each healpix pixel
195  self._ccd_input_pixels_ccd_input_pixels = self.ccd_input_mapccd_input_map.valid_pixels
196 
197  dtype = [(f'v{visit}', 'i4') for visit in self._bits_per_visit_bits_per_visit.keys()]
198 
199  cov = self.config.nside_coverage
200  ns = self.config.nside
201  self._ccd_input_bad_count_map_ccd_input_bad_count_map = hsp.HealSparseMap.make_empty(nside_coverage=cov,
202  nside_sparse=ns,
203  dtype=dtype,
204  primary=dtype[0][0])
205  # Don't set input bad map if there are no ccds which overlap the bbox.
206  if len(self._ccd_input_pixels_ccd_input_pixels) > 0:
207  self._ccd_input_bad_count_map_ccd_input_bad_count_map[self._ccd_input_pixels_ccd_input_pixels] = np.zeros(1, dtype=dtype)
208 
209  def mask_warp_bbox(self, bbox, visit, mask, bit_mask_value):
210  """Mask a subregion from a visit.
211  This must be run after build_ccd_input_map initializes
212  the overall map.
213 
214  Parameters
215  ----------
216  bbox : `lsst.geom.Box2I`
217  Bounding box from region to mask.
218  visit : `int`
219  Visit number corresponding to warp with mask.
220  mask : `lsst.afw.image.MaskX`
221  Mask plane from warp exposure.
222  bit_mask_value : `int`
223  Bit mask to check for bad pixels.
224 
225  Raises
226  ------
227  RuntimeError : Raised if build_ccd_input_map was not run first.
228  """
229  if self.ccd_input_mapccd_input_map is None:
230  raise RuntimeError("Must run build_ccd_input_map before mask_warp_bbox")
231 
232  # Find the bad pixels and convert to healpix
233  bad_pixels = np.where(mask.array & bit_mask_value)
234  if len(bad_pixels[0]) == 0:
235  # No bad pixels
236  return
237 
238  # Bad pixels come from warps which use the overall wcs.
239  bad_ra, bad_dec = self._wcs_wcs.pixelToSkyArray(bad_pixels[1].astype(np.float64),
240  bad_pixels[0].astype(np.float64),
241  degrees=True)
242  bad_hpix = hp.ang2pix(self.config.nside, bad_ra, bad_dec,
243  lonlat=True, nest=True)
244 
245  # Count the number of bad image pixels in each healpix pixel
246  min_bad_hpix = bad_hpix.min()
247  bad_hpix_count = np.zeros(bad_hpix.max() - min_bad_hpix + 1, dtype=np.int32)
248  np.add.at(bad_hpix_count, bad_hpix - min_bad_hpix, 1)
249 
250  # Add these to the accumulator map.
251  # We need to make sure that the "primary" array has valid values for
252  # this pixel to be registered in the accumulator map.
253  pix_to_add, = np.where(bad_hpix_count > 0)
254  count_map_arr = self._ccd_input_bad_count_map_ccd_input_bad_count_map[min_bad_hpix + pix_to_add]
255  primary = self._ccd_input_bad_count_map_ccd_input_bad_count_map.primary
256  count_map_arr[primary] = np.clip(count_map_arr[primary], 0, None)
257 
258  count_map_arr[f'v{visit}'] = np.clip(count_map_arr[f'v{visit}'], 0, None)
259  count_map_arr[f'v{visit}'] += bad_hpix_count[pix_to_add]
260 
261  self._ccd_input_bad_count_map_ccd_input_bad_count_map[min_bad_hpix + pix_to_add] = count_map_arr
262 
264  """Use accumulated mask information to finalize the masking of
265  ccd_input_map.
266 
267  Raises
268  ------
269  RuntimeError : Raised if build_ccd_input_map was not run first.
270  """
271  if self.ccd_input_mapccd_input_map is None:
272  raise RuntimeError("Must run build_ccd_input_map before finalize_ccd_input_map_mask.")
273 
274  count_map_arr = self._ccd_input_bad_count_map_ccd_input_bad_count_map[self._ccd_input_pixels_ccd_input_pixels]
275  for visit in self._bits_per_visit_bits_per_visit:
276  to_mask, = np.where(count_map_arr[f'v{visit}'] > self._min_bad_min_bad)
277  if to_mask.size == 0:
278  continue
279  self.ccd_input_mapccd_input_map.clear_bits_pix(self._ccd_input_pixels_ccd_input_pixels[to_mask],
280  self._bits_per_visit_bits_per_visit[visit])
281 
282  # Clear memory
283  self._ccd_input_bad_count_map_ccd_input_bad_count_map = None
284 
285  def _pixels_to_radec(self, wcs, pixels):
286  """Convert pixels to ra/dec positions using a wcs.
287 
288  Parameters
289  ----------
290  wcs : `lsst.afw.geom.SkyWcs`
291  WCS object.
292  pixels : `list` [`lsst.geom.Point2D`]
293  List of pixels to convert.
294 
295  Returns
296  -------
297  radec : `numpy.ndarray`
298  Nx2 array of ra/dec positions associated with pixels.
299  """
300  sph_pts = wcs.pixelToSky(pixels)
301  return np.array([(sph.getRa().asDegrees(), sph.getDec().asDegrees())
302  for sph in sph_pts])
def mask_warp_bbox(self, bbox, visit, mask, bit_mask_value)