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