lsst.ip.isr  22.0.1-17-gb1dbe93+5cbcee02db
isrTask.py
Go to the documentation of this file.
1 # This file is part of ip_isr.
2 #
3 # Developed for the LSST Data Management System.
4 # This product includes software developed by the LSST Project
5 # (https://www.lsst.org).
6 # See the COPYRIGHT file at the top-level directory of this distribution
7 # for details of code ownership.
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 GNU General Public License
20 # along with this program. If not, see <https://www.gnu.org/licenses/>.
21 
22 import math
23 import numpy
24 
25 import lsst.geom
26 import lsst.afw.image as afwImage
27 import lsst.afw.math as afwMath
28 import lsst.pex.config as pexConfig
29 import lsst.pipe.base as pipeBase
30 import lsst.pipe.base.connectionTypes as cT
31 
32 from contextlib import contextmanager
33 from lsstDebug import getDebugFrame
34 
35 from lsst.afw.cameraGeom import (PIXELS, FOCAL_PLANE, NullLinearityType,
36  ReadoutCorner)
37 from lsst.afw.display import getDisplay
38 from lsst.afw.geom import Polygon
39 from lsst.daf.persistence import ButlerDataRef
40 from lsst.daf.persistence.butler import NoResults
41 from lsst.meas.algorithms.detection import SourceDetectionTask
42 
43 from . import isrFunctions
44 from . import isrQa
45 from . import linearize
46 from .defects import Defects
47 
48 from .assembleCcdTask import AssembleCcdTask
49 from .crosstalk import CrosstalkTask, CrosstalkCalib
50 from .fringe import FringeTask
51 from .isr import maskNans
52 from .masking import MaskingTask
53 from .overscan import OverscanCorrectionTask
54 from .straylight import StrayLightTask
55 from .vignette import VignetteTask
56 from .ampOffset import AmpOffsetTask
57 from lsst.daf.butler import DimensionGraph
58 
59 
60 __all__ = ["IsrTask", "IsrTaskConfig", "RunIsrTask", "RunIsrConfig"]
61 
62 
63 def crosstalkSourceLookup(datasetType, registry, quantumDataId, collections):
64  """Lookup function to identify crosstalkSource entries.
65 
66  This should return an empty list under most circumstances. Only
67  when inter-chip crosstalk has been identified should this be
68  populated.
69 
70  This will be unused until DM-25348 resolves the quantum graph
71  generation issue.
72 
73  Parameters
74  ----------
75  datasetType : `str`
76  Dataset to lookup.
77  registry : `lsst.daf.butler.Registry`
78  Butler registry to query.
79  quantumDataId : `lsst.daf.butler.ExpandedDataCoordinate`
80  Data id to transform to identify crosstalkSources. The
81  ``detector`` entry will be stripped.
82  collections : `lsst.daf.butler.CollectionSearch`
83  Collections to search through.
84 
85  Returns
86  -------
87  results : `list` [`lsst.daf.butler.DatasetRef`]
88  List of datasets that match the query that will be used as
89  crosstalkSources.
90  """
91  newDataId = quantumDataId.subset(DimensionGraph(registry.dimensions, names=["instrument", "exposure"]))
92  results = list(registry.queryDatasets(datasetType,
93  collections=collections,
94  dataId=newDataId,
95  findFirst=True,
96  ).expanded())
97  return results
98 
99 
100 class IsrTaskConnections(pipeBase.PipelineTaskConnections,
101  dimensions={"instrument", "exposure", "detector"},
102  defaultTemplates={}):
103  ccdExposure = cT.Input(
104  name="raw",
105  doc="Input exposure to process.",
106  storageClass="Exposure",
107  dimensions=["instrument", "exposure", "detector"],
108  )
109  camera = cT.PrerequisiteInput(
110  name="camera",
111  storageClass="Camera",
112  doc="Input camera to construct complete exposures.",
113  dimensions=["instrument"],
114  isCalibration=True,
115  )
116 
117  crosstalk = cT.PrerequisiteInput(
118  name="crosstalk",
119  doc="Input crosstalk object",
120  storageClass="CrosstalkCalib",
121  dimensions=["instrument", "detector"],
122  isCalibration=True,
123  minimum=0, # can fall back to cameraGeom
124  )
125  # TODO: DM-25348. This does not work yet to correctly load
126  # possible crosstalk sources.
127  crosstalkSources = cT.PrerequisiteInput(
128  name="isrOverscanCorrected",
129  doc="Overscan corrected input images.",
130  storageClass="Exposure",
131  dimensions=["instrument", "exposure", "detector"],
132  deferLoad=True,
133  multiple=True,
134  lookupFunction=crosstalkSourceLookup,
135  minimum=0, # not needed for all instruments, no config to control this
136  )
137  bias = cT.PrerequisiteInput(
138  name="bias",
139  doc="Input bias calibration.",
140  storageClass="ExposureF",
141  dimensions=["instrument", "detector"],
142  isCalibration=True,
143  )
144  dark = cT.PrerequisiteInput(
145  name='dark',
146  doc="Input dark calibration.",
147  storageClass="ExposureF",
148  dimensions=["instrument", "detector"],
149  isCalibration=True,
150  )
151  flat = cT.PrerequisiteInput(
152  name="flat",
153  doc="Input flat calibration.",
154  storageClass="ExposureF",
155  dimensions=["instrument", "physical_filter", "detector"],
156  isCalibration=True,
157  )
158  ptc = cT.PrerequisiteInput(
159  name="ptc",
160  doc="Input Photon Transfer Curve dataset",
161  storageClass="PhotonTransferCurveDataset",
162  dimensions=["instrument", "detector"],
163  isCalibration=True,
164  )
165  fringes = cT.PrerequisiteInput(
166  name="fringe",
167  doc="Input fringe calibration.",
168  storageClass="ExposureF",
169  dimensions=["instrument", "physical_filter", "detector"],
170  isCalibration=True,
171  minimum=0, # only needed for some bands, even when enabled
172  )
173  strayLightData = cT.PrerequisiteInput(
174  name='yBackground',
175  doc="Input stray light calibration.",
176  storageClass="StrayLightData",
177  dimensions=["instrument", "physical_filter", "detector"],
178  deferLoad=True,
179  isCalibration=True,
180  minimum=0, # only needed for some bands, even when enabled
181  )
182  bfKernel = cT.PrerequisiteInput(
183  name='bfKernel',
184  doc="Input brighter-fatter kernel.",
185  storageClass="NumpyArray",
186  dimensions=["instrument"],
187  isCalibration=True,
188  minimum=0, # can use either bfKernel or newBFKernel
189  )
190  newBFKernel = cT.PrerequisiteInput(
191  name='brighterFatterKernel',
192  doc="Newer complete kernel + gain solutions.",
193  storageClass="BrighterFatterKernel",
194  dimensions=["instrument", "detector"],
195  isCalibration=True,
196  minimum=0, # can use either bfKernel or newBFKernel
197  )
198  defects = cT.PrerequisiteInput(
199  name='defects',
200  doc="Input defect tables.",
201  storageClass="Defects",
202  dimensions=["instrument", "detector"],
203  isCalibration=True,
204  )
205  linearizer = cT.PrerequisiteInput(
206  name='linearizer',
207  storageClass="Linearizer",
208  doc="Linearity correction calibration.",
209  dimensions=["instrument", "detector"],
210  isCalibration=True,
211  minimum=0, # can fall back to cameraGeom
212  )
213  opticsTransmission = cT.PrerequisiteInput(
214  name="transmission_optics",
215  storageClass="TransmissionCurve",
216  doc="Transmission curve due to the optics.",
217  dimensions=["instrument"],
218  isCalibration=True,
219  )
220  filterTransmission = cT.PrerequisiteInput(
221  name="transmission_filter",
222  storageClass="TransmissionCurve",
223  doc="Transmission curve due to the filter.",
224  dimensions=["instrument", "physical_filter"],
225  isCalibration=True,
226  )
227  sensorTransmission = cT.PrerequisiteInput(
228  name="transmission_sensor",
229  storageClass="TransmissionCurve",
230  doc="Transmission curve due to the sensor.",
231  dimensions=["instrument", "detector"],
232  isCalibration=True,
233  )
234  atmosphereTransmission = cT.PrerequisiteInput(
235  name="transmission_atmosphere",
236  storageClass="TransmissionCurve",
237  doc="Transmission curve due to the atmosphere.",
238  dimensions=["instrument"],
239  isCalibration=True,
240  )
241  illumMaskedImage = cT.PrerequisiteInput(
242  name="illum",
243  doc="Input illumination correction.",
244  storageClass="MaskedImageF",
245  dimensions=["instrument", "physical_filter", "detector"],
246  isCalibration=True,
247  )
248 
249  outputExposure = cT.Output(
250  name='postISRCCD',
251  doc="Output ISR processed exposure.",
252  storageClass="Exposure",
253  dimensions=["instrument", "exposure", "detector"],
254  )
255  preInterpExposure = cT.Output(
256  name='preInterpISRCCD',
257  doc="Output ISR processed exposure, with pixels left uninterpolated.",
258  storageClass="ExposureF",
259  dimensions=["instrument", "exposure", "detector"],
260  )
261  outputOssThumbnail = cT.Output(
262  name="OssThumb",
263  doc="Output Overscan-subtracted thumbnail image.",
264  storageClass="Thumbnail",
265  dimensions=["instrument", "exposure", "detector"],
266  )
267  outputFlattenedThumbnail = cT.Output(
268  name="FlattenedThumb",
269  doc="Output flat-corrected thumbnail image.",
270  storageClass="Thumbnail",
271  dimensions=["instrument", "exposure", "detector"],
272  )
273 
274  def __init__(self, *, config=None):
275  super().__init__(config=config)
276 
277  if config.doBias is not True:
278  self.prerequisiteInputs.discard("bias")
279  if config.doLinearize is not True:
280  self.prerequisiteInputs.discard("linearizer")
281  if config.doCrosstalk is not True:
282  self.inputs.discard("crosstalkSources")
283  self.prerequisiteInputs.discard("crosstalk")
284  if config.doBrighterFatter is not True:
285  self.prerequisiteInputs.discard("bfKernel")
286  self.prerequisiteInputs.discard("newBFKernel")
287  if config.doDefect is not True:
288  self.prerequisiteInputs.discard("defects")
289  if config.doDark is not True:
290  self.prerequisiteInputs.discard("dark")
291  if config.doFlat is not True:
292  self.prerequisiteInputs.discard("flat")
293  if config.doFringe is not True:
294  self.prerequisiteInputs.discard("fringe")
295  if config.doStrayLight is not True:
296  self.prerequisiteInputs.discard("strayLightData")
297  if config.usePtcGains is not True and config.usePtcReadNoise is not True:
298  self.prerequisiteInputs.discard("ptc")
299  if config.doAttachTransmissionCurve is not True:
300  self.prerequisiteInputs.discard("opticsTransmission")
301  self.prerequisiteInputs.discard("filterTransmission")
302  self.prerequisiteInputs.discard("sensorTransmission")
303  self.prerequisiteInputs.discard("atmosphereTransmission")
304  if config.doUseOpticsTransmission is not True:
305  self.prerequisiteInputs.discard("opticsTransmission")
306  if config.doUseFilterTransmission is not True:
307  self.prerequisiteInputs.discard("filterTransmission")
308  if config.doUseSensorTransmission is not True:
309  self.prerequisiteInputs.discard("sensorTransmission")
310  if config.doUseAtmosphereTransmission is not True:
311  self.prerequisiteInputs.discard("atmosphereTransmission")
312  if config.doIlluminationCorrection is not True:
313  self.prerequisiteInputs.discard("illumMaskedImage")
314 
315  if config.doWrite is not True:
316  self.outputs.discard("outputExposure")
317  self.outputs.discard("preInterpExposure")
318  self.outputs.discard("outputFlattenedThumbnail")
319  self.outputs.discard("outputOssThumbnail")
320  if config.doSaveInterpPixels is not True:
321  self.outputs.discard("preInterpExposure")
322  if config.qa.doThumbnailOss is not True:
323  self.outputs.discard("outputOssThumbnail")
324  if config.qa.doThumbnailFlattened is not True:
325  self.outputs.discard("outputFlattenedThumbnail")
326 
327 
328 class IsrTaskConfig(pipeBase.PipelineTaskConfig,
329  pipelineConnections=IsrTaskConnections):
330  """Configuration parameters for IsrTask.
331 
332  Items are grouped in the order in which they are executed by the task.
333  """
334  datasetType = pexConfig.Field(
335  dtype=str,
336  doc="Dataset type for input data; users will typically leave this alone, "
337  "but camera-specific ISR tasks will override it",
338  default="raw",
339  )
340 
341  fallbackFilterName = pexConfig.Field(
342  dtype=str,
343  doc="Fallback default filter name for calibrations.",
344  optional=True
345  )
346  useFallbackDate = pexConfig.Field(
347  dtype=bool,
348  doc="Pass observation date when using fallback filter.",
349  default=False,
350  )
351  expectWcs = pexConfig.Field(
352  dtype=bool,
353  default=True,
354  doc="Expect input science images to have a WCS (set False for e.g. spectrographs)."
355  )
356  fwhm = pexConfig.Field(
357  dtype=float,
358  doc="FWHM of PSF in arcseconds.",
359  default=1.0,
360  )
361  qa = pexConfig.ConfigField(
362  dtype=isrQa.IsrQaConfig,
363  doc="QA related configuration options.",
364  )
365 
366  # Image conversion configuration
367  doConvertIntToFloat = pexConfig.Field(
368  dtype=bool,
369  doc="Convert integer raw images to floating point values?",
370  default=True,
371  )
372 
373  # Saturated pixel handling.
374  doSaturation = pexConfig.Field(
375  dtype=bool,
376  doc="Mask saturated pixels? NB: this is totally independent of the"
377  " interpolation option - this is ONLY setting the bits in the mask."
378  " To have them interpolated make sure doSaturationInterpolation=True",
379  default=True,
380  )
381  saturatedMaskName = pexConfig.Field(
382  dtype=str,
383  doc="Name of mask plane to use in saturation detection and interpolation",
384  default="SAT",
385  )
386  saturation = pexConfig.Field(
387  dtype=float,
388  doc="The saturation level to use if no Detector is present in the Exposure (ignored if NaN)",
389  default=float("NaN"),
390  )
391  growSaturationFootprintSize = pexConfig.Field(
392  dtype=int,
393  doc="Number of pixels by which to grow the saturation footprints",
394  default=1,
395  )
396 
397  # Suspect pixel handling.
398  doSuspect = pexConfig.Field(
399  dtype=bool,
400  doc="Mask suspect pixels?",
401  default=False,
402  )
403  suspectMaskName = pexConfig.Field(
404  dtype=str,
405  doc="Name of mask plane to use for suspect pixels",
406  default="SUSPECT",
407  )
408  numEdgeSuspect = pexConfig.Field(
409  dtype=int,
410  doc="Number of edge pixels to be flagged as untrustworthy.",
411  default=0,
412  )
413  edgeMaskLevel = pexConfig.ChoiceField(
414  dtype=str,
415  doc="Mask edge pixels in which coordinate frame: DETECTOR or AMP?",
416  default="DETECTOR",
417  allowed={
418  'DETECTOR': 'Mask only the edges of the full detector.',
419  'AMP': 'Mask edges of each amplifier.',
420  },
421  )
422 
423  # Initial masking options.
424  doSetBadRegions = pexConfig.Field(
425  dtype=bool,
426  doc="Should we set the level of all BAD patches of the chip to the chip's average value?",
427  default=True,
428  )
429  badStatistic = pexConfig.ChoiceField(
430  dtype=str,
431  doc="How to estimate the average value for BAD regions.",
432  default='MEANCLIP',
433  allowed={
434  "MEANCLIP": "Correct using the (clipped) mean of good data",
435  "MEDIAN": "Correct using the median of the good data",
436  },
437  )
438 
439  # Overscan subtraction configuration.
440  doOverscan = pexConfig.Field(
441  dtype=bool,
442  doc="Do overscan subtraction?",
443  default=True,
444  )
445  overscan = pexConfig.ConfigurableField(
446  target=OverscanCorrectionTask,
447  doc="Overscan subtraction task for image segments.",
448  )
449  overscanFitType = pexConfig.ChoiceField(
450  dtype=str,
451  doc="The method for fitting the overscan bias level.",
452  default='MEDIAN',
453  allowed={
454  "POLY": "Fit ordinary polynomial to the longest axis of the overscan region",
455  "CHEB": "Fit Chebyshev polynomial to the longest axis of the overscan region",
456  "LEG": "Fit Legendre polynomial to the longest axis of the overscan region",
457  "NATURAL_SPLINE": "Fit natural spline to the longest axis of the overscan region",
458  "CUBIC_SPLINE": "Fit cubic spline to the longest axis of the overscan region",
459  "AKIMA_SPLINE": "Fit Akima spline to the longest axis of the overscan region",
460  "MEAN": "Correct using the mean of the overscan region",
461  "MEANCLIP": "Correct using a clipped mean of the overscan region",
462  "MEDIAN": "Correct using the median of the overscan region",
463  "MEDIAN_PER_ROW": "Correct using the median per row of the overscan region",
464  },
465  deprecated=("Please configure overscan via the OverscanCorrectionConfig interface."
466  " This option will no longer be used, and will be removed after v20.")
467  )
468  overscanOrder = pexConfig.Field(
469  dtype=int,
470  doc=("Order of polynomial or to fit if overscan fit type is a polynomial, "
471  "or number of spline knots if overscan fit type is a spline."),
472  default=1,
473  deprecated=("Please configure overscan via the OverscanCorrectionConfig interface."
474  " This option will no longer be used, and will be removed after v20.")
475  )
476  overscanNumSigmaClip = pexConfig.Field(
477  dtype=float,
478  doc="Rejection threshold (sigma) for collapsing overscan before fit",
479  default=3.0,
480  deprecated=("Please configure overscan via the OverscanCorrectionConfig interface."
481  " This option will no longer be used, and will be removed after v20.")
482  )
483  overscanIsInt = pexConfig.Field(
484  dtype=bool,
485  doc="Treat overscan as an integer image for purposes of overscan.FitType=MEDIAN"
486  " and overscan.FitType=MEDIAN_PER_ROW.",
487  default=True,
488  deprecated=("Please configure overscan via the OverscanCorrectionConfig interface."
489  " This option will no longer be used, and will be removed after v20.")
490  )
491  # These options do not get deprecated, as they define how we slice up the
492  # image data.
493  overscanNumLeadingColumnsToSkip = pexConfig.Field(
494  dtype=int,
495  doc="Number of columns to skip in overscan, i.e. those closest to amplifier",
496  default=0,
497  )
498  overscanNumTrailingColumnsToSkip = pexConfig.Field(
499  dtype=int,
500  doc="Number of columns to skip in overscan, i.e. those farthest from amplifier",
501  default=0,
502  )
503  overscanMaxDev = pexConfig.Field(
504  dtype=float,
505  doc="Maximum deviation from the median for overscan",
506  default=1000.0, check=lambda x: x > 0
507  )
508  overscanBiasJump = pexConfig.Field(
509  dtype=bool,
510  doc="Fit the overscan in a piecewise-fashion to correct for bias jumps?",
511  default=False,
512  )
513  overscanBiasJumpKeyword = pexConfig.Field(
514  dtype=str,
515  doc="Header keyword containing information about devices.",
516  default="NO_SUCH_KEY",
517  )
518  overscanBiasJumpDevices = pexConfig.ListField(
519  dtype=str,
520  doc="List of devices that need piecewise overscan correction.",
521  default=(),
522  )
523  overscanBiasJumpLocation = pexConfig.Field(
524  dtype=int,
525  doc="Location of bias jump along y-axis.",
526  default=0,
527  )
528 
529  # Amplifier to CCD assembly configuration
530  doAssembleCcd = pexConfig.Field(
531  dtype=bool,
532  default=True,
533  doc="Assemble amp-level exposures into a ccd-level exposure?"
534  )
535  assembleCcd = pexConfig.ConfigurableField(
536  target=AssembleCcdTask,
537  doc="CCD assembly task",
538  )
539 
540  # General calibration configuration.
541  doAssembleIsrExposures = pexConfig.Field(
542  dtype=bool,
543  default=False,
544  doc="Assemble amp-level calibration exposures into ccd-level exposure?"
545  )
546  doTrimToMatchCalib = pexConfig.Field(
547  dtype=bool,
548  default=False,
549  doc="Trim raw data to match calibration bounding boxes?"
550  )
551 
552  # Bias subtraction.
553  doBias = pexConfig.Field(
554  dtype=bool,
555  doc="Apply bias frame correction?",
556  default=True,
557  )
558  biasDataProductName = pexConfig.Field(
559  dtype=str,
560  doc="Name of the bias data product",
561  default="bias",
562  )
563  doBiasBeforeOverscan = pexConfig.Field(
564  dtype=bool,
565  doc="Reverse order of overscan and bias correction.",
566  default=False
567  )
568 
569  # Variance construction
570  doVariance = pexConfig.Field(
571  dtype=bool,
572  doc="Calculate variance?",
573  default=True
574  )
575  gain = pexConfig.Field(
576  dtype=float,
577  doc="The gain to use if no Detector is present in the Exposure (ignored if NaN)",
578  default=float("NaN"),
579  )
580  readNoise = pexConfig.Field(
581  dtype=float,
582  doc="The read noise to use if no Detector is present in the Exposure",
583  default=0.0,
584  )
585  doEmpiricalReadNoise = pexConfig.Field(
586  dtype=bool,
587  default=False,
588  doc="Calculate empirical read noise instead of value from AmpInfo data?"
589  )
590  usePtcReadNoise = pexConfig.Field(
591  dtype=bool,
592  default=False,
593  doc="Use readnoise values from the Photon Transfer Curve?"
594  )
595  maskNegativeVariance = pexConfig.Field(
596  dtype=bool,
597  default=True,
598  doc="Mask pixels that claim a negative variance? This likely indicates a failure "
599  "in the measurement of the overscan at an edge due to the data falling off faster "
600  "than the overscan model can account for it."
601  )
602  negativeVarianceMaskName = pexConfig.Field(
603  dtype=str,
604  default="BAD",
605  doc="Mask plane to use to mark pixels with negative variance, if `maskNegativeVariance` is True.",
606  )
607  # Linearization.
608  doLinearize = pexConfig.Field(
609  dtype=bool,
610  doc="Correct for nonlinearity of the detector's response?",
611  default=True,
612  )
613 
614  # Crosstalk.
615  doCrosstalk = pexConfig.Field(
616  dtype=bool,
617  doc="Apply intra-CCD crosstalk correction?",
618  default=False,
619  )
620  doCrosstalkBeforeAssemble = pexConfig.Field(
621  dtype=bool,
622  doc="Apply crosstalk correction before CCD assembly, and before trimming?",
623  default=False,
624  )
625  crosstalk = pexConfig.ConfigurableField(
626  target=CrosstalkTask,
627  doc="Intra-CCD crosstalk correction",
628  )
629 
630  # Masking options.
631  doDefect = pexConfig.Field(
632  dtype=bool,
633  doc="Apply correction for CCD defects, e.g. hot pixels?",
634  default=True,
635  )
636  doNanMasking = pexConfig.Field(
637  dtype=bool,
638  doc="Mask non-finite (NAN, inf) pixels?",
639  default=True,
640  )
641  doWidenSaturationTrails = pexConfig.Field(
642  dtype=bool,
643  doc="Widen bleed trails based on their width?",
644  default=True
645  )
646 
647  # Brighter-Fatter correction.
648  doBrighterFatter = pexConfig.Field(
649  dtype=bool,
650  default=False,
651  doc="Apply the brighter-fatter correction?"
652  )
653  brighterFatterLevel = pexConfig.ChoiceField(
654  dtype=str,
655  default="DETECTOR",
656  doc="The level at which to correct for brighter-fatter.",
657  allowed={
658  "AMP": "Every amplifier treated separately.",
659  "DETECTOR": "One kernel per detector",
660  }
661  )
662  brighterFatterMaxIter = pexConfig.Field(
663  dtype=int,
664  default=10,
665  doc="Maximum number of iterations for the brighter-fatter correction"
666  )
667  brighterFatterThreshold = pexConfig.Field(
668  dtype=float,
669  default=1000,
670  doc="Threshold used to stop iterating the brighter-fatter correction. It is the "
671  "absolute value of the difference between the current corrected image and the one "
672  "from the previous iteration summed over all the pixels."
673  )
674  brighterFatterApplyGain = pexConfig.Field(
675  dtype=bool,
676  default=True,
677  doc="Should the gain be applied when applying the brighter-fatter correction?"
678  )
679  brighterFatterMaskListToInterpolate = pexConfig.ListField(
680  dtype=str,
681  doc="List of mask planes that should be interpolated over when applying the brighter-fatter "
682  "correction.",
683  default=["SAT", "BAD", "NO_DATA", "UNMASKEDNAN"],
684  )
685  brighterFatterMaskGrowSize = pexConfig.Field(
686  dtype=int,
687  default=0,
688  doc="Number of pixels to grow the masks listed in config.brighterFatterMaskListToInterpolate "
689  "when brighter-fatter correction is applied."
690  )
691 
692  # Dark subtraction.
693  doDark = pexConfig.Field(
694  dtype=bool,
695  doc="Apply dark frame correction?",
696  default=True,
697  )
698  darkDataProductName = pexConfig.Field(
699  dtype=str,
700  doc="Name of the dark data product",
701  default="dark",
702  )
703 
704  # Camera-specific stray light removal.
705  doStrayLight = pexConfig.Field(
706  dtype=bool,
707  doc="Subtract stray light in the y-band (due to encoder LEDs)?",
708  default=False,
709  )
710  strayLight = pexConfig.ConfigurableField(
711  target=StrayLightTask,
712  doc="y-band stray light correction"
713  )
714 
715  # Flat correction.
716  doFlat = pexConfig.Field(
717  dtype=bool,
718  doc="Apply flat field correction?",
719  default=True,
720  )
721  flatDataProductName = pexConfig.Field(
722  dtype=str,
723  doc="Name of the flat data product",
724  default="flat",
725  )
726  flatScalingType = pexConfig.ChoiceField(
727  dtype=str,
728  doc="The method for scaling the flat on the fly.",
729  default='USER',
730  allowed={
731  "USER": "Scale by flatUserScale",
732  "MEAN": "Scale by the inverse of the mean",
733  "MEDIAN": "Scale by the inverse of the median",
734  },
735  )
736  flatUserScale = pexConfig.Field(
737  dtype=float,
738  doc="If flatScalingType is 'USER' then scale flat by this amount; ignored otherwise",
739  default=1.0,
740  )
741  doTweakFlat = pexConfig.Field(
742  dtype=bool,
743  doc="Tweak flats to match observed amplifier ratios?",
744  default=False
745  )
746 
747  # Amplifier normalization based on gains instead of using flats
748  # configuration.
749  doApplyGains = pexConfig.Field(
750  dtype=bool,
751  doc="Correct the amplifiers for their gains instead of applying flat correction",
752  default=False,
753  )
754  usePtcGains = pexConfig.Field(
755  dtype=bool,
756  doc="Use the gain values from the Photon Transfer Curve?",
757  default=False,
758  )
759  normalizeGains = pexConfig.Field(
760  dtype=bool,
761  doc="Normalize all the amplifiers in each CCD to have the same median value.",
762  default=False,
763  )
764 
765  # Fringe correction.
766  doFringe = pexConfig.Field(
767  dtype=bool,
768  doc="Apply fringe correction?",
769  default=True,
770  )
771  fringe = pexConfig.ConfigurableField(
772  target=FringeTask,
773  doc="Fringe subtraction task",
774  )
775  fringeAfterFlat = pexConfig.Field(
776  dtype=bool,
777  doc="Do fringe subtraction after flat-fielding?",
778  default=True,
779  )
780 
781  # Amp offset correction.
782  doAmpOffset = pexConfig.Field(
783  doc="Calculate and apply amp offset corrections?",
784  dtype=bool,
785  default=False,
786  )
787  ampOffset = pexConfig.ConfigurableField(
788  doc="Amp offset correction task.",
789  target=AmpOffsetTask,
790  )
791 
792  # Initial CCD-level background statistics options.
793  doMeasureBackground = pexConfig.Field(
794  dtype=bool,
795  doc="Measure the background level on the reduced image?",
796  default=False,
797  )
798 
799  # Camera-specific masking configuration.
800  doCameraSpecificMasking = pexConfig.Field(
801  dtype=bool,
802  doc="Mask camera-specific bad regions?",
803  default=False,
804  )
805  masking = pexConfig.ConfigurableField(
806  target=MaskingTask,
807  doc="Masking task."
808  )
809 
810  # Interpolation options.
811  doInterpolate = pexConfig.Field(
812  dtype=bool,
813  doc="Interpolate masked pixels?",
814  default=True,
815  )
816  doSaturationInterpolation = pexConfig.Field(
817  dtype=bool,
818  doc="Perform interpolation over pixels masked as saturated?"
819  " NB: This is independent of doSaturation; if that is False this plane"
820  " will likely be blank, resulting in a no-op here.",
821  default=True,
822  )
823  doNanInterpolation = pexConfig.Field(
824  dtype=bool,
825  doc="Perform interpolation over pixels masked as NaN?"
826  " NB: This is independent of doNanMasking; if that is False this plane"
827  " will likely be blank, resulting in a no-op here.",
828  default=True,
829  )
830  doNanInterpAfterFlat = pexConfig.Field(
831  dtype=bool,
832  doc=("If True, ensure we interpolate NaNs after flat-fielding, even if we "
833  "also have to interpolate them before flat-fielding."),
834  default=False,
835  )
836  maskListToInterpolate = pexConfig.ListField(
837  dtype=str,
838  doc="List of mask planes that should be interpolated.",
839  default=['SAT', 'BAD'],
840  )
841  doSaveInterpPixels = pexConfig.Field(
842  dtype=bool,
843  doc="Save a copy of the pre-interpolated pixel values?",
844  default=False,
845  )
846 
847  # Default photometric calibration options.
848  fluxMag0T1 = pexConfig.DictField(
849  keytype=str,
850  itemtype=float,
851  doc="The approximate flux of a zero-magnitude object in a one-second exposure, per filter.",
852  default=dict((f, pow(10.0, 0.4*m)) for f, m in (("Unknown", 28.0),
853  ))
854  )
855  defaultFluxMag0T1 = pexConfig.Field(
856  dtype=float,
857  doc="Default value for fluxMag0T1 (for an unrecognized filter).",
858  default=pow(10.0, 0.4*28.0)
859  )
860 
861  # Vignette correction configuration.
862  doVignette = pexConfig.Field(
863  dtype=bool,
864  doc="Apply vignetting parameters?",
865  default=False,
866  )
867  vignette = pexConfig.ConfigurableField(
868  target=VignetteTask,
869  doc="Vignetting task.",
870  )
871 
872  # Transmission curve configuration.
873  doAttachTransmissionCurve = pexConfig.Field(
874  dtype=bool,
875  default=False,
876  doc="Construct and attach a wavelength-dependent throughput curve for this CCD image?"
877  )
878  doUseOpticsTransmission = pexConfig.Field(
879  dtype=bool,
880  default=True,
881  doc="Load and use transmission_optics (if doAttachTransmissionCurve is True)?"
882  )
883  doUseFilterTransmission = pexConfig.Field(
884  dtype=bool,
885  default=True,
886  doc="Load and use transmission_filter (if doAttachTransmissionCurve is True)?"
887  )
888  doUseSensorTransmission = pexConfig.Field(
889  dtype=bool,
890  default=True,
891  doc="Load and use transmission_sensor (if doAttachTransmissionCurve is True)?"
892  )
893  doUseAtmosphereTransmission = pexConfig.Field(
894  dtype=bool,
895  default=True,
896  doc="Load and use transmission_atmosphere (if doAttachTransmissionCurve is True)?"
897  )
898 
899  # Illumination correction.
900  doIlluminationCorrection = pexConfig.Field(
901  dtype=bool,
902  default=False,
903  doc="Perform illumination correction?"
904  )
905  illuminationCorrectionDataProductName = pexConfig.Field(
906  dtype=str,
907  doc="Name of the illumination correction data product.",
908  default="illumcor",
909  )
910  illumScale = pexConfig.Field(
911  dtype=float,
912  doc="Scale factor for the illumination correction.",
913  default=1.0,
914  )
915  illumFilters = pexConfig.ListField(
916  dtype=str,
917  default=[],
918  doc="Only perform illumination correction for these filters."
919  )
920 
921  # Write the outputs to disk. If ISR is run as a subtask, this may not
922  # be needed.
923  doWrite = pexConfig.Field(
924  dtype=bool,
925  doc="Persist postISRCCD?",
926  default=True,
927  )
928 
929  def validate(self):
930  super().validate()
931  if self.doFlatdoFlat and self.doApplyGainsdoApplyGains:
932  raise ValueError("You may not specify both doFlat and doApplyGains")
933  if self.doBiasBeforeOverscandoBiasBeforeOverscan and self.doTrimToMatchCalibdoTrimToMatchCalib:
934  raise ValueError("You may not specify both doBiasBeforeOverscan and doTrimToMatchCalib")
935  if self.doSaturationInterpolationdoSaturationInterpolation and self.saturatedMaskNamesaturatedMaskName not in self.maskListToInterpolatemaskListToInterpolate:
936  self.maskListToInterpolatemaskListToInterpolate.append(self.saturatedMaskNamesaturatedMaskName)
937  if not self.doSaturationInterpolationdoSaturationInterpolation and self.saturatedMaskNamesaturatedMaskName in self.maskListToInterpolatemaskListToInterpolate:
938  self.maskListToInterpolatemaskListToInterpolate.remove(self.saturatedMaskNamesaturatedMaskName)
939  if self.doNanInterpolationdoNanInterpolation and "UNMASKEDNAN" not in self.maskListToInterpolatemaskListToInterpolate:
940  self.maskListToInterpolatemaskListToInterpolate.append("UNMASKEDNAN")
941 
942 
943 class IsrTask(pipeBase.PipelineTask, pipeBase.CmdLineTask):
944  """Apply common instrument signature correction algorithms to a raw frame.
945 
946  The process for correcting imaging data is very similar from
947  camera to camera. This task provides a vanilla implementation of
948  doing these corrections, including the ability to turn certain
949  corrections off if they are not needed. The inputs to the primary
950  method, `run()`, are a raw exposure to be corrected and the
951  calibration data products. The raw input is a single chip sized
952  mosaic of all amps including overscans and other non-science
953  pixels. The method `runDataRef()` identifies and defines the
954  calibration data products, and is intended for use by a
955  `lsst.pipe.base.cmdLineTask.CmdLineTask` and takes as input only a
956  `daf.persistence.butlerSubset.ButlerDataRef`. This task may be
957  subclassed for different camera, although the most camera specific
958  methods have been split into subtasks that can be redirected
959  appropriately.
960 
961  The __init__ method sets up the subtasks for ISR processing, using
962  the defaults from `lsst.ip.isr`.
963 
964  Parameters
965  ----------
966  args : `list`
967  Positional arguments passed to the Task constructor.
968  None used at this time.
969  kwargs : `dict`, optional
970  Keyword arguments passed on to the Task constructor.
971  None used at this time.
972  """
973  ConfigClass = IsrTaskConfig
974  _DefaultName = "isr"
975 
976  def __init__(self, **kwargs):
977  super().__init__(**kwargs)
978  self.makeSubtask("assembleCcd")
979  self.makeSubtask("crosstalk")
980  self.makeSubtask("strayLight")
981  self.makeSubtask("fringe")
982  self.makeSubtask("masking")
983  self.makeSubtask("overscan")
984  self.makeSubtask("vignette")
985  self.makeSubtask("ampOffset")
986 
987  def runQuantum(self, butlerQC, inputRefs, outputRefs):
988  inputs = butlerQC.get(inputRefs)
989 
990  try:
991  inputs['detectorNum'] = inputRefs.ccdExposure.dataId['detector']
992  except Exception as e:
993  raise ValueError("Failure to find valid detectorNum value for Dataset %s: %s." %
994  (inputRefs, e))
995 
996  inputs['isGen3'] = True
997 
998  detector = inputs['ccdExposure'].getDetector()
999 
1000  if self.config.doCrosstalk is True:
1001  # Crosstalk sources need to be defined by the pipeline
1002  # yaml if they exist.
1003  if 'crosstalk' in inputs and inputs['crosstalk'] is not None:
1004  if not isinstance(inputs['crosstalk'], CrosstalkCalib):
1005  inputs['crosstalk'] = CrosstalkCalib.fromTable(inputs['crosstalk'])
1006  else:
1007  coeffVector = (self.config.crosstalk.crosstalkValues
1008  if self.config.crosstalk.useConfigCoefficients else None)
1009  crosstalkCalib = CrosstalkCalib().fromDetector(detector, coeffVector=coeffVector)
1010  inputs['crosstalk'] = crosstalkCalib
1011  if inputs['crosstalk'].interChip and len(inputs['crosstalk'].interChip) > 0:
1012  if 'crosstalkSources' not in inputs:
1013  self.log.warning("No crosstalkSources found for chip with interChip terms!")
1014 
1015  if self.doLinearizedoLinearize(detector) is True:
1016  if 'linearizer' in inputs:
1017  if isinstance(inputs['linearizer'], dict):
1018  linearizer = linearize.Linearizer(detector=detector, log=self.log)
1019  linearizer.fromYaml(inputs['linearizer'])
1020  self.log.warning("Dictionary linearizers will be deprecated in DM-28741.")
1021  elif isinstance(inputs['linearizer'], numpy.ndarray):
1022  linearizer = linearize.Linearizer(table=inputs.get('linearizer', None),
1023  detector=detector,
1024  log=self.log)
1025  self.log.warning("Bare lookup table linearizers will be deprecated in DM-28741.")
1026  else:
1027  linearizer = inputs['linearizer']
1028  linearizer.log = self.log
1029  inputs['linearizer'] = linearizer
1030  else:
1031  inputs['linearizer'] = linearize.Linearizer(detector=detector, log=self.log)
1032  self.log.warning("Constructing linearizer from cameraGeom information.")
1033 
1034  if self.config.doDefect is True:
1035  if "defects" in inputs and inputs['defects'] is not None:
1036  # defects is loaded as a BaseCatalog with columns
1037  # x0, y0, width, height. Masking expects a list of defects
1038  # defined by their bounding box
1039  if not isinstance(inputs["defects"], Defects):
1040  inputs["defects"] = Defects.fromTable(inputs["defects"])
1041 
1042  # Load the correct style of brighter-fatter kernel, and repack
1043  # the information as a numpy array.
1044  if self.config.doBrighterFatter:
1045  brighterFatterKernel = inputs.pop('newBFKernel', None)
1046  if brighterFatterKernel is None:
1047  brighterFatterKernel = inputs.get('bfKernel', None)
1048 
1049  if brighterFatterKernel is not None and not isinstance(brighterFatterKernel, numpy.ndarray):
1050  # This is a ISR calib kernel
1051  detName = detector.getName()
1052  level = brighterFatterKernel.level
1053 
1054  # This is expected to be a dictionary of amp-wise gains.
1055  inputs['bfGains'] = brighterFatterKernel.gain
1056  if self.config.brighterFatterLevel == 'DETECTOR':
1057  if level == 'DETECTOR':
1058  if detName in brighterFatterKernel.detKernels:
1059  inputs['bfKernel'] = brighterFatterKernel.detKernels[detName]
1060  else:
1061  raise RuntimeError("Failed to extract kernel from new-style BF kernel.")
1062  elif level == 'AMP':
1063  self.log.warning("Making DETECTOR level kernel from AMP based brighter "
1064  "fatter kernels.")
1065  brighterFatterKernel.makeDetectorKernelFromAmpwiseKernels(detName)
1066  inputs['bfKernel'] = brighterFatterKernel.detKernels[detName]
1067  elif self.config.brighterFatterLevel == 'AMP':
1068  raise NotImplementedError("Per-amplifier brighter-fatter correction not implemented")
1069 
1070  if self.config.doFringe is True and self.fringe.checkFilter(inputs['ccdExposure']):
1071  expId = inputs['ccdExposure'].getInfo().getVisitInfo().getExposureId()
1072  inputs['fringes'] = self.fringe.loadFringes(inputs['fringes'],
1073  expId=expId,
1074  assembler=self.assembleCcd
1075  if self.config.doAssembleIsrExposures else None)
1076  else:
1077  inputs['fringes'] = pipeBase.Struct(fringes=None)
1078 
1079  if self.config.doStrayLight is True and self.strayLight.checkFilter(inputs['ccdExposure']):
1080  if 'strayLightData' not in inputs:
1081  inputs['strayLightData'] = None
1082 
1083  outputs = self.runrun(**inputs)
1084  butlerQC.put(outputs, outputRefs)
1085 
1086  def readIsrData(self, dataRef, rawExposure):
1087  """Retrieve necessary frames for instrument signature removal.
1088 
1089  Pre-fetching all required ISR data products limits the IO
1090  required by the ISR. Any conflict between the calibration data
1091  available and that needed for ISR is also detected prior to
1092  doing processing, allowing it to fail quickly.
1093 
1094  Parameters
1095  ----------
1096  dataRef : `daf.persistence.butlerSubset.ButlerDataRef`
1097  Butler reference of the detector data to be processed
1098  rawExposure : `afw.image.Exposure`
1099  The raw exposure that will later be corrected with the
1100  retrieved calibration data; should not be modified in this
1101  method.
1102 
1103  Returns
1104  -------
1105  result : `lsst.pipe.base.Struct`
1106  Result struct with components (which may be `None`):
1107  - ``bias``: bias calibration frame (`afw.image.Exposure`)
1108  - ``linearizer``: functor for linearization
1109  (`ip.isr.linearize.LinearizeBase`)
1110  - ``crosstalkSources``: list of possible crosstalk sources (`list`)
1111  - ``dark``: dark calibration frame (`afw.image.Exposure`)
1112  - ``flat``: flat calibration frame (`afw.image.Exposure`)
1113  - ``bfKernel``: Brighter-Fatter kernel (`numpy.ndarray`)
1114  - ``defects``: list of defects (`lsst.ip.isr.Defects`)
1115  - ``fringes``: `lsst.pipe.base.Struct` with components:
1116  - ``fringes``: fringe calibration frame (`afw.image.Exposure`)
1117  - ``seed``: random seed derived from the ccdExposureId for random
1118  number generator (`uint32`).
1119  - ``opticsTransmission``: `lsst.afw.image.TransmissionCurve`
1120  A ``TransmissionCurve`` that represents the throughput of the
1121  optics, to be evaluated in focal-plane coordinates.
1122  - ``filterTransmission`` : `lsst.afw.image.TransmissionCurve`
1123  A ``TransmissionCurve`` that represents the throughput of the
1124  filter itself, to be evaluated in focal-plane coordinates.
1125  - ``sensorTransmission`` : `lsst.afw.image.TransmissionCurve`
1126  A ``TransmissionCurve`` that represents the throughput of the
1127  sensor itself, to be evaluated in post-assembly trimmed
1128  detector coordinates.
1129  - ``atmosphereTransmission`` : `lsst.afw.image.TransmissionCurve`
1130  A ``TransmissionCurve`` that represents the throughput of the
1131  atmosphere, assumed to be spatially constant.
1132  - ``strayLightData`` : `object`
1133  An opaque object containing calibration information for
1134  stray-light correction. If `None`, no correction will be
1135  performed.
1136  - ``illumMaskedImage`` : illumination correction image
1137  (`lsst.afw.image.MaskedImage`)
1138 
1139  Raises
1140  ------
1141  NotImplementedError :
1142  Raised if a per-amplifier brighter-fatter kernel is requested by
1143  the configuration.
1144  """
1145  try:
1146  dateObs = rawExposure.getInfo().getVisitInfo().getDate()
1147  dateObs = dateObs.toPython().isoformat()
1148  except RuntimeError:
1149  self.log.warning("Unable to identify dateObs for rawExposure.")
1150  dateObs = None
1151 
1152  ccd = rawExposure.getDetector()
1153  filterLabel = rawExposure.getFilterLabel()
1154  physicalFilter = isrFunctions.getPhysicalFilter(filterLabel, self.log)
1155  rawExposure.mask.addMaskPlane("UNMASKEDNAN") # needed to match pre DM-15862 processing.
1156  biasExposure = (self.getIsrExposuregetIsrExposure(dataRef, self.config.biasDataProductName)
1157  if self.config.doBias else None)
1158  # immediate=True required for functors and linearizers are functors
1159  # see ticket DM-6515
1160  linearizer = (dataRef.get("linearizer", immediate=True)
1161  if self.doLinearizedoLinearize(ccd) else None)
1162  if linearizer is not None and not isinstance(linearizer, numpy.ndarray):
1163  linearizer.log = self.log
1164  if isinstance(linearizer, numpy.ndarray):
1165  linearizer = linearize.Linearizer(table=linearizer, detector=ccd)
1166 
1167  crosstalkCalib = None
1168  if self.config.doCrosstalk:
1169  try:
1170  crosstalkCalib = dataRef.get("crosstalk", immediate=True)
1171  except NoResults:
1172  coeffVector = (self.config.crosstalk.crosstalkValues
1173  if self.config.crosstalk.useConfigCoefficients else None)
1174  crosstalkCalib = CrosstalkCalib().fromDetector(ccd, coeffVector=coeffVector)
1175  crosstalkSources = (self.crosstalk.prepCrosstalk(dataRef, crosstalkCalib)
1176  if self.config.doCrosstalk else None)
1177 
1178  darkExposure = (self.getIsrExposuregetIsrExposure(dataRef, self.config.darkDataProductName)
1179  if self.config.doDark else None)
1180  flatExposure = (self.getIsrExposuregetIsrExposure(dataRef, self.config.flatDataProductName,
1181  dateObs=dateObs)
1182  if self.config.doFlat else None)
1183 
1184  brighterFatterKernel = None
1185  brighterFatterGains = None
1186  if self.config.doBrighterFatter is True:
1187  try:
1188  # Use the new-style cp_pipe version of the kernel if it exists
1189  # If using a new-style kernel, always use the self-consistent
1190  # gains, i.e. the ones inside the kernel object itself
1191  brighterFatterKernel = dataRef.get("brighterFatterKernel")
1192  brighterFatterGains = brighterFatterKernel.gain
1193  self.log.info("New style brighter-fatter kernel (brighterFatterKernel) loaded")
1194  except NoResults:
1195  try: # Fall back to the old-style numpy-ndarray style kernel if necessary.
1196  brighterFatterKernel = dataRef.get("bfKernel")
1197  self.log.info("Old style brighter-fatter kernel (bfKernel) loaded")
1198  except NoResults:
1199  brighterFatterKernel = None
1200  if brighterFatterKernel is not None and not isinstance(brighterFatterKernel, numpy.ndarray):
1201  # If the kernel is not an ndarray, it's the cp_pipe version
1202  # so extract the kernel for this detector, or raise an error
1203  if self.config.brighterFatterLevel == 'DETECTOR':
1204  if brighterFatterKernel.detKernels:
1205  brighterFatterKernel = brighterFatterKernel.detKernels[ccd.getName()]
1206  else:
1207  raise RuntimeError("Failed to extract kernel from new-style BF kernel.")
1208  else:
1209  # TODO DM-15631 for implementing this
1210  raise NotImplementedError("Per-amplifier brighter-fatter correction not implemented")
1211 
1212  defectList = (dataRef.get("defects")
1213  if self.config.doDefect else None)
1214  expId = rawExposure.getInfo().getVisitInfo().getExposureId()
1215  fringeStruct = (self.fringe.readFringes(dataRef, expId=expId, assembler=self.assembleCcd
1216  if self.config.doAssembleIsrExposures else None)
1217  if self.config.doFringe and self.fringe.checkFilter(rawExposure)
1218  else pipeBase.Struct(fringes=None))
1219 
1220  if self.config.doAttachTransmissionCurve:
1221  opticsTransmission = (dataRef.get("transmission_optics")
1222  if self.config.doUseOpticsTransmission else None)
1223  filterTransmission = (dataRef.get("transmission_filter")
1224  if self.config.doUseFilterTransmission else None)
1225  sensorTransmission = (dataRef.get("transmission_sensor")
1226  if self.config.doUseSensorTransmission else None)
1227  atmosphereTransmission = (dataRef.get("transmission_atmosphere")
1228  if self.config.doUseAtmosphereTransmission else None)
1229  else:
1230  opticsTransmission = None
1231  filterTransmission = None
1232  sensorTransmission = None
1233  atmosphereTransmission = None
1234 
1235  if self.config.doStrayLight:
1236  strayLightData = self.strayLight.readIsrData(dataRef, rawExposure)
1237  else:
1238  strayLightData = None
1239 
1240  illumMaskedImage = (self.getIsrExposuregetIsrExposure(dataRef,
1241  self.config.illuminationCorrectionDataProductName).getMaskedImage()
1242  if (self.config.doIlluminationCorrection
1243  and physicalFilter in self.config.illumFilters)
1244  else None)
1245 
1246  # Struct should include only kwargs to run()
1247  return pipeBase.Struct(bias=biasExposure,
1248  linearizer=linearizer,
1249  crosstalk=crosstalkCalib,
1250  crosstalkSources=crosstalkSources,
1251  dark=darkExposure,
1252  flat=flatExposure,
1253  bfKernel=brighterFatterKernel,
1254  bfGains=brighterFatterGains,
1255  defects=defectList,
1256  fringes=fringeStruct,
1257  opticsTransmission=opticsTransmission,
1258  filterTransmission=filterTransmission,
1259  sensorTransmission=sensorTransmission,
1260  atmosphereTransmission=atmosphereTransmission,
1261  strayLightData=strayLightData,
1262  illumMaskedImage=illumMaskedImage
1263  )
1264 
1265  @pipeBase.timeMethod
1266  def run(self, ccdExposure, *, camera=None, bias=None, linearizer=None,
1267  crosstalk=None, crosstalkSources=None,
1268  dark=None, flat=None, ptc=None, bfKernel=None, bfGains=None, defects=None,
1269  fringes=pipeBase.Struct(fringes=None), opticsTransmission=None, filterTransmission=None,
1270  sensorTransmission=None, atmosphereTransmission=None,
1271  detectorNum=None, strayLightData=None, illumMaskedImage=None,
1272  isGen3=False,
1273  ):
1274  """Perform instrument signature removal on an exposure.
1275 
1276  Steps included in the ISR processing, in order performed, are:
1277  - saturation and suspect pixel masking
1278  - overscan subtraction
1279  - CCD assembly of individual amplifiers
1280  - bias subtraction
1281  - variance image construction
1282  - linearization of non-linear response
1283  - crosstalk masking
1284  - brighter-fatter correction
1285  - dark subtraction
1286  - fringe correction
1287  - stray light subtraction
1288  - flat correction
1289  - masking of known defects and camera specific features
1290  - vignette calculation
1291  - appending transmission curve and distortion model
1292 
1293  Parameters
1294  ----------
1295  ccdExposure : `lsst.afw.image.Exposure`
1296  The raw exposure that is to be run through ISR. The
1297  exposure is modified by this method.
1298  camera : `lsst.afw.cameraGeom.Camera`, optional
1299  The camera geometry for this exposure. Required if
1300  one or more of ``ccdExposure``, ``bias``, ``dark``, or
1301  ``flat`` does not have an associated detector.
1302  bias : `lsst.afw.image.Exposure`, optional
1303  Bias calibration frame.
1304  linearizer : `lsst.ip.isr.linearize.LinearizeBase`, optional
1305  Functor for linearization.
1306  crosstalk : `lsst.ip.isr.crosstalk.CrosstalkCalib`, optional
1307  Calibration for crosstalk.
1308  crosstalkSources : `list`, optional
1309  List of possible crosstalk sources.
1310  dark : `lsst.afw.image.Exposure`, optional
1311  Dark calibration frame.
1312  flat : `lsst.afw.image.Exposure`, optional
1313  Flat calibration frame.
1314  ptc : `lsst.ip.isr.PhotonTransferCurveDataset`, optional
1315  Photon transfer curve dataset, with, e.g., gains
1316  and read noise.
1317  bfKernel : `numpy.ndarray`, optional
1318  Brighter-fatter kernel.
1319  bfGains : `dict` of `float`, optional
1320  Gains used to override the detector's nominal gains for the
1321  brighter-fatter correction. A dict keyed by amplifier name for
1322  the detector in question.
1323  defects : `lsst.ip.isr.Defects`, optional
1324  List of defects.
1325  fringes : `lsst.pipe.base.Struct`, optional
1326  Struct containing the fringe correction data, with
1327  elements:
1328  - ``fringes``: fringe calibration frame (`afw.image.Exposure`)
1329  - ``seed``: random seed derived from the ccdExposureId for random
1330  number generator (`uint32`)
1331  opticsTransmission: `lsst.afw.image.TransmissionCurve`, optional
1332  A ``TransmissionCurve`` that represents the throughput of the,
1333  optics, to be evaluated in focal-plane coordinates.
1334  filterTransmission : `lsst.afw.image.TransmissionCurve`
1335  A ``TransmissionCurve`` that represents the throughput of the
1336  filter itself, to be evaluated in focal-plane coordinates.
1337  sensorTransmission : `lsst.afw.image.TransmissionCurve`
1338  A ``TransmissionCurve`` that represents the throughput of the
1339  sensor itself, to be evaluated in post-assembly trimmed detector
1340  coordinates.
1341  atmosphereTransmission : `lsst.afw.image.TransmissionCurve`
1342  A ``TransmissionCurve`` that represents the throughput of the
1343  atmosphere, assumed to be spatially constant.
1344  detectorNum : `int`, optional
1345  The integer number for the detector to process.
1346  isGen3 : bool, optional
1347  Flag this call to run() as using the Gen3 butler environment.
1348  strayLightData : `object`, optional
1349  Opaque object containing calibration information for stray-light
1350  correction. If `None`, no correction will be performed.
1351  illumMaskedImage : `lsst.afw.image.MaskedImage`, optional
1352  Illumination correction image.
1353 
1354  Returns
1355  -------
1356  result : `lsst.pipe.base.Struct`
1357  Result struct with component:
1358  - ``exposure`` : `afw.image.Exposure`
1359  The fully ISR corrected exposure.
1360  - ``outputExposure`` : `afw.image.Exposure`
1361  An alias for `exposure`
1362  - ``ossThumb`` : `numpy.ndarray`
1363  Thumbnail image of the exposure after overscan subtraction.
1364  - ``flattenedThumb`` : `numpy.ndarray`
1365  Thumbnail image of the exposure after flat-field correction.
1366 
1367  Raises
1368  ------
1369  RuntimeError
1370  Raised if a configuration option is set to True, but the
1371  required calibration data has not been specified.
1372 
1373  Notes
1374  -----
1375  The current processed exposure can be viewed by setting the
1376  appropriate lsstDebug entries in the `debug.display`
1377  dictionary. The names of these entries correspond to some of
1378  the IsrTaskConfig Boolean options, with the value denoting the
1379  frame to use. The exposure is shown inside the matching
1380  option check and after the processing of that step has
1381  finished. The steps with debug points are:
1382 
1383  doAssembleCcd
1384  doBias
1385  doCrosstalk
1386  doBrighterFatter
1387  doDark
1388  doFringe
1389  doStrayLight
1390  doFlat
1391 
1392  In addition, setting the "postISRCCD" entry displays the
1393  exposure after all ISR processing has finished.
1394 
1395  """
1396 
1397  if isGen3 is True:
1398  # Gen3 currently cannot automatically do configuration overrides.
1399  # DM-15257 looks to discuss this issue.
1400  # Configure input exposures;
1401 
1402  ccdExposure = self.ensureExposureensureExposure(ccdExposure, camera, detectorNum)
1403  bias = self.ensureExposureensureExposure(bias, camera, detectorNum)
1404  dark = self.ensureExposureensureExposure(dark, camera, detectorNum)
1405  flat = self.ensureExposureensureExposure(flat, camera, detectorNum)
1406  else:
1407  if isinstance(ccdExposure, ButlerDataRef):
1408  return self.runDataRefrunDataRef(ccdExposure)
1409 
1410  ccd = ccdExposure.getDetector()
1411  filterLabel = ccdExposure.getFilterLabel()
1412  physicalFilter = isrFunctions.getPhysicalFilter(filterLabel, self.log)
1413 
1414  if not ccd:
1415  assert not self.config.doAssembleCcd, "You need a Detector to run assembleCcd."
1416  ccd = [FakeAmp(ccdExposure, self.config)]
1417 
1418  # Validate Input
1419  if self.config.doBias and bias is None:
1420  raise RuntimeError("Must supply a bias exposure if config.doBias=True.")
1421  if self.doLinearizedoLinearize(ccd) and linearizer is None:
1422  raise RuntimeError("Must supply a linearizer if config.doLinearize=True for this detector.")
1423  if self.config.doBrighterFatter and bfKernel is None:
1424  raise RuntimeError("Must supply a kernel if config.doBrighterFatter=True.")
1425  if self.config.doDark and dark is None:
1426  raise RuntimeError("Must supply a dark exposure if config.doDark=True.")
1427  if self.config.doFlat and flat is None:
1428  raise RuntimeError("Must supply a flat exposure if config.doFlat=True.")
1429  if self.config.doDefect and defects is None:
1430  raise RuntimeError("Must supply defects if config.doDefect=True.")
1431  if (self.config.doFringe and physicalFilter in self.fringe.config.filters
1432  and fringes.fringes is None):
1433  # The `fringes` object needs to be a pipeBase.Struct, as
1434  # we use it as a `dict` for the parameters of
1435  # `FringeTask.run()`. The `fringes.fringes` `list` may
1436  # not be `None` if `doFringe=True`. Otherwise, raise.
1437  raise RuntimeError("Must supply fringe exposure as a pipeBase.Struct.")
1438  if (self.config.doIlluminationCorrection and physicalFilter in self.config.illumFilters
1439  and illumMaskedImage is None):
1440  raise RuntimeError("Must supply an illumcor if config.doIlluminationCorrection=True.")
1441 
1442  # Begin ISR processing.
1443  if self.config.doConvertIntToFloat:
1444  self.log.info("Converting exposure to floating point values.")
1445  ccdExposure = self.convertIntToFloatconvertIntToFloat(ccdExposure)
1446 
1447  if self.config.doBias and self.config.doBiasBeforeOverscan:
1448  self.log.info("Applying bias correction.")
1449  isrFunctions.biasCorrection(ccdExposure.getMaskedImage(), bias.getMaskedImage(),
1450  trimToFit=self.config.doTrimToMatchCalib)
1451  self.debugViewdebugView(ccdExposure, "doBias")
1452 
1453  # Amplifier level processing.
1454  overscans = []
1455  for amp in ccd:
1456  # if ccdExposure is one amp,
1457  # check for coverage to prevent performing ops multiple times
1458  if ccdExposure.getBBox().contains(amp.getBBox()):
1459  # Check for fully masked bad amplifiers,
1460  # and generate masks for SUSPECT and SATURATED values.
1461  badAmp = self.maskAmplifiermaskAmplifier(ccdExposure, amp, defects)
1462 
1463  if self.config.doOverscan and not badAmp:
1464  # Overscan correction on amp-by-amp basis.
1465  overscanResults = self.overscanCorrectionoverscanCorrection(ccdExposure, amp)
1466  self.log.debug("Corrected overscan for amplifier %s.", amp.getName())
1467  if overscanResults is not None and \
1468  self.config.qa is not None and self.config.qa.saveStats is True:
1469  if isinstance(overscanResults.overscanFit, float):
1470  qaMedian = overscanResults.overscanFit
1471  qaStdev = float("NaN")
1472  else:
1473  qaStats = afwMath.makeStatistics(overscanResults.overscanFit,
1474  afwMath.MEDIAN | afwMath.STDEVCLIP)
1475  qaMedian = qaStats.getValue(afwMath.MEDIAN)
1476  qaStdev = qaStats.getValue(afwMath.STDEVCLIP)
1477 
1478  self.metadata.set(f"FIT MEDIAN {amp.getName()}", qaMedian)
1479  self.metadata.set(f"FIT STDEV {amp.getName()}", qaStdev)
1480  self.log.debug(" Overscan stats for amplifer %s: %f +/- %f",
1481  amp.getName(), qaMedian, qaStdev)
1482 
1483  # Residuals after overscan correction
1484  qaStatsAfter = afwMath.makeStatistics(overscanResults.overscanImage,
1485  afwMath.MEDIAN | afwMath.STDEVCLIP)
1486  qaMedianAfter = qaStatsAfter.getValue(afwMath.MEDIAN)
1487  qaStdevAfter = qaStatsAfter.getValue(afwMath.STDEVCLIP)
1488 
1489  self.metadata.set(f"RESIDUAL MEDIAN {amp.getName()}", qaMedianAfter)
1490  self.metadata.set(f"RESIDUAL STDEV {amp.getName()}", qaStdevAfter)
1491  self.log.debug(" Overscan stats for amplifer %s after correction: %f +/- %f",
1492  amp.getName(), qaMedianAfter, qaStdevAfter)
1493 
1494  ccdExposure.getMetadata().set('OVERSCAN', "Overscan corrected")
1495  else:
1496  if badAmp:
1497  self.log.warning("Amplifier %s is bad.", amp.getName())
1498  overscanResults = None
1499 
1500  overscans.append(overscanResults if overscanResults is not None else None)
1501  else:
1502  self.log.info("Skipped OSCAN for %s.", amp.getName())
1503 
1504  if self.config.doCrosstalk and self.config.doCrosstalkBeforeAssemble:
1505  self.log.info("Applying crosstalk correction.")
1506  self.crosstalk.run(ccdExposure, crosstalk=crosstalk,
1507  crosstalkSources=crosstalkSources, camera=camera)
1508  self.debugViewdebugView(ccdExposure, "doCrosstalk")
1509 
1510  if self.config.doAssembleCcd:
1511  self.log.info("Assembling CCD from amplifiers.")
1512  ccdExposure = self.assembleCcd.assembleCcd(ccdExposure)
1513 
1514  if self.config.expectWcs and not ccdExposure.getWcs():
1515  self.log.warning("No WCS found in input exposure.")
1516  self.debugViewdebugView(ccdExposure, "doAssembleCcd")
1517 
1518  ossThumb = None
1519  if self.config.qa.doThumbnailOss:
1520  ossThumb = isrQa.makeThumbnail(ccdExposure, isrQaConfig=self.config.qa)
1521 
1522  if self.config.doBias and not self.config.doBiasBeforeOverscan:
1523  self.log.info("Applying bias correction.")
1524  isrFunctions.biasCorrection(ccdExposure.getMaskedImage(), bias.getMaskedImage(),
1525  trimToFit=self.config.doTrimToMatchCalib)
1526  self.debugViewdebugView(ccdExposure, "doBias")
1527 
1528  if self.config.doVariance:
1529  for amp, overscanResults in zip(ccd, overscans):
1530  if ccdExposure.getBBox().contains(amp.getBBox()):
1531  self.log.debug("Constructing variance map for amplifer %s.", amp.getName())
1532  ampExposure = ccdExposure.Factory(ccdExposure, amp.getBBox())
1533  if overscanResults is not None:
1534  self.updateVarianceupdateVariance(ampExposure, amp,
1535  overscanImage=overscanResults.overscanImage,
1536  ptcDataset=ptc)
1537  else:
1538  self.updateVarianceupdateVariance(ampExposure, amp,
1539  overscanImage=None,
1540  ptcDataset=ptc)
1541  if self.config.qa is not None and self.config.qa.saveStats is True:
1542  qaStats = afwMath.makeStatistics(ampExposure.getVariance(),
1543  afwMath.MEDIAN | afwMath.STDEVCLIP)
1544  self.metadata.set(f"ISR VARIANCE {amp.getName()} MEDIAN",
1545  qaStats.getValue(afwMath.MEDIAN))
1546  self.metadata.set(f"ISR VARIANCE {amp.getName()} STDEV",
1547  qaStats.getValue(afwMath.STDEVCLIP))
1548  self.log.debug(" Variance stats for amplifer %s: %f +/- %f.",
1549  amp.getName(), qaStats.getValue(afwMath.MEDIAN),
1550  qaStats.getValue(afwMath.STDEVCLIP))
1551  if self.config.maskNegativeVariance:
1552  self.maskNegativeVariancemaskNegativeVariance(ccdExposure)
1553 
1554  if self.doLinearizedoLinearize(ccd):
1555  self.log.info("Applying linearizer.")
1556  linearizer.applyLinearity(image=ccdExposure.getMaskedImage().getImage(),
1557  detector=ccd, log=self.log)
1558 
1559  if self.config.doCrosstalk and not self.config.doCrosstalkBeforeAssemble:
1560  self.log.info("Applying crosstalk correction.")
1561  self.crosstalk.run(ccdExposure, crosstalk=crosstalk,
1562  crosstalkSources=crosstalkSources, isTrimmed=True)
1563  self.debugViewdebugView(ccdExposure, "doCrosstalk")
1564 
1565  # Masking block. Optionally mask known defects, NAN/inf pixels,
1566  # widen trails, and do anything else the camera needs. Saturated and
1567  # suspect pixels have already been masked.
1568  if self.config.doDefect:
1569  self.log.info("Masking defects.")
1570  self.maskDefectmaskDefect(ccdExposure, defects)
1571 
1572  if self.config.numEdgeSuspect > 0:
1573  self.log.info("Masking edges as SUSPECT.")
1574  self.maskEdgesmaskEdges(ccdExposure, numEdgePixels=self.config.numEdgeSuspect,
1575  maskPlane="SUSPECT", level=self.config.edgeMaskLevel)
1576 
1577  if self.config.doNanMasking:
1578  self.log.info("Masking non-finite (NAN, inf) value pixels.")
1579  self.maskNanmaskNan(ccdExposure)
1580 
1581  if self.config.doWidenSaturationTrails:
1582  self.log.info("Widening saturation trails.")
1583  isrFunctions.widenSaturationTrails(ccdExposure.getMaskedImage().getMask())
1584 
1585  if self.config.doCameraSpecificMasking:
1586  self.log.info("Masking regions for camera specific reasons.")
1587  self.masking.run(ccdExposure)
1588 
1589  if self.config.doBrighterFatter:
1590  # We need to apply flats and darks before we can interpolate, and
1591  # we need to interpolate before we do B-F, but we do B-F without
1592  # the flats and darks applied so we can work in units of electrons
1593  # or holes. This context manager applies and then removes the darks
1594  # and flats.
1595  #
1596  # We also do not want to interpolate values here, so operate on
1597  # temporary images so we can apply only the BF-correction and roll
1598  # back the interpolation.
1599  interpExp = ccdExposure.clone()
1600  with self.flatContextflatContext(interpExp, flat, dark):
1601  isrFunctions.interpolateFromMask(
1602  maskedImage=interpExp.getMaskedImage(),
1603  fwhm=self.config.fwhm,
1604  growSaturatedFootprints=self.config.growSaturationFootprintSize,
1605  maskNameList=list(self.config.brighterFatterMaskListToInterpolate)
1606  )
1607  bfExp = interpExp.clone()
1608 
1609  self.log.info("Applying brighter-fatter correction using kernel type %s / gains %s.",
1610  type(bfKernel), type(bfGains))
1611  bfResults = isrFunctions.brighterFatterCorrection(bfExp, bfKernel,
1612  self.config.brighterFatterMaxIter,
1613  self.config.brighterFatterThreshold,
1614  self.config.brighterFatterApplyGain,
1615  bfGains)
1616  if bfResults[1] == self.config.brighterFatterMaxIter:
1617  self.log.warning("Brighter-fatter correction did not converge, final difference %f.",
1618  bfResults[0])
1619  else:
1620  self.log.info("Finished brighter-fatter correction in %d iterations.",
1621  bfResults[1])
1622  image = ccdExposure.getMaskedImage().getImage()
1623  bfCorr = bfExp.getMaskedImage().getImage()
1624  bfCorr -= interpExp.getMaskedImage().getImage()
1625  image += bfCorr
1626 
1627  # Applying the brighter-fatter correction applies a
1628  # convolution to the science image. At the edges this
1629  # convolution may not have sufficient valid pixels to
1630  # produce a valid correction. Mark pixels within the size
1631  # of the brighter-fatter kernel as EDGE to warn of this
1632  # fact.
1633  self.log.info("Ensuring image edges are masked as EDGE to the brighter-fatter kernel size.")
1634  self.maskEdgesmaskEdges(ccdExposure, numEdgePixels=numpy.max(bfKernel.shape) // 2,
1635  maskPlane="EDGE")
1636 
1637  if self.config.brighterFatterMaskGrowSize > 0:
1638  self.log.info("Growing masks to account for brighter-fatter kernel convolution.")
1639  for maskPlane in self.config.brighterFatterMaskListToInterpolate:
1640  isrFunctions.growMasks(ccdExposure.getMask(),
1641  radius=self.config.brighterFatterMaskGrowSize,
1642  maskNameList=maskPlane,
1643  maskValue=maskPlane)
1644 
1645  self.debugViewdebugView(ccdExposure, "doBrighterFatter")
1646 
1647  if self.config.doDark:
1648  self.log.info("Applying dark correction.")
1649  self.darkCorrectiondarkCorrection(ccdExposure, dark)
1650  self.debugViewdebugView(ccdExposure, "doDark")
1651 
1652  if self.config.doFringe and not self.config.fringeAfterFlat:
1653  self.log.info("Applying fringe correction before flat.")
1654  self.fringe.run(ccdExposure, **fringes.getDict())
1655  self.debugViewdebugView(ccdExposure, "doFringe")
1656 
1657  if self.config.doStrayLight and self.strayLight.check(ccdExposure):
1658  self.log.info("Checking strayLight correction.")
1659  self.strayLight.run(ccdExposure, strayLightData)
1660  self.debugViewdebugView(ccdExposure, "doStrayLight")
1661 
1662  if self.config.doFlat:
1663  self.log.info("Applying flat correction.")
1664  self.flatCorrectionflatCorrection(ccdExposure, flat)
1665  self.debugViewdebugView(ccdExposure, "doFlat")
1666 
1667  if self.config.doApplyGains:
1668  self.log.info("Applying gain correction instead of flat.")
1669  if self.config.usePtcGains:
1670  self.log.info("Using gains from the Photon Transfer Curve.")
1671  isrFunctions.applyGains(ccdExposure, self.config.normalizeGains,
1672  ptcGains=ptc.gain)
1673  else:
1674  isrFunctions.applyGains(ccdExposure, self.config.normalizeGains)
1675 
1676  if self.config.doFringe and self.config.fringeAfterFlat:
1677  self.log.info("Applying fringe correction after flat.")
1678  self.fringe.run(ccdExposure, **fringes.getDict())
1679 
1680  if self.config.doVignette:
1681  self.log.info("Constructing Vignette polygon.")
1682  self.vignettePolygonvignettePolygon = self.vignette.run(ccdExposure)
1683 
1684  if self.config.vignette.doWriteVignettePolygon:
1685  self.setValidPolygonIntersectsetValidPolygonIntersect(ccdExposure, self.vignettePolygonvignettePolygon)
1686 
1687  if self.config.doAttachTransmissionCurve:
1688  self.log.info("Adding transmission curves.")
1689  isrFunctions.attachTransmissionCurve(ccdExposure, opticsTransmission=opticsTransmission,
1690  filterTransmission=filterTransmission,
1691  sensorTransmission=sensorTransmission,
1692  atmosphereTransmission=atmosphereTransmission)
1693 
1694  flattenedThumb = None
1695  if self.config.qa.doThumbnailFlattened:
1696  flattenedThumb = isrQa.makeThumbnail(ccdExposure, isrQaConfig=self.config.qa)
1697 
1698  if self.config.doIlluminationCorrection and physicalFilter in self.config.illumFilters:
1699  self.log.info("Performing illumination correction.")
1700  isrFunctions.illuminationCorrection(ccdExposure.getMaskedImage(),
1701  illumMaskedImage, illumScale=self.config.illumScale,
1702  trimToFit=self.config.doTrimToMatchCalib)
1703 
1704  preInterpExp = None
1705  if self.config.doSaveInterpPixels:
1706  preInterpExp = ccdExposure.clone()
1707 
1708  # Reset and interpolate bad pixels.
1709  #
1710  # Large contiguous bad regions (which should have the BAD mask
1711  # bit set) should have their values set to the image median.
1712  # This group should include defects and bad amplifiers. As the
1713  # area covered by these defects are large, there's little
1714  # reason to expect that interpolation would provide a more
1715  # useful value.
1716  #
1717  # Smaller defects can be safely interpolated after the larger
1718  # regions have had their pixel values reset. This ensures
1719  # that the remaining defects adjacent to bad amplifiers (as an
1720  # example) do not attempt to interpolate extreme values.
1721  if self.config.doSetBadRegions:
1722  badPixelCount, badPixelValue = isrFunctions.setBadRegions(ccdExposure)
1723  if badPixelCount > 0:
1724  self.log.info("Set %d BAD pixels to %f.", badPixelCount, badPixelValue)
1725 
1726  if self.config.doInterpolate:
1727  self.log.info("Interpolating masked pixels.")
1728  isrFunctions.interpolateFromMask(
1729  maskedImage=ccdExposure.getMaskedImage(),
1730  fwhm=self.config.fwhm,
1731  growSaturatedFootprints=self.config.growSaturationFootprintSize,
1732  maskNameList=list(self.config.maskListToInterpolate)
1733  )
1734 
1735  self.roughZeroPointroughZeroPoint(ccdExposure)
1736 
1737  # correct for amp offsets within the CCD
1738  if self.config.doAmpOffset:
1739  self.log.info("Correcting amp offsets.")
1740  self.ampOffset.run(ccdExposure)
1741 
1742  if self.config.doMeasureBackground:
1743  self.log.info("Measuring background level.")
1744  self.measureBackgroundmeasureBackground(ccdExposure, self.config.qa)
1745 
1746  if self.config.qa is not None and self.config.qa.saveStats is True:
1747  for amp in ccd:
1748  ampExposure = ccdExposure.Factory(ccdExposure, amp.getBBox())
1749  qaStats = afwMath.makeStatistics(ampExposure.getImage(),
1750  afwMath.MEDIAN | afwMath.STDEVCLIP)
1751  self.metadata.set("ISR BACKGROUND {} MEDIAN".format(amp.getName()),
1752  qaStats.getValue(afwMath.MEDIAN))
1753  self.metadata.set("ISR BACKGROUND {} STDEV".format(amp.getName()),
1754  qaStats.getValue(afwMath.STDEVCLIP))
1755  self.log.debug(" Background stats for amplifer %s: %f +/- %f",
1756  amp.getName(), qaStats.getValue(afwMath.MEDIAN),
1757  qaStats.getValue(afwMath.STDEVCLIP))
1758 
1759  self.debugViewdebugView(ccdExposure, "postISRCCD")
1760 
1761  return pipeBase.Struct(
1762  exposure=ccdExposure,
1763  ossThumb=ossThumb,
1764  flattenedThumb=flattenedThumb,
1765 
1766  preInterpExposure=preInterpExp,
1767  outputExposure=ccdExposure,
1768  outputOssThumbnail=ossThumb,
1769  outputFlattenedThumbnail=flattenedThumb,
1770  )
1771 
1772  @pipeBase.timeMethod
1773  def runDataRef(self, sensorRef):
1774  """Perform instrument signature removal on a ButlerDataRef of a Sensor.
1775 
1776  This method contains the `CmdLineTask` interface to the ISR
1777  processing. All IO is handled here, freeing the `run()` method
1778  to manage only pixel-level calculations. The steps performed
1779  are:
1780  - Read in necessary detrending/isr/calibration data.
1781  - Process raw exposure in `run()`.
1782  - Persist the ISR-corrected exposure as "postISRCCD" if
1783  config.doWrite=True.
1784 
1785  Parameters
1786  ----------
1787  sensorRef : `daf.persistence.butlerSubset.ButlerDataRef`
1788  DataRef of the detector data to be processed
1789 
1790  Returns
1791  -------
1792  result : `lsst.pipe.base.Struct`
1793  Result struct with component:
1794  - ``exposure`` : `afw.image.Exposure`
1795  The fully ISR corrected exposure.
1796 
1797  Raises
1798  ------
1799  RuntimeError
1800  Raised if a configuration option is set to True, but the
1801  required calibration data does not exist.
1802 
1803  """
1804  self.log.info("Performing ISR on sensor %s.", sensorRef.dataId)
1805 
1806  ccdExposure = sensorRef.get(self.config.datasetType)
1807 
1808  camera = sensorRef.get("camera")
1809  isrData = self.readIsrDatareadIsrData(sensorRef, ccdExposure)
1810 
1811  result = self.runrun(ccdExposure, camera=camera, **isrData.getDict())
1812 
1813  if self.config.doWrite:
1814  sensorRef.put(result.exposure, "postISRCCD")
1815  if result.preInterpExposure is not None:
1816  sensorRef.put(result.preInterpExposure, "postISRCCD_uninterpolated")
1817  if result.ossThumb is not None:
1818  isrQa.writeThumbnail(sensorRef, result.ossThumb, "ossThumb")
1819  if result.flattenedThumb is not None:
1820  isrQa.writeThumbnail(sensorRef, result.flattenedThumb, "flattenedThumb")
1821 
1822  return result
1823 
1824  def getIsrExposure(self, dataRef, datasetType, dateObs=None, immediate=True):
1825  """Retrieve a calibration dataset for removing instrument signature.
1826 
1827  Parameters
1828  ----------
1829 
1830  dataRef : `daf.persistence.butlerSubset.ButlerDataRef`
1831  DataRef of the detector data to find calibration datasets
1832  for.
1833  datasetType : `str`
1834  Type of dataset to retrieve (e.g. 'bias', 'flat', etc).
1835  dateObs : `str`, optional
1836  Date of the observation. Used to correct butler failures
1837  when using fallback filters.
1838  immediate : `Bool`
1839  If True, disable butler proxies to enable error handling
1840  within this routine.
1841 
1842  Returns
1843  -------
1844  exposure : `lsst.afw.image.Exposure`
1845  Requested calibration frame.
1846 
1847  Raises
1848  ------
1849  RuntimeError
1850  Raised if no matching calibration frame can be found.
1851  """
1852  try:
1853  exp = dataRef.get(datasetType, immediate=immediate)
1854  except Exception as exc1:
1855  if not self.config.fallbackFilterName:
1856  raise RuntimeError("Unable to retrieve %s for %s: %s." % (datasetType, dataRef.dataId, exc1))
1857  try:
1858  if self.config.useFallbackDate and dateObs:
1859  exp = dataRef.get(datasetType, filter=self.config.fallbackFilterName,
1860  dateObs=dateObs, immediate=immediate)
1861  else:
1862  exp = dataRef.get(datasetType, filter=self.config.fallbackFilterName, immediate=immediate)
1863  except Exception as exc2:
1864  raise RuntimeError("Unable to retrieve %s for %s, even with fallback filter %s: %s AND %s." %
1865  (datasetType, dataRef.dataId, self.config.fallbackFilterName, exc1, exc2))
1866  self.log.warning("Using fallback calibration from filter %s.", self.config.fallbackFilterName)
1867 
1868  if self.config.doAssembleIsrExposures:
1869  exp = self.assembleCcd.assembleCcd(exp)
1870  return exp
1871 
1872  def ensureExposure(self, inputExp, camera=None, detectorNum=None):
1873  """Ensure that the data returned by Butler is a fully constructed exp.
1874 
1875  ISR requires exposure-level image data for historical reasons, so if we
1876  did not recieve that from Butler, construct it from what we have,
1877  modifying the input in place.
1878 
1879  Parameters
1880  ----------
1881  inputExp : `lsst.afw.image.Exposure`, `lsst.afw.image.DecoratedImageU`,
1882  or `lsst.afw.image.ImageF`
1883  The input data structure obtained from Butler.
1884  camera : `lsst.afw.cameraGeom.camera`, optional
1885  The camera associated with the image. Used to find the appropriate
1886  detector if detector is not already set.
1887  detectorNum : `int`, optional
1888  The detector in the camera to attach, if the detector is not
1889  already set.
1890 
1891  Returns
1892  -------
1893  inputExp : `lsst.afw.image.Exposure`
1894  The re-constructed exposure, with appropriate detector parameters.
1895 
1896  Raises
1897  ------
1898  TypeError
1899  Raised if the input data cannot be used to construct an exposure.
1900  """
1901  if isinstance(inputExp, afwImage.DecoratedImageU):
1902  inputExp = afwImage.makeExposure(afwImage.makeMaskedImage(inputExp))
1903  elif isinstance(inputExp, afwImage.ImageF):
1904  inputExp = afwImage.makeExposure(afwImage.makeMaskedImage(inputExp))
1905  elif isinstance(inputExp, afwImage.MaskedImageF):
1906  inputExp = afwImage.makeExposure(inputExp)
1907  elif isinstance(inputExp, afwImage.Exposure):
1908  pass
1909  elif inputExp is None:
1910  # Assume this will be caught by the setup if it is a problem.
1911  return inputExp
1912  else:
1913  raise TypeError("Input Exposure is not known type in isrTask.ensureExposure: %s." %
1914  (type(inputExp), ))
1915 
1916  if inputExp.getDetector() is None:
1917  if camera is None or detectorNum is None:
1918  raise RuntimeError('Must supply both a camera and detector number when using exposures '
1919  'without a detector set.')
1920  inputExp.setDetector(camera[detectorNum])
1921 
1922  return inputExp
1923 
1924  def convertIntToFloat(self, exposure):
1925  """Convert exposure image from uint16 to float.
1926 
1927  If the exposure does not need to be converted, the input is
1928  immediately returned. For exposures that are converted to use
1929  floating point pixels, the variance is set to unity and the
1930  mask to zero.
1931 
1932  Parameters
1933  ----------
1934  exposure : `lsst.afw.image.Exposure`
1935  The raw exposure to be converted.
1936 
1937  Returns
1938  -------
1939  newexposure : `lsst.afw.image.Exposure`
1940  The input ``exposure``, converted to floating point pixels.
1941 
1942  Raises
1943  ------
1944  RuntimeError
1945  Raised if the exposure type cannot be converted to float.
1946 
1947  """
1948  if isinstance(exposure, afwImage.ExposureF):
1949  # Nothing to be done
1950  self.log.debug("Exposure already of type float.")
1951  return exposure
1952  if not hasattr(exposure, "convertF"):
1953  raise RuntimeError("Unable to convert exposure (%s) to float." % type(exposure))
1954 
1955  newexposure = exposure.convertF()
1956  newexposure.variance[:] = 1
1957  newexposure.mask[:] = 0x0
1958 
1959  return newexposure
1960 
1961  def maskAmplifier(self, ccdExposure, amp, defects):
1962  """Identify bad amplifiers, saturated and suspect pixels.
1963 
1964  Parameters
1965  ----------
1966  ccdExposure : `lsst.afw.image.Exposure`
1967  Input exposure to be masked.
1968  amp : `lsst.afw.table.AmpInfoCatalog`
1969  Catalog of parameters defining the amplifier on this
1970  exposure to mask.
1971  defects : `lsst.ip.isr.Defects`
1972  List of defects. Used to determine if the entire
1973  amplifier is bad.
1974 
1975  Returns
1976  -------
1977  badAmp : `Bool`
1978  If this is true, the entire amplifier area is covered by
1979  defects and unusable.
1980 
1981  """
1982  maskedImage = ccdExposure.getMaskedImage()
1983 
1984  badAmp = False
1985 
1986  # Check if entire amp region is defined as a defect
1987  # NB: need to use amp.getBBox() for correct comparison with current
1988  # defects definition.
1989  if defects is not None:
1990  badAmp = bool(sum([v.getBBox().contains(amp.getBBox()) for v in defects]))
1991 
1992  # In the case of a bad amp, we will set mask to "BAD"
1993  # (here use amp.getRawBBox() for correct association with pixels in
1994  # current ccdExposure).
1995  if badAmp:
1996  dataView = afwImage.MaskedImageF(maskedImage, amp.getRawBBox(),
1997  afwImage.PARENT)
1998  maskView = dataView.getMask()
1999  maskView |= maskView.getPlaneBitMask("BAD")
2000  del maskView
2001  return badAmp
2002 
2003  # Mask remaining defects after assembleCcd() to allow for defects that
2004  # cross amplifier boundaries. Saturation and suspect pixels can be
2005  # masked now, though.
2006  limits = dict()
2007  if self.config.doSaturation and not badAmp:
2008  limits.update({self.config.saturatedMaskName: amp.getSaturation()})
2009  if self.config.doSuspect and not badAmp:
2010  limits.update({self.config.suspectMaskName: amp.getSuspectLevel()})
2011  if math.isfinite(self.config.saturation):
2012  limits.update({self.config.saturatedMaskName: self.config.saturation})
2013 
2014  for maskName, maskThreshold in limits.items():
2015  if not math.isnan(maskThreshold):
2016  dataView = maskedImage.Factory(maskedImage, amp.getRawBBox())
2017  isrFunctions.makeThresholdMask(
2018  maskedImage=dataView,
2019  threshold=maskThreshold,
2020  growFootprints=0,
2021  maskName=maskName
2022  )
2023 
2024  # Determine if we've fully masked this amplifier with SUSPECT and
2025  # SAT pixels.
2026  maskView = afwImage.Mask(maskedImage.getMask(), amp.getRawDataBBox(),
2027  afwImage.PARENT)
2028  maskVal = maskView.getPlaneBitMask([self.config.saturatedMaskName,
2029  self.config.suspectMaskName])
2030  if numpy.all(maskView.getArray() & maskVal > 0):
2031  badAmp = True
2032  maskView |= maskView.getPlaneBitMask("BAD")
2033 
2034  return badAmp
2035 
2036  def overscanCorrection(self, ccdExposure, amp):
2037  """Apply overscan correction in place.
2038 
2039  This method does initial pixel rejection of the overscan
2040  region. The overscan can also be optionally segmented to
2041  allow for discontinuous overscan responses to be fit
2042  separately. The actual overscan subtraction is performed by
2043  the `lsst.ip.isr.isrFunctions.overscanCorrection` function,
2044  which is called here after the amplifier is preprocessed.
2045 
2046  Parameters
2047  ----------
2048  ccdExposure : `lsst.afw.image.Exposure`
2049  Exposure to have overscan correction performed.
2050  amp : `lsst.afw.cameraGeom.Amplifer`
2051  The amplifier to consider while correcting the overscan.
2052 
2053  Returns
2054  -------
2055  overscanResults : `lsst.pipe.base.Struct`
2056  Result struct with components:
2057  - ``imageFit`` : scalar or `lsst.afw.image.Image`
2058  Value or fit subtracted from the amplifier image data.
2059  - ``overscanFit`` : scalar or `lsst.afw.image.Image`
2060  Value or fit subtracted from the overscan image data.
2061  - ``overscanImage`` : `lsst.afw.image.Image`
2062  Image of the overscan region with the overscan
2063  correction applied. This quantity is used to estimate
2064  the amplifier read noise empirically.
2065 
2066  Raises
2067  ------
2068  RuntimeError
2069  Raised if the ``amp`` does not contain raw pixel information.
2070 
2071  See Also
2072  --------
2073  lsst.ip.isr.isrFunctions.overscanCorrection
2074  """
2075  if amp.getRawHorizontalOverscanBBox().isEmpty():
2076  self.log.info("ISR_OSCAN: No overscan region. Not performing overscan correction.")
2077  return None
2078 
2079  statControl = afwMath.StatisticsControl()
2080  statControl.setAndMask(ccdExposure.mask.getPlaneBitMask("SAT"))
2081 
2082  # Determine the bounding boxes
2083  dataBBox = amp.getRawDataBBox()
2084  oscanBBox = amp.getRawHorizontalOverscanBBox()
2085  dx0 = 0
2086  dx1 = 0
2087 
2088  prescanBBox = amp.getRawPrescanBBox()
2089  if (oscanBBox.getBeginX() > prescanBBox.getBeginX()): # amp is at the right
2090  dx0 += self.config.overscanNumLeadingColumnsToSkip
2091  dx1 -= self.config.overscanNumTrailingColumnsToSkip
2092  else:
2093  dx0 += self.config.overscanNumTrailingColumnsToSkip
2094  dx1 -= self.config.overscanNumLeadingColumnsToSkip
2095 
2096  # Determine if we need to work on subregions of the amplifier
2097  # and overscan.
2098  imageBBoxes = []
2099  overscanBBoxes = []
2100 
2101  if ((self.config.overscanBiasJump
2102  and self.config.overscanBiasJumpLocation)
2103  and (ccdExposure.getMetadata().exists(self.config.overscanBiasJumpKeyword)
2104  and ccdExposure.getMetadata().getScalar(self.config.overscanBiasJumpKeyword) in
2105  self.config.overscanBiasJumpDevices)):
2106  if amp.getReadoutCorner() in (ReadoutCorner.LL, ReadoutCorner.LR):
2107  yLower = self.config.overscanBiasJumpLocation
2108  yUpper = dataBBox.getHeight() - yLower
2109  else:
2110  yUpper = self.config.overscanBiasJumpLocation
2111  yLower = dataBBox.getHeight() - yUpper
2112 
2113  imageBBoxes.append(lsst.geom.Box2I(dataBBox.getBegin(),
2114  lsst.geom.Extent2I(dataBBox.getWidth(), yLower)))
2115  overscanBBoxes.append(lsst.geom.Box2I(oscanBBox.getBegin() + lsst.geom.Extent2I(dx0, 0),
2116  lsst.geom.Extent2I(oscanBBox.getWidth() - dx0 + dx1,
2117  yLower)))
2118 
2119  imageBBoxes.append(lsst.geom.Box2I(dataBBox.getBegin() + lsst.geom.Extent2I(0, yLower),
2120  lsst.geom.Extent2I(dataBBox.getWidth(), yUpper)))
2121  overscanBBoxes.append(lsst.geom.Box2I(oscanBBox.getBegin() + lsst.geom.Extent2I(dx0, yLower),
2122  lsst.geom.Extent2I(oscanBBox.getWidth() - dx0 + dx1,
2123  yUpper)))
2124  else:
2125  imageBBoxes.append(lsst.geom.Box2I(dataBBox.getBegin(),
2126  lsst.geom.Extent2I(dataBBox.getWidth(), dataBBox.getHeight())))
2127  overscanBBoxes.append(lsst.geom.Box2I(oscanBBox.getBegin() + lsst.geom.Extent2I(dx0, 0),
2128  lsst.geom.Extent2I(oscanBBox.getWidth() - dx0 + dx1,
2129  oscanBBox.getHeight())))
2130 
2131  # Perform overscan correction on subregions, ensuring saturated
2132  # pixels are masked.
2133  for imageBBox, overscanBBox in zip(imageBBoxes, overscanBBoxes):
2134  ampImage = ccdExposure.maskedImage[imageBBox]
2135  overscanImage = ccdExposure.maskedImage[overscanBBox]
2136 
2137  overscanArray = overscanImage.image.array
2138  median = numpy.ma.median(numpy.ma.masked_where(overscanImage.mask.array, overscanArray))
2139  bad = numpy.where(numpy.abs(overscanArray - median) > self.config.overscanMaxDev)
2140  overscanImage.mask.array[bad] = overscanImage.mask.getPlaneBitMask("SAT")
2141 
2142  statControl = afwMath.StatisticsControl()
2143  statControl.setAndMask(ccdExposure.mask.getPlaneBitMask("SAT"))
2144 
2145  overscanResults = self.overscan.run(ampImage.getImage(), overscanImage, amp)
2146 
2147  # Measure average overscan levels and record them in the metadata.
2148  levelStat = afwMath.MEDIAN
2149  sigmaStat = afwMath.STDEVCLIP
2150 
2151  sctrl = afwMath.StatisticsControl(self.config.qa.flatness.clipSigma,
2152  self.config.qa.flatness.nIter)
2153  metadata = ccdExposure.getMetadata()
2154  ampNum = amp.getName()
2155  # if self.config.overscanFitType in ("MEDIAN", "MEAN", "MEANCLIP"):
2156  if isinstance(overscanResults.overscanFit, float):
2157  metadata.set("ISR_OSCAN_LEVEL%s" % ampNum, overscanResults.overscanFit)
2158  metadata.set("ISR_OSCAN_SIGMA%s" % ampNum, 0.0)
2159  else:
2160  stats = afwMath.makeStatistics(overscanResults.overscanFit, levelStat | sigmaStat, sctrl)
2161  metadata.set("ISR_OSCAN_LEVEL%s" % ampNum, stats.getValue(levelStat))
2162  metadata.set("ISR_OSCAN_SIGMA%s" % ampNum, stats.getValue(sigmaStat))
2163 
2164  return overscanResults
2165 
2166  def updateVariance(self, ampExposure, amp, overscanImage=None, ptcDataset=None):
2167  """Set the variance plane using the gain and read noise
2168 
2169  The read noise is calculated from the ``overscanImage`` if the
2170  ``doEmpiricalReadNoise`` option is set in the configuration; otherwise
2171  the value from the amplifier data is used.
2172 
2173  Parameters
2174  ----------
2175  ampExposure : `lsst.afw.image.Exposure`
2176  Exposure to process.
2177  amp : `lsst.afw.table.AmpInfoRecord` or `FakeAmp`
2178  Amplifier detector data.
2179  overscanImage : `lsst.afw.image.MaskedImage`, optional.
2180  Image of overscan, required only for empirical read noise.
2181  ptcDataset : `lsst.ip.isr.PhotonTransferCurveDataset`, optional
2182  PTC dataset containing the gains and read noise.
2183 
2184 
2185  Raises
2186  ------
2187  RuntimeError
2188  Raised if either ``usePtcGains`` of ``usePtcReadNoise``
2189  are ``True``, but ptcDataset is not provided.
2190 
2191  Raised if ```doEmpiricalReadNoise`` is ``True`` but
2192  ``overscanImage`` is ``None``.
2193 
2194  See also
2195  --------
2196  lsst.ip.isr.isrFunctions.updateVariance
2197  """
2198  maskPlanes = [self.config.saturatedMaskName, self.config.suspectMaskName]
2199  if self.config.usePtcGains:
2200  if ptcDataset is None:
2201  raise RuntimeError("No ptcDataset provided to use PTC gains.")
2202  else:
2203  gain = ptcDataset.gain[amp.getName()]
2204  self.log.info("Using gain from Photon Transfer Curve.")
2205  else:
2206  gain = amp.getGain()
2207 
2208  if math.isnan(gain):
2209  gain = 1.0
2210  self.log.warning("Gain set to NAN! Updating to 1.0 to generate Poisson variance.")
2211  elif gain <= 0:
2212  patchedGain = 1.0
2213  self.log.warning("Gain for amp %s == %g <= 0; setting to %f.",
2214  amp.getName(), gain, patchedGain)
2215  gain = patchedGain
2216 
2217  if self.config.doEmpiricalReadNoise and overscanImage is None:
2218  raise RuntimeError("Overscan is none for EmpiricalReadNoise.")
2219 
2220  if self.config.doEmpiricalReadNoise and overscanImage is not None:
2221  stats = afwMath.StatisticsControl()
2222  stats.setAndMask(overscanImage.mask.getPlaneBitMask(maskPlanes))
2223  readNoise = afwMath.makeStatistics(overscanImage, afwMath.STDEVCLIP, stats).getValue()
2224  self.log.info("Calculated empirical read noise for amp %s: %f.",
2225  amp.getName(), readNoise)
2226  elif self.config.usePtcReadNoise:
2227  if ptcDataset is None:
2228  raise RuntimeError("No ptcDataset provided to use PTC readnoise.")
2229  else:
2230  readNoise = ptcDataset.noise[amp.getName()]
2231  self.log.info("Using read noise from Photon Transfer Curve.")
2232  else:
2233  readNoise = amp.getReadNoise()
2234 
2235  isrFunctions.updateVariance(
2236  maskedImage=ampExposure.getMaskedImage(),
2237  gain=gain,
2238  readNoise=readNoise,
2239  )
2240 
2241  def maskNegativeVariance(self, exposure):
2242  """Identify and mask pixels with negative variance values.
2243 
2244  Parameters
2245  ----------
2246  exposure : `lsst.afw.image.Exposure`
2247  Exposure to process.
2248 
2249  See Also
2250  --------
2251  lsst.ip.isr.isrFunctions.updateVariance
2252  """
2253  maskPlane = exposure.getMask().getPlaneBitMask(self.config.negativeVarianceMaskName)
2254  bad = numpy.where(exposure.getVariance().getArray() <= 0.0)
2255  exposure.mask.array[bad] |= maskPlane
2256 
2257  def darkCorrection(self, exposure, darkExposure, invert=False):
2258  """Apply dark correction in place.
2259 
2260  Parameters
2261  ----------
2262  exposure : `lsst.afw.image.Exposure`
2263  Exposure to process.
2264  darkExposure : `lsst.afw.image.Exposure`
2265  Dark exposure of the same size as ``exposure``.
2266  invert : `Bool`, optional
2267  If True, re-add the dark to an already corrected image.
2268 
2269  Raises
2270  ------
2271  RuntimeError
2272  Raised if either ``exposure`` or ``darkExposure`` do not
2273  have their dark time defined.
2274 
2275  See Also
2276  --------
2277  lsst.ip.isr.isrFunctions.darkCorrection
2278  """
2279  expScale = exposure.getInfo().getVisitInfo().getDarkTime()
2280  if math.isnan(expScale):
2281  raise RuntimeError("Exposure darktime is NAN.")
2282  if darkExposure.getInfo().getVisitInfo() is not None \
2283  and not math.isnan(darkExposure.getInfo().getVisitInfo().getDarkTime()):
2284  darkScale = darkExposure.getInfo().getVisitInfo().getDarkTime()
2285  else:
2286  # DM-17444: darkExposure.getInfo.getVisitInfo() is None
2287  # so getDarkTime() does not exist.
2288  self.log.warning("darkExposure.getInfo().getVisitInfo() does not exist. Using darkScale = 1.0.")
2289  darkScale = 1.0
2290 
2291  isrFunctions.darkCorrection(
2292  maskedImage=exposure.getMaskedImage(),
2293  darkMaskedImage=darkExposure.getMaskedImage(),
2294  expScale=expScale,
2295  darkScale=darkScale,
2296  invert=invert,
2297  trimToFit=self.config.doTrimToMatchCalib
2298  )
2299 
2300  def doLinearize(self, detector):
2301  """Check if linearization is needed for the detector cameraGeom.
2302 
2303  Checks config.doLinearize and the linearity type of the first
2304  amplifier.
2305 
2306  Parameters
2307  ----------
2308  detector : `lsst.afw.cameraGeom.Detector`
2309  Detector to get linearity type from.
2310 
2311  Returns
2312  -------
2313  doLinearize : `Bool`
2314  If True, linearization should be performed.
2315  """
2316  return self.config.doLinearize and \
2317  detector.getAmplifiers()[0].getLinearityType() != NullLinearityType
2318 
2319  def flatCorrection(self, exposure, flatExposure, invert=False):
2320  """Apply flat correction in place.
2321 
2322  Parameters
2323  ----------
2324  exposure : `lsst.afw.image.Exposure`
2325  Exposure to process.
2326  flatExposure : `lsst.afw.image.Exposure`
2327  Flat exposure of the same size as ``exposure``.
2328  invert : `Bool`, optional
2329  If True, unflatten an already flattened image.
2330 
2331  See Also
2332  --------
2333  lsst.ip.isr.isrFunctions.flatCorrection
2334  """
2335  isrFunctions.flatCorrection(
2336  maskedImage=exposure.getMaskedImage(),
2337  flatMaskedImage=flatExposure.getMaskedImage(),
2338  scalingType=self.config.flatScalingType,
2339  userScale=self.config.flatUserScale,
2340  invert=invert,
2341  trimToFit=self.config.doTrimToMatchCalib
2342  )
2343 
2344  def saturationDetection(self, exposure, amp):
2345  """Detect and mask saturated pixels in config.saturatedMaskName.
2346 
2347  Parameters
2348  ----------
2349  exposure : `lsst.afw.image.Exposure`
2350  Exposure to process. Only the amplifier DataSec is processed.
2351  amp : `lsst.afw.table.AmpInfoCatalog`
2352  Amplifier detector data.
2353 
2354  See Also
2355  --------
2356  lsst.ip.isr.isrFunctions.makeThresholdMask
2357  """
2358  if not math.isnan(amp.getSaturation()):
2359  maskedImage = exposure.getMaskedImage()
2360  dataView = maskedImage.Factory(maskedImage, amp.getRawBBox())
2361  isrFunctions.makeThresholdMask(
2362  maskedImage=dataView,
2363  threshold=amp.getSaturation(),
2364  growFootprints=0,
2365  maskName=self.config.saturatedMaskName,
2366  )
2367 
2368  def saturationInterpolation(self, exposure):
2369  """Interpolate over saturated pixels, in place.
2370 
2371  This method should be called after `saturationDetection`, to
2372  ensure that the saturated pixels have been identified in the
2373  SAT mask. It should also be called after `assembleCcd`, since
2374  saturated regions may cross amplifier boundaries.
2375 
2376  Parameters
2377  ----------
2378  exposure : `lsst.afw.image.Exposure`
2379  Exposure to process.
2380 
2381  See Also
2382  --------
2383  lsst.ip.isr.isrTask.saturationDetection
2384  lsst.ip.isr.isrFunctions.interpolateFromMask
2385  """
2386  isrFunctions.interpolateFromMask(
2387  maskedImage=exposure.getMaskedImage(),
2388  fwhm=self.config.fwhm,
2389  growSaturatedFootprints=self.config.growSaturationFootprintSize,
2390  maskNameList=list(self.config.saturatedMaskName),
2391  )
2392 
2393  def suspectDetection(self, exposure, amp):
2394  """Detect and mask suspect pixels in config.suspectMaskName.
2395 
2396  Parameters
2397  ----------
2398  exposure : `lsst.afw.image.Exposure`
2399  Exposure to process. Only the amplifier DataSec is processed.
2400  amp : `lsst.afw.table.AmpInfoCatalog`
2401  Amplifier detector data.
2402 
2403  See Also
2404  --------
2405  lsst.ip.isr.isrFunctions.makeThresholdMask
2406 
2407  Notes
2408  -----
2409  Suspect pixels are pixels whose value is greater than
2410  amp.getSuspectLevel(). This is intended to indicate pixels that may be
2411  affected by unknown systematics; for example if non-linearity
2412  corrections above a certain level are unstable then that would be a
2413  useful value for suspectLevel. A value of `nan` indicates that no such
2414  level exists and no pixels are to be masked as suspicious.
2415  """
2416  suspectLevel = amp.getSuspectLevel()
2417  if math.isnan(suspectLevel):
2418  return
2419 
2420  maskedImage = exposure.getMaskedImage()
2421  dataView = maskedImage.Factory(maskedImage, amp.getRawBBox())
2422  isrFunctions.makeThresholdMask(
2423  maskedImage=dataView,
2424  threshold=suspectLevel,
2425  growFootprints=0,
2426  maskName=self.config.suspectMaskName,
2427  )
2428 
2429  def maskDefect(self, exposure, defectBaseList):
2430  """Mask defects using mask plane "BAD", in place.
2431 
2432  Parameters
2433  ----------
2434  exposure : `lsst.afw.image.Exposure`
2435  Exposure to process.
2436  defectBaseList : `lsst.ip.isr.Defects` or `list` of
2437  `lsst.afw.image.DefectBase`.
2438  List of defects to mask.
2439 
2440  Notes
2441  -----
2442  Call this after CCD assembly, since defects may cross amplifier
2443  boundaries.
2444  """
2445  maskedImage = exposure.getMaskedImage()
2446  if not isinstance(defectBaseList, Defects):
2447  # Promotes DefectBase to Defect
2448  defectList = Defects(defectBaseList)
2449  else:
2450  defectList = defectBaseList
2451  defectList.maskPixels(maskedImage, maskName="BAD")
2452 
2453  def maskEdges(self, exposure, numEdgePixels=0, maskPlane="SUSPECT", level='DETECTOR'):
2454  """Mask edge pixels with applicable mask plane.
2455 
2456  Parameters
2457  ----------
2458  exposure : `lsst.afw.image.Exposure`
2459  Exposure to process.
2460  numEdgePixels : `int`, optional
2461  Number of edge pixels to mask.
2462  maskPlane : `str`, optional
2463  Mask plane name to use.
2464  level : `str`, optional
2465  Level at which to mask edges.
2466  """
2467  maskedImage = exposure.getMaskedImage()
2468  maskBitMask = maskedImage.getMask().getPlaneBitMask(maskPlane)
2469 
2470  if numEdgePixels > 0:
2471  if level == 'DETECTOR':
2472  boxes = [maskedImage.getBBox()]
2473  elif level == 'AMP':
2474  boxes = [amp.getBBox() for amp in exposure.getDetector()]
2475 
2476  for box in boxes:
2477  # This makes a bbox numEdgeSuspect pixels smaller than the
2478  # image on each side
2479  subImage = maskedImage[box]
2480  box.grow(-numEdgePixels)
2481  # Mask pixels outside box
2482  SourceDetectionTask.setEdgeBits(
2483  subImage,
2484  box,
2485  maskBitMask)
2486 
2487  def maskAndInterpolateDefects(self, exposure, defectBaseList):
2488  """Mask and interpolate defects using mask plane "BAD", in place.
2489 
2490  Parameters
2491  ----------
2492  exposure : `lsst.afw.image.Exposure`
2493  Exposure to process.
2494  defectBaseList : `lsst.ip.isr.Defects` or `list` of
2495  `lsst.afw.image.DefectBase`.
2496  List of defects to mask and interpolate.
2497 
2498  See Also
2499  --------
2500  lsst.ip.isr.isrTask.maskDefect
2501  """
2502  self.maskDefectmaskDefect(exposure, defectBaseList)
2503  self.maskEdgesmaskEdges(exposure, numEdgePixels=self.config.numEdgeSuspect,
2504  maskPlane="SUSPECT", level=self.config.edgeMaskLevel)
2505  isrFunctions.interpolateFromMask(
2506  maskedImage=exposure.getMaskedImage(),
2507  fwhm=self.config.fwhm,
2508  growSaturatedFootprints=0,
2509  maskNameList=["BAD"],
2510  )
2511 
2512  def maskNan(self, exposure):
2513  """Mask NaNs using mask plane "UNMASKEDNAN", in place.
2514 
2515  Parameters
2516  ----------
2517  exposure : `lsst.afw.image.Exposure`
2518  Exposure to process.
2519 
2520  Notes
2521  -----
2522  We mask over all non-finite values (NaN, inf), including those
2523  that are masked with other bits (because those may or may not be
2524  interpolated over later, and we want to remove all NaN/infs).
2525  Despite this behaviour, the "UNMASKEDNAN" mask plane is used to
2526  preserve the historical name.
2527  """
2528  maskedImage = exposure.getMaskedImage()
2529 
2530  # Find and mask NaNs
2531  maskedImage.getMask().addMaskPlane("UNMASKEDNAN")
2532  maskVal = maskedImage.getMask().getPlaneBitMask("UNMASKEDNAN")
2533  numNans = maskNans(maskedImage, maskVal)
2534  self.metadata.set("NUMNANS", numNans)
2535  if numNans > 0:
2536  self.log.warning("There were %d unmasked NaNs.", numNans)
2537 
2538  def maskAndInterpolateNan(self, exposure):
2539  """"Mask and interpolate NaN/infs using mask plane "UNMASKEDNAN",
2540  in place.
2541 
2542  Parameters
2543  ----------
2544  exposure : `lsst.afw.image.Exposure`
2545  Exposure to process.
2546 
2547  See Also
2548  --------
2549  lsst.ip.isr.isrTask.maskNan
2550  """
2551  self.maskNanmaskNan(exposure)
2552  isrFunctions.interpolateFromMask(
2553  maskedImage=exposure.getMaskedImage(),
2554  fwhm=self.config.fwhm,
2555  growSaturatedFootprints=0,
2556  maskNameList=["UNMASKEDNAN"],
2557  )
2558 
2559  def measureBackground(self, exposure, IsrQaConfig=None):
2560  """Measure the image background in subgrids, for quality control.
2561 
2562  Parameters
2563  ----------
2564  exposure : `lsst.afw.image.Exposure`
2565  Exposure to process.
2566  IsrQaConfig : `lsst.ip.isr.isrQa.IsrQaConfig`
2567  Configuration object containing parameters on which background
2568  statistics and subgrids to use.
2569  """
2570  if IsrQaConfig is not None:
2571  statsControl = afwMath.StatisticsControl(IsrQaConfig.flatness.clipSigma,
2572  IsrQaConfig.flatness.nIter)
2573  maskVal = exposure.getMaskedImage().getMask().getPlaneBitMask(["BAD", "SAT", "DETECTED"])
2574  statsControl.setAndMask(maskVal)
2575  maskedImage = exposure.getMaskedImage()
2576  stats = afwMath.makeStatistics(maskedImage, afwMath.MEDIAN | afwMath.STDEVCLIP, statsControl)
2577  skyLevel = stats.getValue(afwMath.MEDIAN)
2578  skySigma = stats.getValue(afwMath.STDEVCLIP)
2579  self.log.info("Flattened sky level: %f +/- %f.", skyLevel, skySigma)
2580  metadata = exposure.getMetadata()
2581  metadata.set('SKYLEVEL', skyLevel)
2582  metadata.set('SKYSIGMA', skySigma)
2583 
2584  # calcluating flatlevel over the subgrids
2585  stat = afwMath.MEANCLIP if IsrQaConfig.flatness.doClip else afwMath.MEAN
2586  meshXHalf = int(IsrQaConfig.flatness.meshX/2.)
2587  meshYHalf = int(IsrQaConfig.flatness.meshY/2.)
2588  nX = int((exposure.getWidth() + meshXHalf) / IsrQaConfig.flatness.meshX)
2589  nY = int((exposure.getHeight() + meshYHalf) / IsrQaConfig.flatness.meshY)
2590  skyLevels = numpy.zeros((nX, nY))
2591 
2592  for j in range(nY):
2593  yc = meshYHalf + j * IsrQaConfig.flatness.meshY
2594  for i in range(nX):
2595  xc = meshXHalf + i * IsrQaConfig.flatness.meshX
2596 
2597  xLLC = xc - meshXHalf
2598  yLLC = yc - meshYHalf
2599  xURC = xc + meshXHalf - 1
2600  yURC = yc + meshYHalf - 1
2601 
2602  bbox = lsst.geom.Box2I(lsst.geom.Point2I(xLLC, yLLC), lsst.geom.Point2I(xURC, yURC))
2603  miMesh = maskedImage.Factory(exposure.getMaskedImage(), bbox, afwImage.LOCAL)
2604 
2605  skyLevels[i, j] = afwMath.makeStatistics(miMesh, stat, statsControl).getValue()
2606 
2607  good = numpy.where(numpy.isfinite(skyLevels))
2608  skyMedian = numpy.median(skyLevels[good])
2609  flatness = (skyLevels[good] - skyMedian) / skyMedian
2610  flatness_rms = numpy.std(flatness)
2611  flatness_pp = flatness.max() - flatness.min() if len(flatness) > 0 else numpy.nan
2612 
2613  self.log.info("Measuring sky levels in %dx%d grids: %f.", nX, nY, skyMedian)
2614  self.log.info("Sky flatness in %dx%d grids - pp: %f rms: %f.",
2615  nX, nY, flatness_pp, flatness_rms)
2616 
2617  metadata.set('FLATNESS_PP', float(flatness_pp))
2618  metadata.set('FLATNESS_RMS', float(flatness_rms))
2619  metadata.set('FLATNESS_NGRIDS', '%dx%d' % (nX, nY))
2620  metadata.set('FLATNESS_MESHX', IsrQaConfig.flatness.meshX)
2621  metadata.set('FLATNESS_MESHY', IsrQaConfig.flatness.meshY)
2622 
2623  def roughZeroPoint(self, exposure):
2624  """Set an approximate magnitude zero point for the exposure.
2625 
2626  Parameters
2627  ----------
2628  exposure : `lsst.afw.image.Exposure`
2629  Exposure to process.
2630  """
2631  filterLabel = exposure.getFilterLabel()
2632  physicalFilter = isrFunctions.getPhysicalFilter(filterLabel, self.log)
2633 
2634  if physicalFilter in self.config.fluxMag0T1:
2635  fluxMag0 = self.config.fluxMag0T1[physicalFilter]
2636  else:
2637  self.log.warning("No rough magnitude zero point defined for filter %s.", physicalFilter)
2638  fluxMag0 = self.config.defaultFluxMag0T1
2639 
2640  expTime = exposure.getInfo().getVisitInfo().getExposureTime()
2641  if not expTime > 0: # handle NaN as well as <= 0
2642  self.log.warning("Non-positive exposure time; skipping rough zero point.")
2643  return
2644 
2645  self.log.info("Setting rough magnitude zero point for filter %s: %f",
2646  physicalFilter, 2.5*math.log10(fluxMag0*expTime))
2647  exposure.setPhotoCalib(afwImage.makePhotoCalibFromCalibZeroPoint(fluxMag0*expTime, 0.0))
2648 
2649  def setValidPolygonIntersect(self, ccdExposure, fpPolygon):
2650  """Set valid polygon as the intersection of fpPolygon and chip corners.
2651 
2652  Parameters
2653  ----------
2654  ccdExposure : `lsst.afw.image.Exposure`
2655  Exposure to process.
2656  fpPolygon : `lsst.afw.geom.Polygon`
2657  Polygon in focal plane coordinates.
2658  """
2659  # Get ccd corners in focal plane coordinates
2660  ccd = ccdExposure.getDetector()
2661  fpCorners = ccd.getCorners(FOCAL_PLANE)
2662  ccdPolygon = Polygon(fpCorners)
2663 
2664  # Get intersection of ccd corners with fpPolygon
2665  intersect = ccdPolygon.intersectionSingle(fpPolygon)
2666 
2667  # Transform back to pixel positions and build new polygon
2668  ccdPoints = ccd.transform(intersect, FOCAL_PLANE, PIXELS)
2669  validPolygon = Polygon(ccdPoints)
2670  ccdExposure.getInfo().setValidPolygon(validPolygon)
2671 
2672  @contextmanager
2673  def flatContext(self, exp, flat, dark=None):
2674  """Context manager that applies and removes flats and darks,
2675  if the task is configured to apply them.
2676 
2677  Parameters
2678  ----------
2679  exp : `lsst.afw.image.Exposure`
2680  Exposure to process.
2681  flat : `lsst.afw.image.Exposure`
2682  Flat exposure the same size as ``exp``.
2683  dark : `lsst.afw.image.Exposure`, optional
2684  Dark exposure the same size as ``exp``.
2685 
2686  Yields
2687  ------
2688  exp : `lsst.afw.image.Exposure`
2689  The flat and dark corrected exposure.
2690  """
2691  if self.config.doDark and dark is not None:
2692  self.darkCorrectiondarkCorrection(exp, dark)
2693  if self.config.doFlat:
2694  self.flatCorrectionflatCorrection(exp, flat)
2695  try:
2696  yield exp
2697  finally:
2698  if self.config.doFlat:
2699  self.flatCorrectionflatCorrection(exp, flat, invert=True)
2700  if self.config.doDark and dark is not None:
2701  self.darkCorrectiondarkCorrection(exp, dark, invert=True)
2702 
2703  def debugView(self, exposure, stepname):
2704  """Utility function to examine ISR exposure at different stages.
2705 
2706  Parameters
2707  ----------
2708  exposure : `lsst.afw.image.Exposure`
2709  Exposure to view.
2710  stepname : `str`
2711  State of processing to view.
2712  """
2713  frame = getDebugFrame(self._display, stepname)
2714  if frame:
2715  display = getDisplay(frame)
2716  display.scale('asinh', 'zscale')
2717  display.mtv(exposure)
2718  prompt = "Press Enter to continue [c]... "
2719  while True:
2720  ans = input(prompt).lower()
2721  if ans in ("", "c",):
2722  break
2723 
2724 
2725 class FakeAmp(object):
2726  """A Detector-like object that supports returning gain and saturation level
2727 
2728  This is used when the input exposure does not have a detector.
2729 
2730  Parameters
2731  ----------
2732  exposure : `lsst.afw.image.Exposure`
2733  Exposure to generate a fake amplifier for.
2734  config : `lsst.ip.isr.isrTaskConfig`
2735  Configuration to apply to the fake amplifier.
2736  """
2737 
2738  def __init__(self, exposure, config):
2739  self._bbox_bbox = exposure.getBBox(afwImage.LOCAL)
2740  self._RawHorizontalOverscanBBox_RawHorizontalOverscanBBox = lsst.geom.Box2I()
2741  self._gain_gain = config.gain
2742  self._readNoise_readNoise = config.readNoise
2743  self._saturation_saturation = config.saturation
2744 
2745  def getBBox(self):
2746  return self._bbox_bbox
2747 
2748  def getRawBBox(self):
2749  return self._bbox_bbox
2750 
2752  return self._RawHorizontalOverscanBBox_RawHorizontalOverscanBBox
2753 
2754  def getGain(self):
2755  return self._gain_gain
2756 
2757  def getReadNoise(self):
2758  return self._readNoise_readNoise
2759 
2760  def getSaturation(self):
2761  return self._saturation_saturation
2762 
2763  def getSuspectLevel(self):
2764  return float("NaN")
2765 
2766 
2767 class RunIsrConfig(pexConfig.Config):
2768  isr = pexConfig.ConfigurableField(target=IsrTask, doc="Instrument signature removal")
2769 
2770 
2771 class RunIsrTask(pipeBase.CmdLineTask):
2772  """Task to wrap the default IsrTask to allow it to be retargeted.
2773 
2774  The standard IsrTask can be called directly from a command line
2775  program, but doing so removes the ability of the task to be
2776  retargeted. As most cameras override some set of the IsrTask
2777  methods, this would remove those data-specific methods in the
2778  output post-ISR images. This wrapping class fixes the issue,
2779  allowing identical post-ISR images to be generated by both the
2780  processCcd and isrTask code.
2781  """
2782  ConfigClass = RunIsrConfig
2783  _DefaultName = "runIsr"
2784 
2785  def __init__(self, *args, **kwargs):
2786  super().__init__(*args, **kwargs)
2787  self.makeSubtask("isr")
2788 
2789  def runDataRef(self, dataRef):
2790  """
2791  Parameters
2792  ----------
2793  dataRef : `lsst.daf.persistence.ButlerDataRef`
2794  data reference of the detector data to be processed
2795 
2796  Returns
2797  -------
2798  result : `pipeBase.Struct`
2799  Result struct with component:
2800 
2801  - exposure : `lsst.afw.image.Exposure`
2802  Post-ISR processed exposure.
2803  """
2804  return self.isr.runDataRef(dataRef)
def getRawHorizontalOverscanBBox(self)
Definition: isrTask.py:2751
def __init__(self, exposure, config)
Definition: isrTask.py:2738
def __init__(self, *config=None)
Definition: isrTask.py:274
def flatCorrection(self, exposure, flatExposure, invert=False)
Definition: isrTask.py:2319
def maskAndInterpolateNan(self, exposure)
Definition: isrTask.py:2538
def saturationInterpolation(self, exposure)
Definition: isrTask.py:2368
def runDataRef(self, sensorRef)
Definition: isrTask.py:1773
def maskNan(self, exposure)
Definition: isrTask.py:2512
def maskAmplifier(self, ccdExposure, amp, defects)
Definition: isrTask.py:1961
def debugView(self, exposure, stepname)
Definition: isrTask.py:2703
def ensureExposure(self, inputExp, camera=None, detectorNum=None)
Definition: isrTask.py:1872
def getIsrExposure(self, dataRef, datasetType, dateObs=None, immediate=True)
Definition: isrTask.py:1824
def maskNegativeVariance(self, exposure)
Definition: isrTask.py:2241
def saturationDetection(self, exposure, amp)
Definition: isrTask.py:2344
def maskDefect(self, exposure, defectBaseList)
Definition: isrTask.py:2429
def __init__(self, **kwargs)
Definition: isrTask.py:976
def runQuantum(self, butlerQC, inputRefs, outputRefs)
Definition: isrTask.py:987
def maskEdges(self, exposure, numEdgePixels=0, maskPlane="SUSPECT", level='DETECTOR')
Definition: isrTask.py:2453
def overscanCorrection(self, ccdExposure, amp)
Definition: isrTask.py:2036
def measureBackground(self, exposure, IsrQaConfig=None)
Definition: isrTask.py:2559
def roughZeroPoint(self, exposure)
Definition: isrTask.py:2623
def maskAndInterpolateDefects(self, exposure, defectBaseList)
Definition: isrTask.py:2487
def setValidPolygonIntersect(self, ccdExposure, fpPolygon)
Definition: isrTask.py:2649
def readIsrData(self, dataRef, rawExposure)
Definition: isrTask.py:1086
def run(self, ccdExposure, *camera=None, bias=None, linearizer=None, crosstalk=None, crosstalkSources=None, dark=None, flat=None, ptc=None, bfKernel=None, bfGains=None, defects=None, fringes=pipeBase.Struct(fringes=None), opticsTransmission=None, filterTransmission=None, sensorTransmission=None, atmosphereTransmission=None, detectorNum=None, strayLightData=None, illumMaskedImage=None, isGen3=False)
Definition: isrTask.py:1273
def doLinearize(self, detector)
Definition: isrTask.py:2300
def flatContext(self, exp, flat, dark=None)
Definition: isrTask.py:2673
def convertIntToFloat(self, exposure)
Definition: isrTask.py:1924
def suspectDetection(self, exposure, amp)
Definition: isrTask.py:2393
def updateVariance(self, ampExposure, amp, overscanImage=None, ptcDataset=None)
Definition: isrTask.py:2166
def darkCorrection(self, exposure, darkExposure, invert=False)
Definition: isrTask.py:2257
def __init__(self, *args, **kwargs)
Definition: isrTask.py:2785
def runDataRef(self, dataRef)
Definition: isrTask.py:2789
def checkFilter(exposure, filterList, log)
def crosstalkSourceLookup(datasetType, registry, quantumDataId, collections)
Definition: isrTask.py:63
size_t maskNans(afw::image::MaskedImage< PixelT > const &mi, afw::image::MaskPixel maskVal, afw::image::MaskPixel allow=0)
Mask NANs in an image.
Definition: Isr.cc:35