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