LSSTApplications  19.0.0-14-gb0260a2+463b159a14,20.0.0+0f5db2c0a9,20.0.0+aad4e2d4eb,20.0.0+bddc4f4cbe,20.0.0+f04b2999c2,20.0.0+f47236d43c,20.0.0-1-g253301a+f47236d43c,20.0.0-1-g2b7511a+821cb29515,20.0.0-1-g5b95a8c+b1532f98d2,20.0.0-1-gedffbd8+e71ed47741,20.0.0-11-gd9dafd18+74fdb5a06a,20.0.0-16-gfab17e72e+e13f79ae9f,20.0.0-19-gb3c8180+37ae35d132,20.0.0-2-g4dae9ad+9ccd17b1a3,20.0.0-2-g61b8584+5d3db074ba,20.0.0-2-gb780d76+d529cf1a41,20.0.0-2-ged6426c+928601e9b9,20.0.0-2-gf072044+f47236d43c,20.0.0-2-gf1f7952+9ccd17b1a3,20.0.0-22-gdf434b7+9ccd17b1a3,20.0.0-23-g10eeb28+96244339af,20.0.0-25-g3dcad98+1dc6abc28c,20.0.0-3-g1653f94+5d3db074ba,20.0.0-3-g4cc78c6+e0676b0dc8,20.0.0-3-g8f21e14+1c0155d006,20.0.0-3-gbd60e8c+187b78b4b8,20.0.0-3-gbecbe05+5dad399003,20.0.0-36-g96d5af2b+2d45f6b8a5,20.0.0-4-g97dc21a+1dc6abc28c,20.0.0-4-gb4befbc+925a411f2f,20.0.0-5-gdfe0fee+60390a8ce6,20.0.0-5-gfbfe500+d529cf1a41,20.0.0-6-g64f541c+d529cf1a41,20.0.0-6-g9a5b7a1+a87ba7c33a,20.0.0-67-g32d6278+9ea4af2944,20.0.0-9-g4aef684+0219fac3b9,w.2020.45
LSSTDataManagementBasePackage
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 
26 __all__ = ("Defects",)
27 
28 import logging
29 import itertools
30 import collections.abc
31 import contextlib
32 import numpy as np
33 import copy
34 import datetime
35 import math
36 import numbers
37 import os.path
38 import warnings
39 import astropy.table
40 
41 import lsst.geom
42 import lsst.afw.table
43 import lsst.afw.detection
44 import lsst.afw.image
45 import lsst.afw.geom
46 from lsst.daf.base import PropertyList
47 
48 from . import Defect
49 
50 log = logging.getLogger(__name__)
51 
52 SCHEMA_NAME_KEY = "DEFECTS_SCHEMA"
53 SCHEMA_VERSION_KEY = "DEFECTS_SCHEMA_VERSION"
54 
55 
56 class Defects(collections.abc.MutableSequence):
57  """Collection of `lsst.meas.algorithms.Defect`.
58 
59  Parameters
60  ----------
61  defectList : iterable of `lsst.meas.algorithms.Defect`
62  or `lsst.geom.BoxI`, optional
63  Collections of defects to apply to the image.
64  metadata : `lsst.daf.base.PropertyList`, optional
65  Metadata to associate with the defects. Will be copied and
66  overwrite existing metadata, if any. If not supplied the existing
67  metadata will be reset.
68  normalize_on_init : `bool`
69  If True, normalization is applied to the defects in ``defectList`` to
70  remove duplicates, eliminate overlaps, etc.
71 
72  Notes
73  -----
74  Defects are stored within this collection in a "reduced" or "normalized"
75  form: rather than simply storing the bounding boxes which are added to the
76  collection, we eliminate overlaps and duplicates. This normalization
77  procedure may introduce overhead when adding many new defects; it may be
78  temporarily disabled using the `Defects.bulk_update` context manager if
79  necessary.
80  """
81 
82  _OBSTYPE = "defects"
83  """The calibration type used for ingest.
84  """
85 
86  def __init__(self, defectList=None, metadata=None, *, normalize_on_init=True):
87  self._defects = []
88 
89  if defectList is not None:
90  self._bulk_update = True
91  for d in defectList:
92  self.append(d)
93  self._bulk_update = False
94 
95  if normalize_on_init:
96  self._normalize()
97 
98  if metadata is not None:
99  self._metadata = metadata
100  else:
101  self.setMetadata()
102 
103  def _check_value(self, value):
104  """Check that the supplied value is a `~lsst.meas.algorithms.Defect`
105  or can be converted to one.
106 
107  Parameters
108  ----------
109  value : `object`
110  Value to check.
111 
112  Returns
113  -------
114  new : `~lsst.meas.algorithms.Defect`
115  Either the supplied value or a new object derived from it.
116 
117  Raises
118  ------
119  ValueError
120  Raised if the supplied value can not be converted to
121  `~lsst.meas.algorithms.Defect`
122  """
123  if isinstance(value, Defect):
124  pass
125  elif isinstance(value, lsst.geom.BoxI):
126  value = Defect(value)
127  elif isinstance(value, lsst.geom.PointI):
128  value = Defect(lsst.geom.Box2I(value, lsst.geom.Extent2I(1, 1)))
129  elif isinstance(value, lsst.afw.image.DefectBase):
130  value = Defect(value.getBBox())
131  else:
132  raise ValueError(f"Defects must be of type Defect, BoxI, or PointI, not '{value!r}'")
133  return value
134 
135  def __len__(self):
136  return len(self._defects)
137 
138  def __getitem__(self, index):
139  return self._defects[index]
140 
141  def __setitem__(self, index, value):
142  """Can be given a `~lsst.meas.algorithms.Defect` or a `lsst.geom.BoxI`
143  """
144  self._defects[index] = self._check_value(value)
145  self._normalize()
146 
147  def __iter__(self):
148  return iter(self._defects)
149 
150  def __delitem__(self, index):
151  del self._defects[index]
152 
153  def __eq__(self, other):
154  """Compare if two `Defects` are equal.
155 
156  Two `Defects` are equal if their bounding boxes are equal and in
157  the same order. Metadata content is ignored.
158  """
159  if not isinstance(other, self.__class__):
160  return False
161 
162  # checking the bboxes with zip() only works if same length
163  if len(self) != len(other):
164  return False
165 
166  # Assume equal if bounding boxes are equal
167  for d1, d2 in zip(self, other):
168  if d1.getBBox() != d2.getBBox():
169  return False
170 
171  return True
172 
173  def __str__(self):
174  return "Defects(" + ",".join(str(d.getBBox()) for d in self) + ")"
175 
176  def _normalize(self):
177  """Recalculate defect bounding boxes for efficiency.
178 
179  Notes
180  -----
181  Ideally, this would generate the provably-minimal set of bounding
182  boxes necessary to represent the defects. At present, however, that
183  doesn't happen: see DM-24781. In the cases of substantial overlaps or
184  duplication, though, this will produce a much reduced set.
185  """
186  # In bulk-update mode, normalization is a no-op.
187  if self._bulk_update:
188  return
189 
190  # work out the minimum and maximum bounds from all defect regions.
191  minX, minY, maxX, maxY = float('inf'), float('inf'), float('-inf'), float('-inf')
192  for defect in self:
193  bbox = defect.getBBox()
194  minX = min(minX, bbox.getMinX())
195  minY = min(minY, bbox.getMinY())
196  maxX = max(maxX, bbox.getMaxX())
197  maxY = max(maxY, bbox.getMaxY())
198 
199  region = lsst.geom.Box2I(lsst.geom.Point2I(minX, minY),
200  lsst.geom.Point2I(maxX, maxY))
201 
202  mi = lsst.afw.image.MaskedImageF(region)
203  self.maskPixels(mi, maskName="BAD")
204  self._defects = Defects.fromMask(mi, "BAD")._defects
205 
206  @contextlib.contextmanager
207  def bulk_update(self):
208  """Temporarily suspend normalization of the defect list.
209  """
210  self._bulk_update = True
211  try:
212  yield
213  finally:
214  self._bulk_update = False
215  self._normalize()
216 
217  def insert(self, index, value):
218  self._defects.insert(index, self._check_value(value))
219  self._normalize()
220 
221  def getMetadata(self):
222  """Retrieve metadata associated with these `Defects`.
223 
224  Returns
225  -------
226  meta : `lsst.daf.base.PropertyList`
227  Metadata. The returned `~lsst.daf.base.PropertyList` can be
228  modified by the caller and the changes will be written to
229  external files.
230  """
231  return self._metadata
232 
233  def setMetadata(self, metadata=None):
234  """Store a copy of the supplied metadata with the defects.
235 
236  Parameters
237  ----------
238  metadata : `lsst.daf.base.PropertyList`, optional
239  Metadata to associate with the defects. Will be copied and
240  overwrite existing metadata. If not supplied the existing
241  metadata will be reset.
242  """
243  if metadata is None:
244  self._metadata = PropertyList()
245  else:
246  self._metadata = copy.copy(metadata)
247 
248  # Ensure that we have the obs type required by calibration ingest
249  self._metadata["OBSTYPE"] = self._OBSTYPE
250 
251  def copy(self):
252  """Copy the defects to a new list, creating new defects from the
253  bounding boxes.
254 
255  Returns
256  -------
257  new : `Defects`
258  New list with new `Defect` entries.
259 
260  Notes
261  -----
262  This is not a shallow copy in that new `Defect` instances are
263  created from the original bounding boxes. It's also not a deep
264  copy since the bounding boxes are not recreated.
265  """
266  return self.__class__(d.getBBox() for d in self)
267 
268  def transpose(self):
269  """Make a transposed copy of this defect list.
270 
271  Returns
272  -------
273  retDefectList : `Defects`
274  Transposed list of defects.
275  """
276  retDefectList = self.__class__()
277  for defect in self:
278  bbox = defect.getBBox()
279  dimensions = bbox.getDimensions()
280  nbbox = lsst.geom.Box2I(lsst.geom.Point2I(bbox.getMinY(), bbox.getMinX()),
281  lsst.geom.Extent2I(dimensions[1], dimensions[0]))
282  retDefectList.append(nbbox)
283  return retDefectList
284 
285  def maskPixels(self, maskedImage, maskName="BAD"):
286  """Set mask plane based on these defects.
287 
288  Parameters
289  ----------
290  maskedImage : `lsst.afw.image.MaskedImage`
291  Image to process. Only the mask plane is updated.
292  maskName : str, optional
293  Mask plane name to use.
294  """
295  # mask bad pixels
296  mask = maskedImage.getMask()
297  bitmask = mask.getPlaneBitMask(maskName)
298  for defect in self:
299  bbox = defect.getBBox()
300  lsst.afw.geom.SpanSet(bbox).clippedTo(mask.getBBox()).setMask(mask, bitmask)
301 
302  def toFitsRegionTable(self):
303  """Convert defect list to `~lsst.afw.table.BaseCatalog` using the
304  FITS region standard.
305 
306  Returns
307  -------
308  table : `lsst.afw.table.BaseCatalog`
309  Defects in tabular form.
310 
311  Notes
312  -----
313  The table created uses the
314  `FITS regions <https://fits.gsfc.nasa.gov/registry/region.html>`_
315  definition tabular format. The ``X`` and ``Y`` coordinates are
316  converted to FITS Physical coordinates that have origin pixel (1, 1)
317  rather than the (0, 0) used in LSST software.
318  """
319 
320  nrows = len(self._defects)
321 
322  schema = lsst.afw.table.Schema()
323  x = schema.addField("X", type="D", units="pix", doc="X coordinate of center of shape")
324  y = schema.addField("Y", type="D", units="pix", doc="Y coordinate of center of shape")
325  shape = schema.addField("SHAPE", type="String", size=16, doc="Shape defined by these values")
326  r = schema.addField("R", type="ArrayD", size=2, units="pix", doc="Extents")
327  rotang = schema.addField("ROTANG", type="D", units="deg", doc="Rotation angle")
328  component = schema.addField("COMPONENT", type="I", doc="Index of this region")
329  table = lsst.afw.table.BaseCatalog(schema)
330  table.resize(nrows)
331 
332  if nrows:
333  # Adding entire columns is more efficient than adding
334  # each element separately
335  xCol = []
336  yCol = []
337  rCol = []
338 
339  for i, defect in enumerate(self._defects):
340  box = defect.getBBox()
341  center = box.getCenter()
342  # Correct for the FITS 1-based offset
343  xCol.append(center.getX() + 1.0)
344  yCol.append(center.getY() + 1.0)
345 
346  width = box.width
347  height = box.height
348 
349  if width == 1 and height == 1:
350  # Call this a point
351  shapeType = "POINT"
352  else:
353  shapeType = "BOX"
354 
355  # Strings have to be added per row
356  table[i][shape] = shapeType
357 
358  rCol.append(np.array([width, height], dtype=np.float64))
359 
360  # Assign the columns
361  table[x] = np.array(xCol, dtype=np.float64)
362  table[y] = np.array(yCol, dtype=np.float64)
363 
364  table[r] = np.array(rCol)
365  table[rotang] = np.zeros(nrows, dtype=np.float64)
366  table[component] = np.arange(nrows)
367 
368  # Set some metadata in the table (force OBSTYPE to exist)
369  metadata = copy.copy(self.getMetadata())
370  metadata["OBSTYPE"] = self._OBSTYPE
371  metadata[SCHEMA_NAME_KEY] = "FITS Region"
372  metadata[SCHEMA_VERSION_KEY] = 1
373  table.setMetadata(metadata)
374 
375  return table
376 
377  def writeFits(self, *args):
378  """Write defect list to FITS.
379 
380  Parameters
381  ----------
382  *args
383  Arguments to be forwarded to
384  `lsst.afw.table.BaseCatalog.writeFits`.
385  """
386  table = self.toFitsRegionTable()
387 
388  # Add some additional headers useful for tracking purposes
389  metadata = table.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.writeFits(*args)
396 
397  def toSimpleTable(self):
398  """Convert defects to a simple table form that we use to write
399  to text files.
400 
401  Returns
402  -------
403  table : `lsst.afw.table.BaseCatalog`
404  Defects in simple tabular form.
405 
406  Notes
407  -----
408  These defect tables are used as the human readable definitions
409  of defects in calibration data definition repositories. The format
410  is to use four columns defined as follows:
411 
412  x0 : `int`
413  X coordinate of bottom left corner of box.
414  y0 : `int`
415  Y coordinate of bottom left corner of box.
416  width : `int`
417  X extent of the box.
418  height : `int`
419  Y extent of the box.
420  """
421  schema = lsst.afw.table.Schema()
422  x = schema.addField("x0", type="I", units="pix",
423  doc="X coordinate of bottom left corner of box")
424  y = schema.addField("y0", type="I", units="pix",
425  doc="Y coordinate of bottom left corner of box")
426  width = schema.addField("width", type="I", units="pix",
427  doc="X extent of box")
428  height = schema.addField("height", type="I", units="pix",
429  doc="Y extent of box")
430  table = lsst.afw.table.BaseCatalog(schema)
431 
432  nrows = len(self._defects)
433  table.resize(nrows)
434 
435  if nrows:
436 
437  xCol = []
438  yCol = []
439  widthCol = []
440  heightCol = []
441 
442  for defect in self._defects:
443  box = defect.getBBox()
444  xCol.append(box.getBeginX())
445  yCol.append(box.getBeginY())
446  widthCol.append(box.getWidth())
447  heightCol.append(box.getHeight())
448 
449  table[x] = np.array(xCol, dtype=np.int64)
450  table[y] = np.array(yCol, dtype=np.int64)
451  table[width] = np.array(widthCol, dtype=np.int64)
452  table[height] = np.array(heightCol, dtype=np.int64)
453 
454  # Set some metadata in the table (force OBSTYPE to exist)
455  metadata = copy.copy(self.getMetadata())
456  metadata["OBSTYPE"] = self._OBSTYPE
457  metadata[SCHEMA_NAME_KEY] = "Simple"
458  metadata[SCHEMA_VERSION_KEY] = 1
459  table.setMetadata(metadata)
460 
461  return table
462 
463  def writeText(self, filename):
464  """Write the defects out to a text file with the specified name.
465 
466  Parameters
467  ----------
468  filename : `str`
469  Name of the file to write. The file extension ".ecsv" will
470  always be used.
471 
472  Returns
473  -------
474  used : `str`
475  The name of the file used to write the data (which may be
476  different from the supplied name given the change to file
477  extension).
478 
479  Notes
480  -----
481  The file is written to ECSV format and will include any metadata
482  associated with the `Defects`.
483  """
484 
485  # Using astropy table is the easiest way to serialize to ecsv
486  afwTable = self.toSimpleTable()
487  table = afwTable.asAstropy()
488 
489  metadata = afwTable.getMetadata()
490  now = datetime.datetime.utcnow()
491  metadata["DATE"] = now.isoformat()
492  metadata["CALIB_CREATION_DATE"] = now.strftime("%Y-%m-%d")
493  metadata["CALIB_CREATION_TIME"] = now.strftime("%T %Z").strip()
494 
495  table.meta = metadata.toDict()
496 
497  # Force file extension to .ecsv
498  path, ext = os.path.splitext(filename)
499  filename = path + ".ecsv"
500  table.write(filename, format="ascii.ecsv")
501  return filename
502 
503  @staticmethod
504  def _get_values(values, n=1):
505  """Retrieve N values from the supplied values.
506 
507  Parameters
508  ----------
509  values : `numbers.Number` or `list` or `np.array`
510  Input values.
511  n : `int`
512  Number of values to retrieve.
513 
514  Returns
515  -------
516  vals : `list` or `np.array` or `numbers.Number`
517  Single value from supplied list if ``n`` is 1, or `list`
518  containing first ``n`` values from supplied values.
519 
520  Notes
521  -----
522  Some supplied tables have vectors in some columns that can also
523  be scalars. This method can be used to get the first number as
524  a scalar or the first N items from a vector as a vector.
525  """
526  if n == 1:
527  if isinstance(values, numbers.Number):
528  return values
529  else:
530  return values[0]
531 
532  return values[:n]
533 
534  @classmethod
535  def fromTable(cls, table, normalize_on_init=True):
536  """Construct a `Defects` from the contents of a
537  `~lsst.afw.table.BaseCatalog`.
538 
539  Parameters
540  ----------
541  table : `lsst.afw.table.BaseCatalog`
542  Table with one row per defect.
543  normalize_on_init : `bool`, optional
544  If `True`, normalization is applied to the defects listed in the
545  table to remove duplicates, eliminate overlaps, etc. Otherwise
546  the defects in the returned object exactly match those in the
547  table.
548 
549  Returns
550  -------
551  defects : `Defects`
552  A `Defects` list.
553 
554  Notes
555  -----
556  Two table formats are recognized. The first is the
557  `FITS regions <https://fits.gsfc.nasa.gov/registry/region.html>`_
558  definition tabular format written by `toFitsRegionTable` where the
559  pixel origin is corrected from FITS 1-based to a 0-based origin.
560  The second is the legacy defects format using columns ``x0``, ``y0``
561  (bottom left hand pixel of box in 0-based coordinates), ``width``
562  and ``height``.
563 
564  The FITS standard regions can only read BOX, POINT, or ROTBOX with
565  a zero degree rotation.
566  """
567 
568  defectList = []
569 
570  schema = table.getSchema()
571 
572  # Check schema to see which definitions we have
573  if "X" in schema and "Y" in schema and "R" in schema and "SHAPE" in schema:
574  # This is a FITS region style table
575  isFitsRegion = True
576 
577  # Preselect the keys
578  xKey = schema["X"].asKey()
579  yKey = schema["Y"].asKey()
580  shapeKey = schema["SHAPE"].asKey()
581  rKey = schema["R"].asKey()
582  rotangKey = schema["ROTANG"].asKey()
583 
584  elif "x0" in schema and "y0" in schema and "width" in schema and "height" in schema:
585  # This is a classic LSST-style defect table
586  isFitsRegion = False
587 
588  # Preselect the keys
589  xKey = schema["x0"].asKey()
590  yKey = schema["y0"].asKey()
591  widthKey = schema["width"].asKey()
592  heightKey = schema["height"].asKey()
593 
594  else:
595  raise ValueError("Unsupported schema for defects extraction")
596 
597  for record in table:
598 
599  if isFitsRegion:
600  # Coordinates can be arrays (some shapes in the standard
601  # require this)
602  # Correct for FITS 1-based origin
603  xcen = cls._get_values(record[xKey]) - 1.0
604  ycen = cls._get_values(record[yKey]) - 1.0
605  shape = record[shapeKey].upper()
606  if shape == "BOX":
608  lsst.geom.Extent2I(cls._get_values(record[rKey],
609  n=2)))
610  elif shape == "POINT":
611  # Handle the case where we have an externally created
612  # FITS file.
613  box = lsst.geom.Point2I(xcen, ycen)
614  elif shape == "ROTBOX":
615  # Astropy regions always writes ROTBOX
616  rotang = cls._get_values(record[rotangKey])
617  # We can support 0 or 90 deg
618  if math.isclose(rotang % 90.0, 0.0):
619  # Two values required
620  r = cls._get_values(record[rKey], n=2)
621  if math.isclose(rotang % 180.0, 0.0):
622  width = r[0]
623  height = r[1]
624  else:
625  width = r[1]
626  height = r[0]
628  lsst.geom.Extent2I(width, height))
629  else:
630  log.warning("Defect can not be defined using ROTBOX with non-aligned rotation angle")
631  continue
632  else:
633  log.warning("Defect lists can only be defined using BOX or POINT not %s", shape)
634  continue
635 
636  else:
637  # This is a classic LSST-style defect table
638  box = lsst.geom.Box2I(lsst.geom.Point2I(record[xKey], record[yKey]),
639  lsst.geom.Extent2I(record[widthKey], record[heightKey]))
640 
641  defectList.append(box)
642 
643  defects = cls(defectList, normalize_on_init=normalize_on_init)
644  defects.setMetadata(table.getMetadata())
645 
646  # Once read, the schema headers are irrelevant
647  metadata = defects.getMetadata()
648  for k in (SCHEMA_NAME_KEY, SCHEMA_VERSION_KEY):
649  if k in metadata:
650  del metadata[k]
651 
652  return defects
653 
654  @classmethod
655  def readFits(cls, *args, normalize_on_init=False):
656  """Read defect list from FITS table.
657 
658  Parameters
659  ----------
660  *args
661  Arguments to be forwarded to
662  `lsst.afw.table.BaseCatalog.readFits`.
663  normalize_on_init : `bool`, optional
664  If `True`, normalization is applied to the defects read fom the
665  file to remove duplicates, eliminate overlaps, etc. Otherwise
666  the defects in the returned object exactly match those in the
667  file.
668 
669  Returns
670  -------
671  defects : `Defects`
672  Defects read from a FITS table.
673  """
674  table = lsst.afw.table.BaseCatalog.readFits(*args)
675  return cls.fromTable(table, normalize_on_init=normalize_on_init)
676 
677  @classmethod
678  def readText(cls, filename, normalize_on_init=False):
679  """Read defect list from standard format text table file.
680 
681  Parameters
682  ----------
683  filename : `str`
684  Name of the file containing the defects definitions.
685  normalize_on_init : `bool`, optional
686  If `True`, normalization is applied to the defects read fom the
687  file to remove duplicates, eliminate overlaps, etc. Otherwise
688  the defects in the returned object exactly match those in the
689  file.
690 
691  Returns
692  -------
693  defects : `Defects`
694  Defects read from a FITS table.
695  """
696  with warnings.catch_warnings():
697  # Squash warnings due to astropy failure to close files; we think
698  # this is a real problem, but the warnings are even worse.
699  # https://github.com/astropy/astropy/issues/8673
700  warnings.filterwarnings("ignore", category=ResourceWarning, module="astropy.io.ascii")
701  table = astropy.table.Table.read(filename)
702 
703  # Need to convert the Astropy table to afw table
704  schema = lsst.afw.table.Schema()
705  for colName in table.columns:
706  schema.addField(colName, units=str(table[colName].unit),
707  type=table[colName].dtype.type)
708 
709  # Create AFW table that is required by fromTable()
710  afwTable = lsst.afw.table.BaseCatalog(schema)
711 
712  afwTable.resize(len(table))
713  for colName in table.columns:
714  # String columns will fail -- currently we do not expect any
715  afwTable[colName] = table[colName]
716 
717  # Copy in the metadata from the astropy table
718  metadata = PropertyList()
719  for k, v in table.meta.items():
720  metadata[k] = v
721  afwTable.setMetadata(metadata)
722 
723  # Extract defect information from the table itself
724  return cls.fromTable(afwTable, normalize_on_init=normalize_on_init)
725 
726  @classmethod
727  def readLsstDefectsFile(cls, filename, normalize_on_init=False):
728  """Read defects information from a legacy LSST format text file.
729 
730  Parameters
731  ----------
732  filename : `str`
733  Name of text file containing the defect information.
734  normalize_on_init : `bool`, optional
735  If `True`, normalization is applied to the defects read fom the
736  file to remove duplicates, eliminate overlaps, etc. Otherwise
737  the defects in the returned object exactly match those in the
738  file.
739 
740  Returns
741  -------
742  defects : `Defects`
743  The defects.
744 
745  Notes
746  -----
747  These defect text files are used as the human readable definitions
748  of defects in calibration data definition repositories. The format
749  is to use four columns defined as follows:
750 
751  x0 : `int`
752  X coordinate of bottom left corner of box.
753  y0 : `int`
754  Y coordinate of bottom left corner of box.
755  width : `int`
756  X extent of the box.
757  height : `int`
758  Y extent of the box.
759 
760  Files of this format were used historically to represent defects
761  in simple text form. Use `Defects.readText` and `Defects.writeText`
762  to use the more modern format.
763  """
764  # Use loadtxt so that ValueError is thrown if the file contains a
765  # non-integer value. genfromtxt converts bad values to -1.
766  defect_array = np.loadtxt(filename,
767  dtype=[("x0", "int"), ("y0", "int"),
768  ("x_extent", "int"), ("y_extent", "int")])
769 
770  defects = (lsst.geom.Box2I(lsst.geom.Point2I(row["x0"], row["y0"]),
771  lsst.geom.Extent2I(row["x_extent"], row["y_extent"]))
772  for row in defect_array)
773 
774  return cls(defects, normalize_on_init=normalize_on_init)
775 
776  @classmethod
777  def fromFootprintList(cls, fpList):
778  """Compute a defect list from a footprint list, optionally growing
779  the footprints.
780 
781  Parameters
782  ----------
783  fpList : `list` of `lsst.afw.detection.Footprint`
784  Footprint list to process.
785 
786  Returns
787  -------
788  defects : `Defects`
789  List of defects.
790  """
791  # normalize_on_init is set to False to avoid recursively calling
792  # fromMask/fromFootprintList in Defects.__init__.
793  return cls(itertools.chain.from_iterable(lsst.afw.detection.footprintToBBoxList(fp)
794  for fp in fpList), normalize_on_init=False)
795 
796  @classmethod
797  def fromMask(cls, maskedImage, maskName):
798  """Compute a defect list from a specified mask plane.
799 
800  Parameters
801  ----------
802  maskedImage : `lsst.afw.image.MaskedImage`
803  Image to process.
804  maskName : `str` or `list`
805  Mask plane name, or list of names to convert.
806 
807  Returns
808  -------
809  defects : `Defects`
810  Defect list constructed from masked pixels.
811  """
812  mask = maskedImage.getMask()
813  thresh = lsst.afw.detection.Threshold(mask.getPlaneBitMask(maskName),
814  lsst.afw.detection.Threshold.BITMASK)
815  fpList = lsst.afw.detection.FootprintSet(mask, thresh).getFootprints()
816  return cls.fromFootprintList(fpList)
lsst::meas::algorithms.defects.Defects.__setitem__
def __setitem__(self, index, value)
Definition: defects.py:141
lsst::afw::image
Backwards-compatibility support for depersisting the old Calib (FluxMag0/FluxMag0Err) objects.
Definition: imageAlgorithm.dox:1
lsst::meas::algorithms.defects.Defects._get_values
def _get_values(values, n=1)
Definition: defects.py:504
lsst::meas::algorithms.defects.Defects
Definition: defects.py:56
lsst::meas::algorithms.defects.Defects.__delitem__
def __delitem__(self, index)
Definition: defects.py:150
lsst::meas::algorithms.defects.Defects.readText
def readText(cls, filename, normalize_on_init=False)
Definition: defects.py:678
lsst::meas::algorithms.defects.Defects.__init__
def __init__(self, defectList=None, metadata=None, *normalize_on_init=True)
Definition: defects.py:86
lsst::meas::algorithms.defects.Defects.fromTable
def fromTable(cls, table, normalize_on_init=True)
Definition: defects.py:535
lsst::geom::Box2I::makeCenteredBox
static Box2I makeCenteredBox(Point2D const &center, Extent const &size)
Create a box centered as closely as possible on a particular point.
Definition: Box.cc:97
lsst::meas::algorithms.defects.Defects._check_value
def _check_value(self, value)
Definition: defects.py:103
lsst::meas::algorithms.defects.Defects.bulk_update
def bulk_update(self)
Definition: defects.py:207
lsst::daf::base::PropertyList
Class for storing ordered metadata with comments.
Definition: PropertyList.h:68
lsst::meas::algorithms.defects.Defects.getMetadata
def getMetadata(self)
Definition: defects.py:221
lsst::meas::algorithms.defects.Defects._normalize
def _normalize(self)
Definition: defects.py:176
lsst::afw::table::Schema
Defines the fields and offsets for a table.
Definition: Schema.h:50
strip
bool strip
Definition: fits.cc:911
lsst::meas::algorithms.defects.Defects.__str__
def __str__(self)
Definition: defects.py:173
lsst::afw::geom.transform.transformContinued.cls
cls
Definition: transformContinued.py:33
lsst::meas::algorithms.defects.Defects.fromFootprintList
def fromFootprintList(cls, fpList)
Definition: defects.py:777
lsst::meas::algorithms.defects.Defects._OBSTYPE
string _OBSTYPE
Definition: defects.py:82
lsst::meas::algorithms.defects.Defects.toSimpleTable
def toSimpleTable(self)
Definition: defects.py:397
lsst::afw::detection::FootprintSet
A set of Footprints, associated with a MaskedImage.
Definition: FootprintSet.h:53
lsst::afw::image::DefectBase
Encapsulate information about a bad portion of a detector.
Definition: Defect.h:41
lsst::meas::algorithms.defects.Defects.writeFits
def writeFits(self, *args)
Definition: defects.py:377
lsst::afw::geom::SpanSet
A compact representation of a collection of pixels.
Definition: SpanSet.h:78
lsst::afw::detection::footprintToBBoxList
std::vector< lsst::geom::Box2I > footprintToBBoxList(Footprint const &footprint)
Return a list of BBoxs, whose union contains exactly the pixels in the footprint, neither more nor le...
Definition: Footprint.cc:352
lsst::meas::algorithms.defects.Defects.__iter__
def __iter__(self)
Definition: defects.py:147
lsst::meas::algorithms.defects.Defects._defects
_defects
Definition: defects.py:87
lsst::afw::detection::Threshold
A Threshold is used to pass a threshold value to detection algorithms.
Definition: Threshold.h:43
lsst::meas::algorithms.defects.Defects.readFits
def readFits(cls, *args, normalize_on_init=False)
Definition: defects.py:655
lsst::meas::algorithms.defects.Defects.__eq__
def __eq__(self, other)
Definition: defects.py:153
lsst::meas::algorithms.defects.Defects.readLsstDefectsFile
def readLsstDefectsFile(cls, filename, normalize_on_init=False)
Definition: defects.py:727
max
int max
Definition: BoundedField.cc:104
lsst::afw::table
Definition: table.dox:3
lsst::meas::algorithms.defects.Defects.setMetadata
def setMetadata(self, metadata=None)
Definition: defects.py:233
lsst::meas::algorithms.defects.Defects.fromMask
def fromMask(cls, maskedImage, maskName)
Definition: defects.py:797
lsst::afw::detection
Definition: Footprint.h:50
lsst::meas::algorithms.defects.Defects.transpose
def transpose(self)
Definition: defects.py:268
lsst::geom
Definition: AffineTransform.h:36
lsst::meas::algorithms.defects.Defects.__getitem__
def __getitem__(self, index)
Definition: defects.py:138
lsst::daf::base
Definition: Utils.h:47
lsst::meas::algorithms.defects.Defects.maskPixels
def maskPixels(self, maskedImage, maskName="BAD")
Definition: defects.py:285
min
int min
Definition: BoundedField.cc:103
lsst::geom::Point< int, 2 >
lsst::meas::algorithms.defects.Defects.copy
def copy(self)
Definition: defects.py:251
lsst::geom::Box2I
An integer coordinate rectangle.
Definition: Box.h:55
lsst::meas::algorithms.defects.Defects._metadata
_metadata
Definition: defects.py:99
lsst::meas::algorithms.defects.Defects.__len__
def __len__(self)
Definition: defects.py:135
lsst::meas::algorithms.defects.Defects._bulk_update
_bulk_update
Definition: defects.py:90
lsst::meas::algorithms.defects.Defects.toFitsRegionTable
def toFitsRegionTable(self)
Definition: defects.py:302
lsst::afw::table::CatalogT< BaseRecord >
lsst::geom::Extent< int, 2 >
astshim.fitsChanContinued.iter
def iter(self)
Definition: fitsChanContinued.py:88
lsst::meas::algorithms.defects.Defects.insert
def insert(self, index, value)
Definition: defects.py:217
lsst::meas::algorithms.defects.Defects.writeText
def writeText(self, filename)
Definition: defects.py:463
lsst::afw::geom
Definition: frameSetUtils.h:40