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