lsst.meas.algorithms  17.0.1-14-gd9e67d69+2
defects.py
Go to the documentation of this file.
1 #
2 # LSST Data Management System
3 #
4 # Copyright 2008-2017 AURA/LSST.
5 #
6 # This product includes software developed by the
7 # LSST Project (http://www.lsst.org/).
8 #
9 # This program is free software: you can redistribute it and/or modify
10 # it under the terms of the GNU General Public License as published by
11 # the Free Software Foundation, either version 3 of the License, or
12 # (at your option) any later version.
13 #
14 # This program is distributed in the hope that it will be useful,
15 # but WITHOUT ANY WARRANTY; without even the implied warranty of
16 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 # GNU General Public License for more details.
18 #
19 # You should have received a copy of the LSST License Statement and
20 # the GNU General Public License along with this program. If not,
21 # see <https://www.lsstcorp.org/LegalNotices/>.
22 #
23 """Support for image defects"""
24 
25 __all__ = ("Defects",)
26 
27 import logging
28 import itertools
29 import collections.abc
30 from deprecated.sphinx import deprecated
31 import numpy as np
32 import copy
33 import datetime
34 import math
35 import numbers
36 import os.path
37 import astropy.table
38 
39 import lsst.geom
40 import lsst.pex.policy as policy
41 import lsst.afw.table
42 import lsst.afw.detection
43 import lsst.afw.image
44 import lsst.afw.geom
45 from lsst.daf.base import PropertyList
46 
47 from . import Defect
48 
49 log = logging.getLogger(__name__)
50 
51 SCHEMA_NAME_KEY = "DEFECTS_SCHEMA"
52 SCHEMA_VERSION_KEY = "DEFECTS_SCHEMA_VERSION"
53 
54 
55 @deprecated(reason="Policy defect files no longer supported (will be removed after v18)",
56  category=FutureWarning)
57 def policyToBadRegionList(policyFile):
58  """Given a Policy file describing a CCD's bad pixels, return a vector of BadRegion::Ptr"""
59 
60  badPixelsPolicy = policy.Policy.createPolicy(policyFile)
61  badPixels = []
62 
63  if badPixelsPolicy.exists("Defects"):
64  d = badPixelsPolicy.getArray("Defects")
65  for reg in d:
66  x0 = reg.get("x0")
67  width = reg.get("width")
68  if not width:
69  x1 = reg.get("x1")
70  width = x1 - x0 - 1
71 
72  y0 = reg.get("y0")
73  if reg.exists("height"):
74  height = reg.get("height")
75  else:
76  y1 = reg.get("y1")
77  height = y1 - y0 - 1
78 
79  bbox = lsst.geom.BoxI(lsst.geom.PointI(x0, y0), lsst.geom.ExtentI(width, height))
80  badPixels.append(Defect(bbox))
81 
82  return badPixels
83 
84 
85 class Defects(collections.abc.MutableSequence):
86  """Collection of `lsst.meas.algorithms.Defect`.
87 
88  Parameters
89  ----------
90  defectList : iterable of `lsst.meas.algorithms.Defect`
91  or `lsst.geom.BoxI`, optional
92  Collections of defects to apply to the image.
93  """
94 
95  _OBSTYPE = "defects"
96  """The calibration type used for ingest."""
97 
98  def __init__(self, defectList=None, metadata=None):
99  self._defects = []
100 
101  if metadata is not None:
102  self._metadata = metadata
103  else:
104  self.setMetadata()
105 
106  if defectList is None:
107  return
108 
109  # Ensure that type checking
110  for d in defectList:
111  self.append(d)
112 
113  def _check_value(self, value):
114  """Check that the supplied value is a `~lsst.meas.algorithms.Defect`
115  or can be converted to one.
116 
117  Parameters
118  ----------
119  value : `object`
120  Value to check.
121 
122  Returns
123  -------
124  new : `~lsst.meas.algorithms.Defect`
125  Either the supplied value or a new object derived from it.
126 
127  Raises
128  ------
129  ValueError
130  Raised if the supplied value can not be converted to
131  `~lsst.meas.algorithms.Defect`
132  """
133  if isinstance(value, Defect):
134  pass
135  elif isinstance(value, lsst.geom.BoxI):
136  value = Defect(value)
137  elif isinstance(value, lsst.geom.PointI):
138  value = Defect(lsst.geom.Box2I(value, lsst.geom.Extent2I(1, 1)))
139  elif isinstance(value, lsst.afw.image.DefectBase):
140  value = Defect(value.getBBox())
141  else:
142  raise ValueError(f"Defects must be of type Defect, BoxI, or PointI, not '{value!r}'")
143  return value
144 
145  def __len__(self):
146  return len(self._defects)
147 
148  def __getitem__(self, index):
149  return self._defects[index]
150 
151  def __setitem__(self, index, value):
152  """Can be given a `~lsst.meas.algorithms.Defect` or a `lsst.geom.BoxI`
153  """
154  self._defects[index] = self._check_value(value)
155 
156  def __iter__(self):
157  return iter(self._defects)
158 
159  def __delitem__(self, index):
160  del self._defects[index]
161 
162  def __eq__(self, other):
163  """Compare if two `Defects` are equal.
164 
165  Two `Defects` are equal if their bounding boxes are equal and in
166  the same order. Metadata content is ignored.
167  """
168  if not isinstance(other, self.__class__):
169  return False
170 
171  # Assume equal if bounding boxes are equal
172  for d1, d2 in zip(self, other):
173  if d1.getBBox() != d2.getBBox():
174  return False
175 
176  return True
177 
178  def __str__(self):
179  return "Defects(" + ",".join(str(d.getBBox()) for d in self) + ")"
180 
181  def insert(self, index, value):
182  self._defects.insert(index, self._check_value(value))
183 
184  def getMetadata(self):
185  """Retrieve metadata associated with these `Defects`.
186 
187  Returns
188  -------
189  meta : `lsst.daf.base.PropertyList`
190  Metadata. The returned `~lsst.daf.base.PropertyList` can be
191  modified by the caller and the changes will be written to
192  external files.
193  """
194  return self._metadata
195 
196  def setMetadata(self, metadata=None):
197  """Store a copy of the supplied metadata with the defects.
198 
199  Parameters
200  ----------
201  metadata : `lsst.daf.base.PropertyList`, optional
202  Metadata to associate with the defects. Will be copied and
203  overwrite existing metadata. If not supplied the existing
204  metadata will be reset.
205  """
206  if metadata is None:
207  self._metadata = PropertyList()
208  else:
209  self._metadata = copy.copy(metadata)
210 
211  # Ensure that we have the obs type required by calibration ingest
212  self._metadata["OBSTYPE"] = self._OBSTYPE
213 
214  def copy(self):
215  """Copy the defects to a new list, creating new defects from the
216  bounding boxes.
217 
218  Returns
219  -------
220  new : `Defects`
221  New list with new `Defect` entries.
222 
223  Notes
224  -----
225  This is not a shallow copy in that new `Defect` instances are
226  created from the original bounding boxes. It's also not a deep
227  copy since the bounding boxes are not recreated.
228  """
229  return self.__class__(d.getBBox() for d in self)
230 
231  def transpose(self):
232  """Make a transposed copy of this defect list.
233 
234  Returns
235  -------
236  retDefectList : `Defects`
237  Transposed list of defects.
238  """
239  retDefectList = self.__class__()
240  for defect in self:
241  bbox = defect.getBBox()
242  dimensions = bbox.getDimensions()
243  nbbox = lsst.geom.Box2I(lsst.geom.Point2I(bbox.getMinY(), bbox.getMinX()),
244  lsst.geom.Extent2I(dimensions[1], dimensions[0]))
245  retDefectList.append(nbbox)
246  return retDefectList
247 
248  def maskPixels(self, maskedImage, maskName="BAD"):
249  """Set mask plane based on these defects.
250 
251  Parameters
252  ----------
253  maskedImage : `lsst.afw.image.MaskedImage`
254  Image to process. Only the mask plane is updated.
255  maskName : str, optional
256  Mask plane name to use.
257  """
258  # mask bad pixels
259  mask = maskedImage.getMask()
260  bitmask = mask.getPlaneBitMask(maskName)
261  for defect in self:
262  bbox = defect.getBBox()
263  lsst.afw.geom.SpanSet(bbox).clippedTo(mask.getBBox()).setMask(mask, bitmask)
264 
265  def toFitsRegionTable(self):
266  """Convert defect list to `~lsst.afw.table.BaseCatalog` using the
267  FITS region standard.
268 
269  Returns
270  -------
271  table : `lsst.afw.table.BaseCatalog`
272  Defects in tabular form.
273 
274  Notes
275  -----
276  The table created uses the
277  `FITS regions <https://fits.gsfc.nasa.gov/registry/region.html>`_
278  definition tabular format. The ``X`` and ``Y`` coordinates are
279  converted to FITS Physical coordinates that have origin pixel (1, 1)
280  rather than the (0, 0) used in LSST software.
281  """
282  schema = lsst.afw.table.Schema()
283  x = schema.addField("X", type="D", units="pix", doc="X coordinate of center of shape")
284  y = schema.addField("Y", type="D", units="pix", doc="Y coordinate of center of shape")
285  shape = schema.addField("SHAPE", type="String", size=16, doc="Shape defined by these values")
286  r = schema.addField("R", type="ArrayD", size=2, units="pix", doc="Extents")
287  rotang = schema.addField("ROTANG", type="D", units="deg", doc="Rotation angle")
288  component = schema.addField("COMPONENT", type="I", doc="Index of this region")
289  table = lsst.afw.table.BaseCatalog(schema)
290  table.resize(len(self._defects))
291 
292  for i, defect in enumerate(self._defects):
293  box = defect.getBBox()
294  # Correct for the FITS 1-based offset
295  table[i][x] = box.getCenterX() + 1.0
296  table[i][y] = box.getCenterY() + 1.0
297  width = box.getWidth()
298  height = box.getHeight()
299 
300  if width == 1 and height == 1:
301  # Call this a point
302  shapeType = "POINT"
303  else:
304  shapeType = "BOX"
305  table[i][shape] = shapeType
306  table[i][r] = np.array([width, height], dtype=np.float64)
307  table[i][rotang] = 0.0
308  table[i][component] = i
309 
310  # Set some metadata in the table (force OBSTYPE to exist)
311  metadata = copy.copy(self.getMetadata())
312  metadata["OBSTYPE"] = self._OBSTYPE
313  metadata[SCHEMA_NAME_KEY] = "FITS Region"
314  metadata[SCHEMA_VERSION_KEY] = 1
315  table.setMetadata(metadata)
316 
317  return table
318 
319  def writeFits(self, *args):
320  """Write defect list to FITS.
321 
322  Parameters
323  ----------
324  *args
325  Arguments to be forwarded to
326  `lsst.afw.table.BaseCatalog.writeFits`.
327  """
328  table = self.toFitsRegionTable()
329 
330  # Add some additional headers useful for tracking purposes
331  metadata = table.getMetadata()
332  now = datetime.datetime.utcnow()
333  metadata["DATE"] = now.isoformat()
334  metadata["CALIB_CREATION_DATE"] = now.strftime("%Y-%m-%d")
335  metadata["CALIB_CREATION_TIME"] = now.strftime("%T %Z").strip()
336 
337  table.writeFits(*args)
338 
339  def toSimpleTable(self):
340  """Convert defects to a simple table form that we use to write
341  to text files.
342 
343  Returns
344  -------
345  table : `lsst.afw.table.BaseCatalog`
346  Defects in simple tabular form.
347 
348  Notes
349  -----
350  These defect tables are used as the human readable definitions
351  of defects in calibration data definition repositories. The format
352  is to use four columns defined as follows:
353 
354  x0 : `int`
355  X coordinate of bottom left corner of box.
356  y0 : `int`
357  Y coordinate of bottom left corner of box.
358  width : `int`
359  X extent of the box.
360  height : `int`
361  Y extent of the box.
362  """
363  schema = lsst.afw.table.Schema()
364  x = schema.addField("x0", type="I", units="pix",
365  doc="X coordinate of bottom left corner of box")
366  y = schema.addField("y0", type="I", units="pix",
367  doc="Y coordinate of bottom left corner of box")
368  width = schema.addField("width", type="I", units="pix",
369  doc="X extent of box")
370  height = schema.addField("height", type="I", units="pix",
371  doc="Y extent of box")
372  table = lsst.afw.table.BaseCatalog(schema)
373  table.resize(len(self._defects))
374 
375  for i, defect in enumerate(self._defects):
376  box = defect.getBBox()
377  table[i][x] = box.getBeginX()
378  table[i][y] = box.getBeginY()
379  table[i][width] = box.getWidth()
380  table[i][height] = box.getHeight()
381 
382  # Set some metadata in the table (force OBSTYPE to exist)
383  metadata = copy.copy(self.getMetadata())
384  metadata["OBSTYPE"] = self._OBSTYPE
385  metadata[SCHEMA_NAME_KEY] = "Simple"
386  metadata[SCHEMA_VERSION_KEY] = 1
387  table.setMetadata(metadata)
388 
389  return table
390 
391  def writeText(self, filename):
392  """Write the defects out to a text file with the specified name.
393 
394  Parameters
395  ----------
396  filename : `str`
397  Name of the file to write. The file extension ".ecsv" will
398  always be used.
399 
400  Returns
401  -------
402  used : `str`
403  The name of the file used to write the data (which may be
404  different from the supplied name given the change to file
405  extension).
406 
407  Notes
408  -----
409  The file is written to ECSV format and will include any metadata
410  associated with the `Defects`.
411  """
412 
413  # Using astropy table is the easiest way to serialize to ecsv
414  afwTable = self.toSimpleTable()
415  table = afwTable.asAstropy()
416 
417  metadata = afwTable.getMetadata()
418  now = datetime.datetime.utcnow()
419  metadata["DATE"] = now.isoformat()
420  metadata["CALIB_CREATION_DATE"] = now.strftime("%Y-%m-%d")
421  metadata["CALIB_CREATION_TIME"] = now.strftime("%T %Z").strip()
422 
423  table.meta = metadata.toDict()
424 
425  # Force file extension to .ecsv
426  path, ext = os.path.splitext(filename)
427  filename = path + ".ecsv"
428  table.write(filename, format="ascii.ecsv")
429  return filename
430 
431  @staticmethod
432  def _get_values(values, n=1):
433  """Retrieve N values from the supplied values.
434 
435  Parameters
436  ----------
437  values : `numbers.Number` or `list` or `np.array`
438  Input values.
439  n : `int`
440  Number of values to retrieve.
441 
442  Returns
443  -------
444  vals : `list` or `np.array` or `numbers.Number`
445  Single value from supplied list if ``n`` is 1, or `list`
446  containing first ``n`` values from supplied values.
447 
448  Notes
449  -----
450  Some supplied tables have vectors in some columns that can also
451  be scalars. This method can be used to get the first number as
452  a scalar or the first N items from a vector as a vector.
453  """
454  if n == 1:
455  if isinstance(values, numbers.Number):
456  return values
457  else:
458  return values[0]
459 
460  return values[:n]
461 
462  @classmethod
463  def fromTable(cls, table):
464  """Construct a `Defects` from the contents of a
465  `~lsst.afw.table.BaseCatalog`.
466 
467  Parameters
468  ----------
469  table : `lsst.afw.table.BaseCatalog`
470  Table with one row per defect.
471 
472  Returns
473  -------
474  defects : `Defects`
475  A `Defects` list.
476 
477  Notes
478  -----
479  Two table formats are recognized. The first is the
480  `FITS regions <https://fits.gsfc.nasa.gov/registry/region.html>`_
481  definition tabular format written by `toFitsRegionTable` where the
482  pixel origin is corrected from FITS 1-based to a 0-based origin.
483  The second is the legacy defects format using columns ``x0``, ``y0``
484  (bottom left hand pixel of box in 0-based coordinates), ``width``
485  and ``height``.
486 
487  The FITS standard regions can only read BOX, POINT, or ROTBOX with
488  a zero degree rotation.
489  """
490 
491  defectList = []
492 
493  schema = table.getSchema()
494 
495  # Check schema to see which definitions we have
496  if "X" in schema and "Y" in schema and "R" in schema and "SHAPE" in schema:
497  # This is a FITS region style table
498  isFitsRegion = True
499 
500  elif "x0" in schema and "y0" in schema and "width" in schema and "height" in schema:
501  # This is a classic LSST-style defect table
502  isFitsRegion = False
503 
504  else:
505  raise ValueError("Unsupported schema for defects extraction")
506 
507  for r in table:
508  record = r.extract("*")
509 
510  if isFitsRegion:
511  # Coordinates can be arrays (some shapes in the standard
512  # require this)
513  # Correct for FITS 1-based origin
514  xcen = cls._get_values(record["X"]) - 1.0
515  ycen = cls._get_values(record["Y"]) - 1.0
516  shape = record["SHAPE"].upper()
517  if shape == "BOX":
519  lsst.geom.Extent2I(cls._get_values(record["R"],
520  n=2)))
521  elif shape == "POINT":
522  # Handle the case where we have an externally created
523  # FITS file.
524  box = lsst.geom.Point2I(xcen, ycen)
525  elif shape == "ROTBOX":
526  # Astropy regions always writes ROTBOX
527  rotang = cls._get_values(record["ROTANG"])
528  # We can support 0 or 90 deg
529  if math.isclose(rotang % 90.0, 0.0):
530  # Two values required
531  r = cls._get_values(record["R"], n=2)
532  if math.isclose(rotang % 180.0, 0.0):
533  width = r[0]
534  height = r[1]
535  else:
536  width = r[1]
537  height = r[0]
539  lsst.geom.Extent2I(width, height))
540  else:
541  log.warning("Defect can not be defined using ROTBOX with non-aligned rotation angle")
542  continue
543  else:
544  log.warning("Defect lists can only be defined using BOX or POINT not %s", shape)
545  continue
546 
547  elif "x0" in record and "y0" in record and "width" in record and "height" in record:
548  # This is a classic LSST-style defect table
549  box = lsst.geom.Box2I(lsst.geom.Point2I(record["x0"], record["y0"]),
550  lsst.geom.Extent2I(record["width"], record["height"]))
551 
552  defectList.append(box)
553 
554  defects = cls(defectList)
555  defects.setMetadata(table.getMetadata())
556 
557  # Once read, the schema headers are irrelevant
558  metadata = defects.getMetadata()
559  for k in (SCHEMA_NAME_KEY, SCHEMA_VERSION_KEY):
560  if k in metadata:
561  del metadata[k]
562 
563  return defects
564 
565  @classmethod
566  def readFits(cls, *args):
567  """Read defect list from FITS table.
568 
569  Parameters
570  ----------
571  *args
572  Arguments to be forwarded to
573  `lsst.afw.table.BaseCatalog.writeFits`.
574 
575  Returns
576  -------
577  defects : `Defects`
578  Defects read from a FITS table.
579  """
580  table = lsst.afw.table.BaseCatalog.readFits(*args)
581  return cls.fromTable(table)
582 
583  @classmethod
584  def readText(cls, filename):
585  """Read defect list from standard format text table file.
586 
587  Parameters
588  ----------
589  filename : `str`
590  Name of the file containing the defects definitions.
591 
592  Returns
593  -------
594  defects : `Defects`
595  Defects read from a FITS table.
596  """
597  table = astropy.table.Table.read(filename)
598 
599  # Need to convert the Astropy table to afw table
600  schema = lsst.afw.table.Schema()
601  for colName in table.columns:
602  schema.addField(colName, units=str(table[colName].unit),
603  type=table[colName].dtype.type)
604 
605  # Create AFW table that is required by fromTable()
606  afwTable = lsst.afw.table.BaseCatalog(schema)
607 
608  afwTable.resize(len(table))
609  for colName in table.columns:
610  # String columns will fail -- currently we do not expect any
611  afwTable[colName] = table[colName]
612 
613  # Copy in the metadata from the astropy table
614  metadata = PropertyList()
615  for k, v in table.meta.items():
616  metadata[k] = v
617  afwTable.setMetadata(metadata)
618 
619  # Extract defect information from the table itself
620  return cls.fromTable(afwTable)
621 
622  @classmethod
623  def readLsstDefectsFile(cls, filename):
624  """Read defects information from a legacy LSST format text file.
625 
626  Parameters
627  ----------
628  filename : `str`
629  Name of text file containing the defect information.
630 
631  Returns
632  -------
633  defects : `Defects`
634  The defects.
635 
636  Notes
637  -----
638  These defect text files are used as the human readable definitions
639  of defects in calibration data definition repositories. The format
640  is to use four columns defined as follows:
641 
642  x0 : `int`
643  X coordinate of bottom left corner of box.
644  y0 : `int`
645  Y coordinate of bottom left corner of box.
646  width : `int`
647  X extent of the box.
648  height : `int`
649  Y extent of the box.
650 
651  Files of this format were used historically to represent defects
652  in simple text form. Use `Defects.readText` and `Defects.writeText`
653  to use the more modern format.
654  """
655  # Use loadtxt so that ValueError is thrown if the file contains a
656  # non-integer value. genfromtxt converts bad values to -1.
657  defect_array = np.loadtxt(filename,
658  dtype=[("x0", "int"), ("y0", "int"),
659  ("x_extent", "int"), ("y_extent", "int")])
660 
661  return cls(lsst.geom.Box2I(lsst.geom.Point2I(row["x0"], row["y0"]),
662  lsst.geom.Extent2I(row["x_extent"], row["y_extent"]))
663  for row in defect_array)
664 
665  @classmethod
666  def fromFootprintList(cls, fpList):
667  """Compute a defect list from a footprint list, optionally growing
668  the footprints.
669 
670  Parameters
671  ----------
672  fpList : `list` of `lsst.afw.detection.Footprint`
673  Footprint list to process.
674 
675  Returns
676  -------
677  defects : `Defects`
678  List of defects.
679  """
680  return cls(itertools.chain.from_iterable(lsst.afw.detection.footprintToBBoxList(fp)
681  for fp in fpList))
682 
683  @classmethod
684  def fromMask(cls, maskedImage, maskName):
685  """Compute a defect list from a specified mask plane.
686 
687  Parameters
688  ----------
689  maskedImage : `lsst.afw.image.MaskedImage`
690  Image to process.
691  maskName : `str` or `list`
692  Mask plane name, or list of names to convert.
693 
694  Returns
695  -------
696  defects : `Defects`
697  Defect list constructed from masked pixels.
698  """
699  mask = maskedImage.getMask()
700  thresh = lsst.afw.detection.Threshold(mask.getPlaneBitMask(maskName),
701  lsst.afw.detection.Threshold.BITMASK)
702  fpList = lsst.afw.detection.FootprintSet(mask, thresh).getFootprints()
703  return cls.fromFootprintList(fpList)
def maskPixels(self, maskedImage, maskName="BAD")
Definition: defects.py:248
def setMetadata(self, metadata=None)
Definition: defects.py:196
def __setitem__(self, index, value)
Definition: defects.py:151
def __init__(self, defectList=None, metadata=None)
Definition: defects.py:98
def policyToBadRegionList(policyFile)
Definition: defects.py:57
def readLsstDefectsFile(cls, filename)
Definition: defects.py:623
std::vector< lsst::geom::Box2I > footprintToBBoxList(Footprint const &footprint)
def insert(self, index, value)
Definition: defects.py:181
static Box2I makeCenteredBox(Point2D const &center, Extent const &size)
def fromMask(cls, maskedImage, maskName)
Definition: defects.py:684