lsst.ip.isr  21.0.0-7-gde99fe0+bda8df9b93
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  # TODO DM-28093: change this to: rawExposure.getFilterLabel().physicalLabel
1080  filterName = afwImage.Filter(rawExposure.getFilter().getId()).getName() # Canonical name for filter
1081  rawExposure.mask.addMaskPlane("UNMASKEDNAN") # needed to match pre DM-15862 processing.
1082  biasExposure = (self.getIsrExposuregetIsrExposure(dataRef, self.config.biasDataProductName)
1083  if self.config.doBias else None)
1084  # immediate=True required for functors and linearizers are functors; see ticket DM-6515
1085  linearizer = (dataRef.get("linearizer", immediate=True)
1086  if self.doLinearizedoLinearize(ccd) else None)
1087  if linearizer is not None and not isinstance(linearizer, numpy.ndarray):
1088  linearizer.log = self.log
1089  if isinstance(linearizer, numpy.ndarray):
1090  linearizer = linearize.Linearizer(table=linearizer, detector=ccd)
1091 
1092  crosstalkCalib = None
1093  if self.config.doCrosstalk:
1094  try:
1095  crosstalkCalib = dataRef.get("crosstalk", immediate=True)
1096  except NoResults:
1097  coeffVector = (self.config.crosstalk.crosstalkValues
1098  if self.config.crosstalk.useConfigCoefficients else None)
1099  crosstalkCalib = CrosstalkCalib().fromDetector(ccd, coeffVector=coeffVector)
1100  crosstalkSources = (self.crosstalk.prepCrosstalk(dataRef, crosstalkCalib)
1101  if self.config.doCrosstalk else None)
1102 
1103  darkExposure = (self.getIsrExposuregetIsrExposure(dataRef, self.config.darkDataProductName)
1104  if self.config.doDark else None)
1105  flatExposure = (self.getIsrExposuregetIsrExposure(dataRef, self.config.flatDataProductName,
1106  dateObs=dateObs)
1107  if self.config.doFlat else None)
1108 
1109  brighterFatterKernel = None
1110  brighterFatterGains = None
1111  if self.config.doBrighterFatter is True:
1112  try:
1113  # Use the new-style cp_pipe version of the kernel if it exists
1114  # If using a new-style kernel, always use the self-consistent
1115  # gains, i.e. the ones inside the kernel object itself
1116  brighterFatterKernel = dataRef.get("brighterFatterKernel")
1117  brighterFatterGains = brighterFatterKernel.gain
1118  self.log.info("New style bright-fatter kernel (brighterFatterKernel) loaded")
1119  except NoResults:
1120  try: # Fall back to the old-style numpy-ndarray style kernel if necessary.
1121  brighterFatterKernel = dataRef.get("bfKernel")
1122  self.log.info("Old style bright-fatter kernel (np.array) loaded")
1123  except NoResults:
1124  brighterFatterKernel = None
1125  if brighterFatterKernel is not None and not isinstance(brighterFatterKernel, numpy.ndarray):
1126  # If the kernel is not an ndarray, it's the cp_pipe version
1127  # so extract the kernel for this detector, or raise an error
1128  if self.config.brighterFatterLevel == 'DETECTOR':
1129  if brighterFatterKernel.detectorKernel:
1130  brighterFatterKernel = brighterFatterKernel.detectorKernel[ccd.getId()]
1131  elif brighterFatterKernel.detectorKernelFromAmpKernels:
1132  brighterFatterKernel = brighterFatterKernel.detectorKernelFromAmpKernels[ccd.getId()]
1133  else:
1134  raise RuntimeError("Failed to extract kernel from new-style BF kernel.")
1135  else:
1136  # TODO DM-15631 for implementing this
1137  raise NotImplementedError("Per-amplifier brighter-fatter correction not implemented")
1138 
1139  defectList = (dataRef.get("defects")
1140  if self.config.doDefect else None)
1141  fringeStruct = (self.fringe.readFringes(dataRef, assembler=self.assembleCcd
1142  if self.config.doAssembleIsrExposures else None)
1143  if self.config.doFringe and self.fringe.checkFilter(rawExposure)
1144  else pipeBase.Struct(fringes=None))
1145 
1146  if self.config.doAttachTransmissionCurve:
1147  opticsTransmission = (dataRef.get("transmission_optics")
1148  if self.config.doUseOpticsTransmission else None)
1149  filterTransmission = (dataRef.get("transmission_filter")
1150  if self.config.doUseFilterTransmission else None)
1151  sensorTransmission = (dataRef.get("transmission_sensor")
1152  if self.config.doUseSensorTransmission else None)
1153  atmosphereTransmission = (dataRef.get("transmission_atmosphere")
1154  if self.config.doUseAtmosphereTransmission else None)
1155  else:
1156  opticsTransmission = None
1157  filterTransmission = None
1158  sensorTransmission = None
1159  atmosphereTransmission = None
1160 
1161  if self.config.doStrayLight:
1162  strayLightData = self.strayLight.readIsrData(dataRef, rawExposure)
1163  else:
1164  strayLightData = None
1165 
1166  illumMaskedImage = (self.getIsrExposuregetIsrExposure(dataRef,
1167  self.config.illuminationCorrectionDataProductName).getMaskedImage()
1168  if (self.config.doIlluminationCorrection
1169  and filterName in self.config.illumFilters)
1170  else None)
1171 
1172  # Struct should include only kwargs to run()
1173  return pipeBase.Struct(bias=biasExposure,
1174  linearizer=linearizer,
1175  crosstalk=crosstalkCalib,
1176  crosstalkSources=crosstalkSources,
1177  dark=darkExposure,
1178  flat=flatExposure,
1179  bfKernel=brighterFatterKernel,
1180  bfGains=brighterFatterGains,
1181  defects=defectList,
1182  fringes=fringeStruct,
1183  opticsTransmission=opticsTransmission,
1184  filterTransmission=filterTransmission,
1185  sensorTransmission=sensorTransmission,
1186  atmosphereTransmission=atmosphereTransmission,
1187  strayLightData=strayLightData,
1188  illumMaskedImage=illumMaskedImage
1189  )
1190 
1191  @pipeBase.timeMethod
1192  def run(self, ccdExposure, camera=None, bias=None, linearizer=None,
1193  crosstalk=None, crosstalkSources=None,
1194  dark=None, flat=None, bfKernel=None, bfGains=None, defects=None,
1195  fringes=pipeBase.Struct(fringes=None), opticsTransmission=None, filterTransmission=None,
1196  sensorTransmission=None, atmosphereTransmission=None,
1197  detectorNum=None, strayLightData=None, illumMaskedImage=None,
1198  isGen3=False,
1199  ):
1200  """Perform instrument signature removal on an exposure.
1201 
1202  Steps included in the ISR processing, in order performed, are:
1203  - saturation and suspect pixel masking
1204  - overscan subtraction
1205  - CCD assembly of individual amplifiers
1206  - bias subtraction
1207  - variance image construction
1208  - linearization of non-linear response
1209  - crosstalk masking
1210  - brighter-fatter correction
1211  - dark subtraction
1212  - fringe correction
1213  - stray light subtraction
1214  - flat correction
1215  - masking of known defects and camera specific features
1216  - vignette calculation
1217  - appending transmission curve and distortion model
1218 
1219  Parameters
1220  ----------
1221  ccdExposure : `lsst.afw.image.Exposure`
1222  The raw exposure that is to be run through ISR. The
1223  exposure is modified by this method.
1224  camera : `lsst.afw.cameraGeom.Camera`, optional
1225  The camera geometry for this exposure. Required if ``isGen3`` is
1226  `True` and one or more of ``ccdExposure``, ``bias``, ``dark``, or
1227  ``flat`` does not have an associated detector.
1228  bias : `lsst.afw.image.Exposure`, optional
1229  Bias calibration frame.
1230  linearizer : `lsst.ip.isr.linearize.LinearizeBase`, optional
1231  Functor for linearization.
1232  crosstalk : `lsst.ip.isr.crosstalk.CrosstalkCalib`, optional
1233  Calibration for crosstalk.
1234  crosstalkSources : `list`, optional
1235  List of possible crosstalk sources.
1236  dark : `lsst.afw.image.Exposure`, optional
1237  Dark calibration frame.
1238  flat : `lsst.afw.image.Exposure`, optional
1239  Flat calibration frame.
1240  bfKernel : `numpy.ndarray`, optional
1241  Brighter-fatter kernel.
1242  bfGains : `dict` of `float`, optional
1243  Gains used to override the detector's nominal gains for the
1244  brighter-fatter correction. A dict keyed by amplifier name for
1245  the detector in question.
1246  defects : `lsst.ip.isr.Defects`, optional
1247  List of defects.
1248  fringes : `lsst.pipe.base.Struct`, optional
1249  Struct containing the fringe correction data, with
1250  elements:
1251  - ``fringes``: fringe calibration frame (`afw.image.Exposure`)
1252  - ``seed``: random seed derived from the ccdExposureId for random
1253  number generator (`uint32`)
1254  opticsTransmission: `lsst.afw.image.TransmissionCurve`, optional
1255  A ``TransmissionCurve`` that represents the throughput of the optics,
1256  to be evaluated in focal-plane coordinates.
1257  filterTransmission : `lsst.afw.image.TransmissionCurve`
1258  A ``TransmissionCurve`` that represents the throughput of the filter
1259  itself, to be evaluated in focal-plane coordinates.
1260  sensorTransmission : `lsst.afw.image.TransmissionCurve`
1261  A ``TransmissionCurve`` that represents the throughput of the sensor
1262  itself, to be evaluated in post-assembly trimmed detector coordinates.
1263  atmosphereTransmission : `lsst.afw.image.TransmissionCurve`
1264  A ``TransmissionCurve`` that represents the throughput of the
1265  atmosphere, assumed to be spatially constant.
1266  detectorNum : `int`, optional
1267  The integer number for the detector to process.
1268  isGen3 : bool, optional
1269  Flag this call to run() as using the Gen3 butler environment.
1270  strayLightData : `object`, optional
1271  Opaque object containing calibration information for stray-light
1272  correction. If `None`, no correction will be performed.
1273  illumMaskedImage : `lsst.afw.image.MaskedImage`, optional
1274  Illumination correction image.
1275 
1276  Returns
1277  -------
1278  result : `lsst.pipe.base.Struct`
1279  Result struct with component:
1280  - ``exposure`` : `afw.image.Exposure`
1281  The fully ISR corrected exposure.
1282  - ``outputExposure`` : `afw.image.Exposure`
1283  An alias for `exposure`
1284  - ``ossThumb`` : `numpy.ndarray`
1285  Thumbnail image of the exposure after overscan subtraction.
1286  - ``flattenedThumb`` : `numpy.ndarray`
1287  Thumbnail image of the exposure after flat-field correction.
1288 
1289  Raises
1290  ------
1291  RuntimeError
1292  Raised if a configuration option is set to True, but the
1293  required calibration data has not been specified.
1294 
1295  Notes
1296  -----
1297  The current processed exposure can be viewed by setting the
1298  appropriate lsstDebug entries in the `debug.display`
1299  dictionary. The names of these entries correspond to some of
1300  the IsrTaskConfig Boolean options, with the value denoting the
1301  frame to use. The exposure is shown inside the matching
1302  option check and after the processing of that step has
1303  finished. The steps with debug points are:
1304 
1305  doAssembleCcd
1306  doBias
1307  doCrosstalk
1308  doBrighterFatter
1309  doDark
1310  doFringe
1311  doStrayLight
1312  doFlat
1313 
1314  In addition, setting the "postISRCCD" entry displays the
1315  exposure after all ISR processing has finished.
1316 
1317  """
1318 
1319  if isGen3 is True:
1320  # Gen3 currently cannot automatically do configuration overrides.
1321  # DM-15257 looks to discuss this issue.
1322  # Configure input exposures;
1323  if detectorNum is None:
1324  raise RuntimeError("Must supply the detectorNum if running as Gen3.")
1325 
1326  ccdExposure = self.ensureExposureensureExposure(ccdExposure, camera, detectorNum)
1327  bias = self.ensureExposureensureExposure(bias, camera, detectorNum)
1328  dark = self.ensureExposureensureExposure(dark, camera, detectorNum)
1329  flat = self.ensureExposureensureExposure(flat, camera, detectorNum)
1330  else:
1331  if isinstance(ccdExposure, ButlerDataRef):
1332  return self.runDataRefrunDataRef(ccdExposure)
1333 
1334  ccd = ccdExposure.getDetector()
1335  # TODO DM-28093: change this to: ccdExposure.getFilterLabel().physicalLabel
1336  filterName = afwImage.Filter(ccdExposure.getFilter().getId()).getName() # Canonical name for filter
1337 
1338  if not ccd:
1339  assert not self.config.doAssembleCcd, "You need a Detector to run assembleCcd."
1340  ccd = [FakeAmp(ccdExposure, self.config)]
1341 
1342  # Validate Input
1343  if self.config.doBias and bias is None:
1344  raise RuntimeError("Must supply a bias exposure if config.doBias=True.")
1345  if self.doLinearizedoLinearize(ccd) and linearizer is None:
1346  raise RuntimeError("Must supply a linearizer if config.doLinearize=True for this detector.")
1347  if self.config.doBrighterFatter and bfKernel is None:
1348  raise RuntimeError("Must supply a kernel if config.doBrighterFatter=True.")
1349  if self.config.doDark and dark is None:
1350  raise RuntimeError("Must supply a dark exposure if config.doDark=True.")
1351  if self.config.doFlat and flat is None:
1352  raise RuntimeError("Must supply a flat exposure if config.doFlat=True.")
1353  if self.config.doDefect and defects is None:
1354  raise RuntimeError("Must supply defects if config.doDefect=True.")
1355  if (self.config.doFringe and filterName in self.fringe.config.filters
1356  and fringes.fringes is None):
1357  # The `fringes` object needs to be a pipeBase.Struct, as
1358  # we use it as a `dict` for the parameters of
1359  # `FringeTask.run()`. The `fringes.fringes` `list` may
1360  # not be `None` if `doFringe=True`. Otherwise, raise.
1361  raise RuntimeError("Must supply fringe exposure as a pipeBase.Struct.")
1362  if (self.config.doIlluminationCorrection and filterName in self.config.illumFilters
1363  and illumMaskedImage is None):
1364  raise RuntimeError("Must supply an illumcor if config.doIlluminationCorrection=True.")
1365 
1366  # Begin ISR processing.
1367  if self.config.doConvertIntToFloat:
1368  self.log.info("Converting exposure to floating point values.")
1369  ccdExposure = self.convertIntToFloatconvertIntToFloat(ccdExposure)
1370 
1371  if self.config.doBias and self.config.doBiasBeforeOverscan:
1372  self.log.info("Applying bias correction.")
1373  isrFunctions.biasCorrection(ccdExposure.getMaskedImage(), bias.getMaskedImage(),
1374  trimToFit=self.config.doTrimToMatchCalib)
1375  self.debugViewdebugView(ccdExposure, "doBias")
1376 
1377  # Amplifier level processing.
1378  overscans = []
1379  for amp in ccd:
1380  # if ccdExposure is one amp, check for coverage to prevent performing ops multiple times
1381  if ccdExposure.getBBox().contains(amp.getBBox()):
1382  # Check for fully masked bad amplifiers, and generate masks for SUSPECT and SATURATED values.
1383  badAmp = self.maskAmplifiermaskAmplifier(ccdExposure, amp, defects)
1384 
1385  if self.config.doOverscan and not badAmp:
1386  # Overscan correction on amp-by-amp basis.
1387  overscanResults = self.overscanCorrectionoverscanCorrection(ccdExposure, amp)
1388  self.log.debug("Corrected overscan for amplifier %s.", amp.getName())
1389  if overscanResults is not None and \
1390  self.config.qa is not None and self.config.qa.saveStats is True:
1391  if isinstance(overscanResults.overscanFit, float):
1392  qaMedian = overscanResults.overscanFit
1393  qaStdev = float("NaN")
1394  else:
1395  qaStats = afwMath.makeStatistics(overscanResults.overscanFit,
1396  afwMath.MEDIAN | afwMath.STDEVCLIP)
1397  qaMedian = qaStats.getValue(afwMath.MEDIAN)
1398  qaStdev = qaStats.getValue(afwMath.STDEVCLIP)
1399 
1400  self.metadata.set(f"FIT MEDIAN {amp.getName()}", qaMedian)
1401  self.metadata.set(f"FIT STDEV {amp.getName()}", qaStdev)
1402  self.log.debug(" Overscan stats for amplifer %s: %f +/- %f",
1403  amp.getName(), qaMedian, qaStdev)
1404 
1405  # Residuals after overscan correction
1406  qaStatsAfter = afwMath.makeStatistics(overscanResults.overscanImage,
1407  afwMath.MEDIAN | afwMath.STDEVCLIP)
1408  qaMedianAfter = qaStatsAfter.getValue(afwMath.MEDIAN)
1409  qaStdevAfter = qaStatsAfter.getValue(afwMath.STDEVCLIP)
1410 
1411  self.metadata.set(f"RESIDUAL MEDIAN {amp.getName()}", qaMedianAfter)
1412  self.metadata.set(f"RESIDUAL STDEV {amp.getName()}", qaStdevAfter)
1413  self.log.debug(" Overscan stats for amplifer %s after correction: %f +/- %f",
1414  amp.getName(), qaMedianAfter, qaStdevAfter)
1415 
1416  ccdExposure.getMetadata().set('OVERSCAN', "Overscan corrected")
1417  else:
1418  if badAmp:
1419  self.log.warn("Amplifier %s is bad.", amp.getName())
1420  overscanResults = None
1421 
1422  overscans.append(overscanResults if overscanResults is not None else None)
1423  else:
1424  self.log.info("Skipped OSCAN for %s.", amp.getName())
1425 
1426  if self.config.doCrosstalk and self.config.doCrosstalkBeforeAssemble:
1427  self.log.info("Applying crosstalk correction.")
1428  self.crosstalk.run(ccdExposure, crosstalk=crosstalk,
1429  crosstalkSources=crosstalkSources)
1430  self.debugViewdebugView(ccdExposure, "doCrosstalk")
1431 
1432  if self.config.doAssembleCcd:
1433  self.log.info("Assembling CCD from amplifiers.")
1434  ccdExposure = self.assembleCcd.assembleCcd(ccdExposure)
1435 
1436  if self.config.expectWcs and not ccdExposure.getWcs():
1437  self.log.warn("No WCS found in input exposure.")
1438  self.debugViewdebugView(ccdExposure, "doAssembleCcd")
1439 
1440  ossThumb = None
1441  if self.config.qa.doThumbnailOss:
1442  ossThumb = isrQa.makeThumbnail(ccdExposure, isrQaConfig=self.config.qa)
1443 
1444  if self.config.doBias and not self.config.doBiasBeforeOverscan:
1445  self.log.info("Applying bias correction.")
1446  isrFunctions.biasCorrection(ccdExposure.getMaskedImage(), bias.getMaskedImage(),
1447  trimToFit=self.config.doTrimToMatchCalib)
1448  self.debugViewdebugView(ccdExposure, "doBias")
1449 
1450  if self.config.doVariance:
1451  for amp, overscanResults in zip(ccd, overscans):
1452  if ccdExposure.getBBox().contains(amp.getBBox()):
1453  self.log.debug("Constructing variance map for amplifer %s.", amp.getName())
1454  ampExposure = ccdExposure.Factory(ccdExposure, amp.getBBox())
1455  if overscanResults is not None:
1456  self.updateVarianceupdateVariance(ampExposure, amp,
1457  overscanImage=overscanResults.overscanImage)
1458  else:
1459  self.updateVarianceupdateVariance(ampExposure, amp,
1460  overscanImage=None)
1461  if self.config.qa is not None and self.config.qa.saveStats is True:
1462  qaStats = afwMath.makeStatistics(ampExposure.getVariance(),
1463  afwMath.MEDIAN | afwMath.STDEVCLIP)
1464  self.metadata.set(f"ISR VARIANCE {amp.getName()} MEDIAN",
1465  qaStats.getValue(afwMath.MEDIAN))
1466  self.metadata.set(f"ISR VARIANCE {amp.getName()} STDEV",
1467  qaStats.getValue(afwMath.STDEVCLIP))
1468  self.log.debug(" Variance stats for amplifer %s: %f +/- %f.",
1469  amp.getName(), qaStats.getValue(afwMath.MEDIAN),
1470  qaStats.getValue(afwMath.STDEVCLIP))
1471 
1472  if self.doLinearizedoLinearize(ccd):
1473  self.log.info("Applying linearizer.")
1474  linearizer.applyLinearity(image=ccdExposure.getMaskedImage().getImage(),
1475  detector=ccd, log=self.log)
1476 
1477  if self.config.doCrosstalk and not self.config.doCrosstalkBeforeAssemble:
1478  self.log.info("Applying crosstalk correction.")
1479  self.crosstalk.run(ccdExposure, crosstalk=crosstalk,
1480  crosstalkSources=crosstalkSources, isTrimmed=True)
1481  self.debugViewdebugView(ccdExposure, "doCrosstalk")
1482 
1483  # Masking block. Optionally mask known defects, NAN pixels, widen trails, and do
1484  # anything else the camera needs. Saturated and suspect pixels have already been masked.
1485  if self.config.doDefect:
1486  self.log.info("Masking defects.")
1487  self.maskDefectmaskDefect(ccdExposure, defects)
1488 
1489  if self.config.numEdgeSuspect > 0:
1490  self.log.info("Masking edges as SUSPECT.")
1491  self.maskEdgesmaskEdges(ccdExposure, numEdgePixels=self.config.numEdgeSuspect,
1492  maskPlane="SUSPECT", level=self.config.edgeMaskLevel)
1493 
1494  if self.config.doNanMasking:
1495  self.log.info("Masking NAN value pixels.")
1496  self.maskNanmaskNan(ccdExposure)
1497 
1498  if self.config.doWidenSaturationTrails:
1499  self.log.info("Widening saturation trails.")
1500  isrFunctions.widenSaturationTrails(ccdExposure.getMaskedImage().getMask())
1501 
1502  if self.config.doCameraSpecificMasking:
1503  self.log.info("Masking regions for camera specific reasons.")
1504  self.masking.run(ccdExposure)
1505 
1506  if self.config.doBrighterFatter:
1507  # We need to apply flats and darks before we can interpolate, and we
1508  # need to interpolate before we do B-F, but we do B-F without the
1509  # flats and darks applied so we can work in units of electrons or holes.
1510  # This context manager applies and then removes the darks and flats.
1511  #
1512  # We also do not want to interpolate values here, so operate on temporary
1513  # images so we can apply only the BF-correction and roll back the
1514  # interpolation.
1515  interpExp = ccdExposure.clone()
1516  with self.flatContextflatContext(interpExp, flat, dark):
1517  isrFunctions.interpolateFromMask(
1518  maskedImage=interpExp.getMaskedImage(),
1519  fwhm=self.config.fwhm,
1520  growSaturatedFootprints=self.config.growSaturationFootprintSize,
1521  maskNameList=self.config.maskListToInterpolate
1522  )
1523  bfExp = interpExp.clone()
1524 
1525  self.log.info("Applying brighter fatter correction using kernel type %s / gains %s.",
1526  type(bfKernel), type(bfGains))
1527  bfResults = isrFunctions.brighterFatterCorrection(bfExp, bfKernel,
1528  self.config.brighterFatterMaxIter,
1529  self.config.brighterFatterThreshold,
1530  self.config.brighterFatterApplyGain,
1531  bfGains)
1532  if bfResults[1] == self.config.brighterFatterMaxIter:
1533  self.log.warn("Brighter fatter correction did not converge, final difference %f.",
1534  bfResults[0])
1535  else:
1536  self.log.info("Finished brighter fatter correction in %d iterations.",
1537  bfResults[1])
1538  image = ccdExposure.getMaskedImage().getImage()
1539  bfCorr = bfExp.getMaskedImage().getImage()
1540  bfCorr -= interpExp.getMaskedImage().getImage()
1541  image += bfCorr
1542 
1543  # Applying the brighter-fatter correction applies a
1544  # convolution to the science image. At the edges this
1545  # convolution may not have sufficient valid pixels to
1546  # produce a valid correction. Mark pixels within the size
1547  # of the brighter-fatter kernel as EDGE to warn of this
1548  # fact.
1549  self.log.info("Ensuring image edges are masked as SUSPECT to the brighter-fatter kernel size.")
1550  self.maskEdgesmaskEdges(ccdExposure, numEdgePixels=numpy.max(bfKernel.shape) // 2,
1551  maskPlane="EDGE")
1552 
1553  if self.config.brighterFatterMaskGrowSize > 0:
1554  self.log.info("Growing masks to account for brighter-fatter kernel convolution.")
1555  for maskPlane in self.config.maskListToInterpolate:
1556  isrFunctions.growMasks(ccdExposure.getMask(),
1557  radius=self.config.brighterFatterMaskGrowSize,
1558  maskNameList=maskPlane,
1559  maskValue=maskPlane)
1560 
1561  self.debugViewdebugView(ccdExposure, "doBrighterFatter")
1562 
1563  if self.config.doDark:
1564  self.log.info("Applying dark correction.")
1565  self.darkCorrectiondarkCorrection(ccdExposure, dark)
1566  self.debugViewdebugView(ccdExposure, "doDark")
1567 
1568  if self.config.doFringe and not self.config.fringeAfterFlat:
1569  self.log.info("Applying fringe correction before flat.")
1570  self.fringe.run(ccdExposure, **fringes.getDict())
1571  self.debugViewdebugView(ccdExposure, "doFringe")
1572 
1573  if self.config.doStrayLight and self.strayLight.check(ccdExposure):
1574  self.log.info("Checking strayLight correction.")
1575  self.strayLight.run(ccdExposure, strayLightData)
1576  self.debugViewdebugView(ccdExposure, "doStrayLight")
1577 
1578  if self.config.doFlat:
1579  self.log.info("Applying flat correction.")
1580  self.flatCorrectionflatCorrection(ccdExposure, flat)
1581  self.debugViewdebugView(ccdExposure, "doFlat")
1582 
1583  if self.config.doApplyGains:
1584  self.log.info("Applying gain correction instead of flat.")
1585  isrFunctions.applyGains(ccdExposure, self.config.normalizeGains)
1586 
1587  if self.config.doFringe and self.config.fringeAfterFlat:
1588  self.log.info("Applying fringe correction after flat.")
1589  self.fringe.run(ccdExposure, **fringes.getDict())
1590 
1591  if self.config.doVignette:
1592  self.log.info("Constructing Vignette polygon.")
1593  self.vignettePolygonvignettePolygon = self.vignette.run(ccdExposure)
1594 
1595  if self.config.vignette.doWriteVignettePolygon:
1596  self.setValidPolygonIntersectsetValidPolygonIntersect(ccdExposure, self.vignettePolygonvignettePolygon)
1597 
1598  if self.config.doAttachTransmissionCurve:
1599  self.log.info("Adding transmission curves.")
1600  isrFunctions.attachTransmissionCurve(ccdExposure, opticsTransmission=opticsTransmission,
1601  filterTransmission=filterTransmission,
1602  sensorTransmission=sensorTransmission,
1603  atmosphereTransmission=atmosphereTransmission)
1604 
1605  flattenedThumb = None
1606  if self.config.qa.doThumbnailFlattened:
1607  flattenedThumb = isrQa.makeThumbnail(ccdExposure, isrQaConfig=self.config.qa)
1608 
1609  if self.config.doIlluminationCorrection and filterName in self.config.illumFilters:
1610  self.log.info("Performing illumination correction.")
1611  isrFunctions.illuminationCorrection(ccdExposure.getMaskedImage(),
1612  illumMaskedImage, illumScale=self.config.illumScale,
1613  trimToFit=self.config.doTrimToMatchCalib)
1614 
1615  preInterpExp = None
1616  if self.config.doSaveInterpPixels:
1617  preInterpExp = ccdExposure.clone()
1618 
1619  # Reset and interpolate bad pixels.
1620  #
1621  # Large contiguous bad regions (which should have the BAD mask
1622  # bit set) should have their values set to the image median.
1623  # This group should include defects and bad amplifiers. As the
1624  # area covered by these defects are large, there's little
1625  # reason to expect that interpolation would provide a more
1626  # useful value.
1627  #
1628  # Smaller defects can be safely interpolated after the larger
1629  # regions have had their pixel values reset. This ensures
1630  # that the remaining defects adjacent to bad amplifiers (as an
1631  # example) do not attempt to interpolate extreme values.
1632  if self.config.doSetBadRegions:
1633  badPixelCount, badPixelValue = isrFunctions.setBadRegions(ccdExposure)
1634  if badPixelCount > 0:
1635  self.log.info("Set %d BAD pixels to %f.", badPixelCount, badPixelValue)
1636 
1637  if self.config.doInterpolate:
1638  self.log.info("Interpolating masked pixels.")
1639  isrFunctions.interpolateFromMask(
1640  maskedImage=ccdExposure.getMaskedImage(),
1641  fwhm=self.config.fwhm,
1642  growSaturatedFootprints=self.config.growSaturationFootprintSize,
1643  maskNameList=list(self.config.maskListToInterpolate)
1644  )
1645 
1646  self.roughZeroPointroughZeroPoint(ccdExposure)
1647 
1648  if self.config.doMeasureBackground:
1649  self.log.info("Measuring background level.")
1650  self.measureBackgroundmeasureBackground(ccdExposure, self.config.qa)
1651 
1652  if self.config.qa is not None and self.config.qa.saveStats is True:
1653  for amp in ccd:
1654  ampExposure = ccdExposure.Factory(ccdExposure, amp.getBBox())
1655  qaStats = afwMath.makeStatistics(ampExposure.getImage(),
1656  afwMath.MEDIAN | afwMath.STDEVCLIP)
1657  self.metadata.set("ISR BACKGROUND {} MEDIAN".format(amp.getName()),
1658  qaStats.getValue(afwMath.MEDIAN))
1659  self.metadata.set("ISR BACKGROUND {} STDEV".format(amp.getName()),
1660  qaStats.getValue(afwMath.STDEVCLIP))
1661  self.log.debug(" Background stats for amplifer %s: %f +/- %f",
1662  amp.getName(), qaStats.getValue(afwMath.MEDIAN),
1663  qaStats.getValue(afwMath.STDEVCLIP))
1664 
1665  self.debugViewdebugView(ccdExposure, "postISRCCD")
1666 
1667  return pipeBase.Struct(
1668  exposure=ccdExposure,
1669  ossThumb=ossThumb,
1670  flattenedThumb=flattenedThumb,
1671 
1672  preInterpolatedExposure=preInterpExp,
1673  outputExposure=ccdExposure,
1674  outputOssThumbnail=ossThumb,
1675  outputFlattenedThumbnail=flattenedThumb,
1676  )
1677 
1678  @pipeBase.timeMethod
1679  def runDataRef(self, sensorRef):
1680  """Perform instrument signature removal on a ButlerDataRef of a Sensor.
1681 
1682  This method contains the `CmdLineTask` interface to the ISR
1683  processing. All IO is handled here, freeing the `run()` method
1684  to manage only pixel-level calculations. The steps performed
1685  are:
1686  - Read in necessary detrending/isr/calibration data.
1687  - Process raw exposure in `run()`.
1688  - Persist the ISR-corrected exposure as "postISRCCD" if
1689  config.doWrite=True.
1690 
1691  Parameters
1692  ----------
1693  sensorRef : `daf.persistence.butlerSubset.ButlerDataRef`
1694  DataRef of the detector data to be processed
1695 
1696  Returns
1697  -------
1698  result : `lsst.pipe.base.Struct`
1699  Result struct with component:
1700  - ``exposure`` : `afw.image.Exposure`
1701  The fully ISR corrected exposure.
1702 
1703  Raises
1704  ------
1705  RuntimeError
1706  Raised if a configuration option is set to True, but the
1707  required calibration data does not exist.
1708 
1709  """
1710  self.log.info("Performing ISR on sensor %s.", sensorRef.dataId)
1711 
1712  ccdExposure = sensorRef.get(self.config.datasetType)
1713 
1714  camera = sensorRef.get("camera")
1715  isrData = self.readIsrDatareadIsrData(sensorRef, ccdExposure)
1716 
1717  result = self.runrun(ccdExposure, camera=camera, **isrData.getDict())
1718 
1719  if self.config.doWrite:
1720  sensorRef.put(result.exposure, "postISRCCD")
1721  if result.preInterpolatedExposure is not None:
1722  sensorRef.put(result.preInterpolatedExposure, "postISRCCD_uninterpolated")
1723  if result.ossThumb is not None:
1724  isrQa.writeThumbnail(sensorRef, result.ossThumb, "ossThumb")
1725  if result.flattenedThumb is not None:
1726  isrQa.writeThumbnail(sensorRef, result.flattenedThumb, "flattenedThumb")
1727 
1728  return result
1729 
1730  def getIsrExposure(self, dataRef, datasetType, dateObs=None, immediate=True):
1731  """Retrieve a calibration dataset for removing instrument signature.
1732 
1733  Parameters
1734  ----------
1735 
1736  dataRef : `daf.persistence.butlerSubset.ButlerDataRef`
1737  DataRef of the detector data to find calibration datasets
1738  for.
1739  datasetType : `str`
1740  Type of dataset to retrieve (e.g. 'bias', 'flat', etc).
1741  dateObs : `str`, optional
1742  Date of the observation. Used to correct butler failures
1743  when using fallback filters.
1744  immediate : `Bool`
1745  If True, disable butler proxies to enable error handling
1746  within this routine.
1747 
1748  Returns
1749  -------
1750  exposure : `lsst.afw.image.Exposure`
1751  Requested calibration frame.
1752 
1753  Raises
1754  ------
1755  RuntimeError
1756  Raised if no matching calibration frame can be found.
1757  """
1758  try:
1759  exp = dataRef.get(datasetType, immediate=immediate)
1760  except Exception as exc1:
1761  if not self.config.fallbackFilterName:
1762  raise RuntimeError("Unable to retrieve %s for %s: %s." % (datasetType, dataRef.dataId, exc1))
1763  try:
1764  if self.config.useFallbackDate and dateObs:
1765  exp = dataRef.get(datasetType, filter=self.config.fallbackFilterName,
1766  dateObs=dateObs, immediate=immediate)
1767  else:
1768  exp = dataRef.get(datasetType, filter=self.config.fallbackFilterName, immediate=immediate)
1769  except Exception as exc2:
1770  raise RuntimeError("Unable to retrieve %s for %s, even with fallback filter %s: %s AND %s." %
1771  (datasetType, dataRef.dataId, self.config.fallbackFilterName, exc1, exc2))
1772  self.log.warn("Using fallback calibration from filter %s.", self.config.fallbackFilterName)
1773 
1774  if self.config.doAssembleIsrExposures:
1775  exp = self.assembleCcd.assembleCcd(exp)
1776  return exp
1777 
1778  def ensureExposure(self, inputExp, camera, detectorNum):
1779  """Ensure that the data returned by Butler is a fully constructed exposure.
1780 
1781  ISR requires exposure-level image data for historical reasons, so if we did
1782  not recieve that from Butler, construct it from what we have, modifying the
1783  input in place.
1784 
1785  Parameters
1786  ----------
1787  inputExp : `lsst.afw.image.Exposure`, `lsst.afw.image.DecoratedImageU`, or
1788  `lsst.afw.image.ImageF`
1789  The input data structure obtained from Butler.
1790  camera : `lsst.afw.cameraGeom.camera`
1791  The camera associated with the image. Used to find the appropriate
1792  detector.
1793  detectorNum : `int`
1794  The detector this exposure should match.
1795 
1796  Returns
1797  -------
1798  inputExp : `lsst.afw.image.Exposure`
1799  The re-constructed exposure, with appropriate detector parameters.
1800 
1801  Raises
1802  ------
1803  TypeError
1804  Raised if the input data cannot be used to construct an exposure.
1805  """
1806  if isinstance(inputExp, afwImage.DecoratedImageU):
1807  inputExp = afwImage.makeExposure(afwImage.makeMaskedImage(inputExp))
1808  elif isinstance(inputExp, afwImage.ImageF):
1809  inputExp = afwImage.makeExposure(afwImage.makeMaskedImage(inputExp))
1810  elif isinstance(inputExp, afwImage.MaskedImageF):
1811  inputExp = afwImage.makeExposure(inputExp)
1812  elif isinstance(inputExp, afwImage.Exposure):
1813  pass
1814  elif inputExp is None:
1815  # Assume this will be caught by the setup if it is a problem.
1816  return inputExp
1817  else:
1818  raise TypeError("Input Exposure is not known type in isrTask.ensureExposure: %s." %
1819  (type(inputExp), ))
1820 
1821  if inputExp.getDetector() is None:
1822  inputExp.setDetector(camera[detectorNum])
1823 
1824  return inputExp
1825 
1826  def convertIntToFloat(self, exposure):
1827  """Convert exposure image from uint16 to float.
1828 
1829  If the exposure does not need to be converted, the input is
1830  immediately returned. For exposures that are converted to use
1831  floating point pixels, the variance is set to unity and the
1832  mask to zero.
1833 
1834  Parameters
1835  ----------
1836  exposure : `lsst.afw.image.Exposure`
1837  The raw exposure to be converted.
1838 
1839  Returns
1840  -------
1841  newexposure : `lsst.afw.image.Exposure`
1842  The input ``exposure``, converted to floating point pixels.
1843 
1844  Raises
1845  ------
1846  RuntimeError
1847  Raised if the exposure type cannot be converted to float.
1848 
1849  """
1850  if isinstance(exposure, afwImage.ExposureF):
1851  # Nothing to be done
1852  self.log.debug("Exposure already of type float.")
1853  return exposure
1854  if not hasattr(exposure, "convertF"):
1855  raise RuntimeError("Unable to convert exposure (%s) to float." % type(exposure))
1856 
1857  newexposure = exposure.convertF()
1858  newexposure.variance[:] = 1
1859  newexposure.mask[:] = 0x0
1860 
1861  return newexposure
1862 
1863  def maskAmplifier(self, ccdExposure, amp, defects):
1864  """Identify bad amplifiers, saturated and suspect pixels.
1865 
1866  Parameters
1867  ----------
1868  ccdExposure : `lsst.afw.image.Exposure`
1869  Input exposure to be masked.
1870  amp : `lsst.afw.table.AmpInfoCatalog`
1871  Catalog of parameters defining the amplifier on this
1872  exposure to mask.
1873  defects : `lsst.ip.isr.Defects`
1874  List of defects. Used to determine if the entire
1875  amplifier is bad.
1876 
1877  Returns
1878  -------
1879  badAmp : `Bool`
1880  If this is true, the entire amplifier area is covered by
1881  defects and unusable.
1882 
1883  """
1884  maskedImage = ccdExposure.getMaskedImage()
1885 
1886  badAmp = False
1887 
1888  # Check if entire amp region is defined as a defect (need to use amp.getBBox() for correct
1889  # comparison with current defects definition.
1890  if defects is not None:
1891  badAmp = bool(sum([v.getBBox().contains(amp.getBBox()) for v in defects]))
1892 
1893  # In the case of a bad amp, we will set mask to "BAD" (here use amp.getRawBBox() for correct
1894  # association with pixels in current ccdExposure).
1895  if badAmp:
1896  dataView = afwImage.MaskedImageF(maskedImage, amp.getRawBBox(),
1897  afwImage.PARENT)
1898  maskView = dataView.getMask()
1899  maskView |= maskView.getPlaneBitMask("BAD")
1900  del maskView
1901  return badAmp
1902 
1903  # Mask remaining defects after assembleCcd() to allow for defects that cross amplifier boundaries.
1904  # Saturation and suspect pixels can be masked now, though.
1905  limits = dict()
1906  if self.config.doSaturation and not badAmp:
1907  limits.update({self.config.saturatedMaskName: amp.getSaturation()})
1908  if self.config.doSuspect and not badAmp:
1909  limits.update({self.config.suspectMaskName: amp.getSuspectLevel()})
1910  if math.isfinite(self.config.saturation):
1911  limits.update({self.config.saturatedMaskName: self.config.saturation})
1912 
1913  for maskName, maskThreshold in limits.items():
1914  if not math.isnan(maskThreshold):
1915  dataView = maskedImage.Factory(maskedImage, amp.getRawBBox())
1916  isrFunctions.makeThresholdMask(
1917  maskedImage=dataView,
1918  threshold=maskThreshold,
1919  growFootprints=0,
1920  maskName=maskName
1921  )
1922 
1923  # Determine if we've fully masked this amplifier with SUSPECT and SAT pixels.
1924  maskView = afwImage.Mask(maskedImage.getMask(), amp.getRawDataBBox(),
1925  afwImage.PARENT)
1926  maskVal = maskView.getPlaneBitMask([self.config.saturatedMaskName,
1927  self.config.suspectMaskName])
1928  if numpy.all(maskView.getArray() & maskVal > 0):
1929  badAmp = True
1930  maskView |= maskView.getPlaneBitMask("BAD")
1931 
1932  return badAmp
1933 
1934  def overscanCorrection(self, ccdExposure, amp):
1935  """Apply overscan correction in place.
1936 
1937  This method does initial pixel rejection of the overscan
1938  region. The overscan can also be optionally segmented to
1939  allow for discontinuous overscan responses to be fit
1940  separately. The actual overscan subtraction is performed by
1941  the `lsst.ip.isr.isrFunctions.overscanCorrection` function,
1942  which is called here after the amplifier is preprocessed.
1943 
1944  Parameters
1945  ----------
1946  ccdExposure : `lsst.afw.image.Exposure`
1947  Exposure to have overscan correction performed.
1948  amp : `lsst.afw.cameraGeom.Amplifer`
1949  The amplifier to consider while correcting the overscan.
1950 
1951  Returns
1952  -------
1953  overscanResults : `lsst.pipe.base.Struct`
1954  Result struct with components:
1955  - ``imageFit`` : scalar or `lsst.afw.image.Image`
1956  Value or fit subtracted from the amplifier image data.
1957  - ``overscanFit`` : scalar or `lsst.afw.image.Image`
1958  Value or fit subtracted from the overscan image data.
1959  - ``overscanImage`` : `lsst.afw.image.Image`
1960  Image of the overscan region with the overscan
1961  correction applied. This quantity is used to estimate
1962  the amplifier read noise empirically.
1963 
1964  Raises
1965  ------
1966  RuntimeError
1967  Raised if the ``amp`` does not contain raw pixel information.
1968 
1969  See Also
1970  --------
1971  lsst.ip.isr.isrFunctions.overscanCorrection
1972  """
1973  if amp.getRawHorizontalOverscanBBox().isEmpty():
1974  self.log.info("ISR_OSCAN: No overscan region. Not performing overscan correction.")
1975  return None
1976 
1977  statControl = afwMath.StatisticsControl()
1978  statControl.setAndMask(ccdExposure.mask.getPlaneBitMask("SAT"))
1979 
1980  # Determine the bounding boxes
1981  dataBBox = amp.getRawDataBBox()
1982  oscanBBox = amp.getRawHorizontalOverscanBBox()
1983  dx0 = 0
1984  dx1 = 0
1985 
1986  prescanBBox = amp.getRawPrescanBBox()
1987  if (oscanBBox.getBeginX() > prescanBBox.getBeginX()): # amp is at the right
1988  dx0 += self.config.overscanNumLeadingColumnsToSkip
1989  dx1 -= self.config.overscanNumTrailingColumnsToSkip
1990  else:
1991  dx0 += self.config.overscanNumTrailingColumnsToSkip
1992  dx1 -= self.config.overscanNumLeadingColumnsToSkip
1993 
1994  # Determine if we need to work on subregions of the amplifier and overscan.
1995  imageBBoxes = []
1996  overscanBBoxes = []
1997 
1998  if ((self.config.overscanBiasJump
1999  and self.config.overscanBiasJumpLocation)
2000  and (ccdExposure.getMetadata().exists(self.config.overscanBiasJumpKeyword)
2001  and ccdExposure.getMetadata().getScalar(self.config.overscanBiasJumpKeyword) in
2002  self.config.overscanBiasJumpDevices)):
2003  if amp.getReadoutCorner() in (ReadoutCorner.LL, ReadoutCorner.LR):
2004  yLower = self.config.overscanBiasJumpLocation
2005  yUpper = dataBBox.getHeight() - yLower
2006  else:
2007  yUpper = self.config.overscanBiasJumpLocation
2008  yLower = dataBBox.getHeight() - yUpper
2009 
2010  imageBBoxes.append(lsst.geom.Box2I(dataBBox.getBegin(),
2011  lsst.geom.Extent2I(dataBBox.getWidth(), yLower)))
2012  overscanBBoxes.append(lsst.geom.Box2I(oscanBBox.getBegin() + lsst.geom.Extent2I(dx0, 0),
2013  lsst.geom.Extent2I(oscanBBox.getWidth() - dx0 + dx1,
2014  yLower)))
2015 
2016  imageBBoxes.append(lsst.geom.Box2I(dataBBox.getBegin() + lsst.geom.Extent2I(0, yLower),
2017  lsst.geom.Extent2I(dataBBox.getWidth(), yUpper)))
2018  overscanBBoxes.append(lsst.geom.Box2I(oscanBBox.getBegin() + lsst.geom.Extent2I(dx0, yLower),
2019  lsst.geom.Extent2I(oscanBBox.getWidth() - dx0 + dx1,
2020  yUpper)))
2021  else:
2022  imageBBoxes.append(lsst.geom.Box2I(dataBBox.getBegin(),
2023  lsst.geom.Extent2I(dataBBox.getWidth(), dataBBox.getHeight())))
2024  overscanBBoxes.append(lsst.geom.Box2I(oscanBBox.getBegin() + lsst.geom.Extent2I(dx0, 0),
2025  lsst.geom.Extent2I(oscanBBox.getWidth() - dx0 + dx1,
2026  oscanBBox.getHeight())))
2027 
2028  # Perform overscan correction on subregions, ensuring saturated pixels are masked.
2029  for imageBBox, overscanBBox in zip(imageBBoxes, overscanBBoxes):
2030  ampImage = ccdExposure.maskedImage[imageBBox]
2031  overscanImage = ccdExposure.maskedImage[overscanBBox]
2032 
2033  overscanArray = overscanImage.image.array
2034  median = numpy.ma.median(numpy.ma.masked_where(overscanImage.mask.array, overscanArray))
2035  bad = numpy.where(numpy.abs(overscanArray - median) > self.config.overscanMaxDev)
2036  overscanImage.mask.array[bad] = overscanImage.mask.getPlaneBitMask("SAT")
2037 
2038  statControl = afwMath.StatisticsControl()
2039  statControl.setAndMask(ccdExposure.mask.getPlaneBitMask("SAT"))
2040 
2041  overscanResults = self.overscan.run(ampImage.getImage(), overscanImage, amp)
2042 
2043  # Measure average overscan levels and record them in the metadata.
2044  levelStat = afwMath.MEDIAN
2045  sigmaStat = afwMath.STDEVCLIP
2046 
2047  sctrl = afwMath.StatisticsControl(self.config.qa.flatness.clipSigma,
2048  self.config.qa.flatness.nIter)
2049  metadata = ccdExposure.getMetadata()
2050  ampNum = amp.getName()
2051  # if self.config.overscanFitType in ("MEDIAN", "MEAN", "MEANCLIP"):
2052  if isinstance(overscanResults.overscanFit, float):
2053  metadata.set("ISR_OSCAN_LEVEL%s" % ampNum, overscanResults.overscanFit)
2054  metadata.set("ISR_OSCAN_SIGMA%s" % ampNum, 0.0)
2055  else:
2056  stats = afwMath.makeStatistics(overscanResults.overscanFit, levelStat | sigmaStat, sctrl)
2057  metadata.set("ISR_OSCAN_LEVEL%s" % ampNum, stats.getValue(levelStat))
2058  metadata.set("ISR_OSCAN_SIGMA%s" % ampNum, stats.getValue(sigmaStat))
2059 
2060  return overscanResults
2061 
2062  def updateVariance(self, ampExposure, amp, overscanImage=None):
2063  """Set the variance plane using the amplifier gain and read noise
2064 
2065  The read noise is calculated from the ``overscanImage`` if the
2066  ``doEmpiricalReadNoise`` option is set in the configuration; otherwise
2067  the value from the amplifier data is used.
2068 
2069  Parameters
2070  ----------
2071  ampExposure : `lsst.afw.image.Exposure`
2072  Exposure to process.
2073  amp : `lsst.afw.table.AmpInfoRecord` or `FakeAmp`
2074  Amplifier detector data.
2075  overscanImage : `lsst.afw.image.MaskedImage`, optional.
2076  Image of overscan, required only for empirical read noise.
2077 
2078  See also
2079  --------
2080  lsst.ip.isr.isrFunctions.updateVariance
2081  """
2082  maskPlanes = [self.config.saturatedMaskName, self.config.suspectMaskName]
2083  gain = amp.getGain()
2084 
2085  if math.isnan(gain):
2086  gain = 1.0
2087  self.log.warn("Gain set to NAN! Updating to 1.0 to generate Poisson variance.")
2088  elif gain <= 0:
2089  patchedGain = 1.0
2090  self.log.warn("Gain for amp %s == %g <= 0; setting to %f.",
2091  amp.getName(), gain, patchedGain)
2092  gain = patchedGain
2093 
2094  if self.config.doEmpiricalReadNoise and overscanImage is None:
2095  self.log.info("Overscan is none for EmpiricalReadNoise.")
2096 
2097  if self.config.doEmpiricalReadNoise and overscanImage is not None:
2098  stats = afwMath.StatisticsControl()
2099  stats.setAndMask(overscanImage.mask.getPlaneBitMask(maskPlanes))
2100  readNoise = afwMath.makeStatistics(overscanImage, afwMath.STDEVCLIP, stats).getValue()
2101  self.log.info("Calculated empirical read noise for amp %s: %f.",
2102  amp.getName(), readNoise)
2103  else:
2104  readNoise = amp.getReadNoise()
2105 
2106  isrFunctions.updateVariance(
2107  maskedImage=ampExposure.getMaskedImage(),
2108  gain=gain,
2109  readNoise=readNoise,
2110  )
2111 
2112  def darkCorrection(self, exposure, darkExposure, invert=False):
2113  """Apply dark correction in place.
2114 
2115  Parameters
2116  ----------
2117  exposure : `lsst.afw.image.Exposure`
2118  Exposure to process.
2119  darkExposure : `lsst.afw.image.Exposure`
2120  Dark exposure of the same size as ``exposure``.
2121  invert : `Bool`, optional
2122  If True, re-add the dark to an already corrected image.
2123 
2124  Raises
2125  ------
2126  RuntimeError
2127  Raised if either ``exposure`` or ``darkExposure`` do not
2128  have their dark time defined.
2129 
2130  See Also
2131  --------
2132  lsst.ip.isr.isrFunctions.darkCorrection
2133  """
2134  expScale = exposure.getInfo().getVisitInfo().getDarkTime()
2135  if math.isnan(expScale):
2136  raise RuntimeError("Exposure darktime is NAN.")
2137  if darkExposure.getInfo().getVisitInfo() is not None \
2138  and not math.isnan(darkExposure.getInfo().getVisitInfo().getDarkTime()):
2139  darkScale = darkExposure.getInfo().getVisitInfo().getDarkTime()
2140  else:
2141  # DM-17444: darkExposure.getInfo.getVisitInfo() is None
2142  # so getDarkTime() does not exist.
2143  self.log.warn("darkExposure.getInfo().getVisitInfo() does not exist. Using darkScale = 1.0.")
2144  darkScale = 1.0
2145 
2146  isrFunctions.darkCorrection(
2147  maskedImage=exposure.getMaskedImage(),
2148  darkMaskedImage=darkExposure.getMaskedImage(),
2149  expScale=expScale,
2150  darkScale=darkScale,
2151  invert=invert,
2152  trimToFit=self.config.doTrimToMatchCalib
2153  )
2154 
2155  def doLinearize(self, detector):
2156  """Check if linearization is needed for the detector cameraGeom.
2157 
2158  Checks config.doLinearize and the linearity type of the first
2159  amplifier.
2160 
2161  Parameters
2162  ----------
2163  detector : `lsst.afw.cameraGeom.Detector`
2164  Detector to get linearity type from.
2165 
2166  Returns
2167  -------
2168  doLinearize : `Bool`
2169  If True, linearization should be performed.
2170  """
2171  return self.config.doLinearize and \
2172  detector.getAmplifiers()[0].getLinearityType() != NullLinearityType
2173 
2174  def flatCorrection(self, exposure, flatExposure, invert=False):
2175  """Apply flat correction in place.
2176 
2177  Parameters
2178  ----------
2179  exposure : `lsst.afw.image.Exposure`
2180  Exposure to process.
2181  flatExposure : `lsst.afw.image.Exposure`
2182  Flat exposure of the same size as ``exposure``.
2183  invert : `Bool`, optional
2184  If True, unflatten an already flattened image.
2185 
2186  See Also
2187  --------
2188  lsst.ip.isr.isrFunctions.flatCorrection
2189  """
2190  isrFunctions.flatCorrection(
2191  maskedImage=exposure.getMaskedImage(),
2192  flatMaskedImage=flatExposure.getMaskedImage(),
2193  scalingType=self.config.flatScalingType,
2194  userScale=self.config.flatUserScale,
2195  invert=invert,
2196  trimToFit=self.config.doTrimToMatchCalib
2197  )
2198 
2199  def saturationDetection(self, exposure, amp):
2200  """Detect saturated pixels and mask them using mask plane config.saturatedMaskName, in place.
2201 
2202  Parameters
2203  ----------
2204  exposure : `lsst.afw.image.Exposure`
2205  Exposure to process. Only the amplifier DataSec is processed.
2206  amp : `lsst.afw.table.AmpInfoCatalog`
2207  Amplifier detector data.
2208 
2209  See Also
2210  --------
2211  lsst.ip.isr.isrFunctions.makeThresholdMask
2212  """
2213  if not math.isnan(amp.getSaturation()):
2214  maskedImage = exposure.getMaskedImage()
2215  dataView = maskedImage.Factory(maskedImage, amp.getRawBBox())
2216  isrFunctions.makeThresholdMask(
2217  maskedImage=dataView,
2218  threshold=amp.getSaturation(),
2219  growFootprints=0,
2220  maskName=self.config.saturatedMaskName,
2221  )
2222 
2223  def saturationInterpolation(self, exposure):
2224  """Interpolate over saturated pixels, in place.
2225 
2226  This method should be called after `saturationDetection`, to
2227  ensure that the saturated pixels have been identified in the
2228  SAT mask. It should also be called after `assembleCcd`, since
2229  saturated regions may cross amplifier boundaries.
2230 
2231  Parameters
2232  ----------
2233  exposure : `lsst.afw.image.Exposure`
2234  Exposure to process.
2235 
2236  See Also
2237  --------
2238  lsst.ip.isr.isrTask.saturationDetection
2239  lsst.ip.isr.isrFunctions.interpolateFromMask
2240  """
2241  isrFunctions.interpolateFromMask(
2242  maskedImage=exposure.getMaskedImage(),
2243  fwhm=self.config.fwhm,
2244  growSaturatedFootprints=self.config.growSaturationFootprintSize,
2245  maskNameList=list(self.config.saturatedMaskName),
2246  )
2247 
2248  def suspectDetection(self, exposure, amp):
2249  """Detect suspect pixels and mask them using mask plane config.suspectMaskName, in place.
2250 
2251  Parameters
2252  ----------
2253  exposure : `lsst.afw.image.Exposure`
2254  Exposure to process. Only the amplifier DataSec is processed.
2255  amp : `lsst.afw.table.AmpInfoCatalog`
2256  Amplifier detector data.
2257 
2258  See Also
2259  --------
2260  lsst.ip.isr.isrFunctions.makeThresholdMask
2261 
2262  Notes
2263  -----
2264  Suspect pixels are pixels whose value is greater than amp.getSuspectLevel().
2265  This is intended to indicate pixels that may be affected by unknown systematics;
2266  for example if non-linearity corrections above a certain level are unstable
2267  then that would be a useful value for suspectLevel. A value of `nan` indicates
2268  that no such level exists and no pixels are to be masked as suspicious.
2269  """
2270  suspectLevel = amp.getSuspectLevel()
2271  if math.isnan(suspectLevel):
2272  return
2273 
2274  maskedImage = exposure.getMaskedImage()
2275  dataView = maskedImage.Factory(maskedImage, amp.getRawBBox())
2276  isrFunctions.makeThresholdMask(
2277  maskedImage=dataView,
2278  threshold=suspectLevel,
2279  growFootprints=0,
2280  maskName=self.config.suspectMaskName,
2281  )
2282 
2283  def maskDefect(self, exposure, defectBaseList):
2284  """Mask defects using mask plane "BAD", in place.
2285 
2286  Parameters
2287  ----------
2288  exposure : `lsst.afw.image.Exposure`
2289  Exposure to process.
2290  defectBaseList : `lsst.ip.isr.Defects` or `list` of
2291  `lsst.afw.image.DefectBase`.
2292  List of defects to mask.
2293 
2294  Notes
2295  -----
2296  Call this after CCD assembly, since defects may cross amplifier boundaries.
2297  """
2298  maskedImage = exposure.getMaskedImage()
2299  if not isinstance(defectBaseList, Defects):
2300  # Promotes DefectBase to Defect
2301  defectList = Defects(defectBaseList)
2302  else:
2303  defectList = defectBaseList
2304  defectList.maskPixels(maskedImage, maskName="BAD")
2305 
2306  def maskEdges(self, exposure, numEdgePixels=0, maskPlane="SUSPECT", level='DETECTOR'):
2307  """Mask edge pixels with applicable mask plane.
2308 
2309  Parameters
2310  ----------
2311  exposure : `lsst.afw.image.Exposure`
2312  Exposure to process.
2313  numEdgePixels : `int`, optional
2314  Number of edge pixels to mask.
2315  maskPlane : `str`, optional
2316  Mask plane name to use.
2317  level : `str`, optional
2318  Level at which to mask edges.
2319  """
2320  maskedImage = exposure.getMaskedImage()
2321  maskBitMask = maskedImage.getMask().getPlaneBitMask(maskPlane)
2322 
2323  if numEdgePixels > 0:
2324  if level == 'DETECTOR':
2325  boxes = [maskedImage.getBBox()]
2326  elif level == 'AMP':
2327  boxes = [amp.getBBox() for amp in exposure.getDetector()]
2328 
2329  for box in boxes:
2330  # This makes a bbox numEdgeSuspect pixels smaller than the image on each side
2331  subImage = maskedImage[box]
2332  box.grow(-numEdgePixels)
2333  # Mask pixels outside box
2334  SourceDetectionTask.setEdgeBits(
2335  subImage,
2336  box,
2337  maskBitMask)
2338 
2339  def maskAndInterpolateDefects(self, exposure, defectBaseList):
2340  """Mask and interpolate defects using mask plane "BAD", in place.
2341 
2342  Parameters
2343  ----------
2344  exposure : `lsst.afw.image.Exposure`
2345  Exposure to process.
2346  defectBaseList : `lsst.ip.isr.Defects` or `list` of
2347  `lsst.afw.image.DefectBase`.
2348  List of defects to mask and interpolate.
2349 
2350  See Also
2351  --------
2352  lsst.ip.isr.isrTask.maskDefect
2353  """
2354  self.maskDefectmaskDefect(exposure, defectBaseList)
2355  self.maskEdgesmaskEdges(exposure, numEdgePixels=self.config.numEdgeSuspect,
2356  maskPlane="SUSPECT", level=self.config.edgeMaskLevel)
2357  isrFunctions.interpolateFromMask(
2358  maskedImage=exposure.getMaskedImage(),
2359  fwhm=self.config.fwhm,
2360  growSaturatedFootprints=0,
2361  maskNameList=["BAD"],
2362  )
2363 
2364  def maskNan(self, exposure):
2365  """Mask NaNs using mask plane "UNMASKEDNAN", in place.
2366 
2367  Parameters
2368  ----------
2369  exposure : `lsst.afw.image.Exposure`
2370  Exposure to process.
2371 
2372  Notes
2373  -----
2374  We mask over all NaNs, including those that are masked with
2375  other bits (because those may or may not be interpolated over
2376  later, and we want to remove all NaNs). Despite this
2377  behaviour, the "UNMASKEDNAN" mask plane is used to preserve
2378  the historical name.
2379  """
2380  maskedImage = exposure.getMaskedImage()
2381 
2382  # Find and mask NaNs
2383  maskedImage.getMask().addMaskPlane("UNMASKEDNAN")
2384  maskVal = maskedImage.getMask().getPlaneBitMask("UNMASKEDNAN")
2385  numNans = maskNans(maskedImage, maskVal)
2386  self.metadata.set("NUMNANS", numNans)
2387  if numNans > 0:
2388  self.log.warn("There were %d unmasked NaNs.", numNans)
2389 
2390  def maskAndInterpolateNan(self, exposure):
2391  """"Mask and interpolate NaNs using mask plane "UNMASKEDNAN", in place.
2392 
2393  Parameters
2394  ----------
2395  exposure : `lsst.afw.image.Exposure`
2396  Exposure to process.
2397 
2398  See Also
2399  --------
2400  lsst.ip.isr.isrTask.maskNan
2401  """
2402  self.maskNanmaskNan(exposure)
2403  isrFunctions.interpolateFromMask(
2404  maskedImage=exposure.getMaskedImage(),
2405  fwhm=self.config.fwhm,
2406  growSaturatedFootprints=0,
2407  maskNameList=["UNMASKEDNAN"],
2408  )
2409 
2410  def measureBackground(self, exposure, IsrQaConfig=None):
2411  """Measure the image background in subgrids, for quality control purposes.
2412 
2413  Parameters
2414  ----------
2415  exposure : `lsst.afw.image.Exposure`
2416  Exposure to process.
2417  IsrQaConfig : `lsst.ip.isr.isrQa.IsrQaConfig`
2418  Configuration object containing parameters on which background
2419  statistics and subgrids to use.
2420  """
2421  if IsrQaConfig is not None:
2422  statsControl = afwMath.StatisticsControl(IsrQaConfig.flatness.clipSigma,
2423  IsrQaConfig.flatness.nIter)
2424  maskVal = exposure.getMaskedImage().getMask().getPlaneBitMask(["BAD", "SAT", "DETECTED"])
2425  statsControl.setAndMask(maskVal)
2426  maskedImage = exposure.getMaskedImage()
2427  stats = afwMath.makeStatistics(maskedImage, afwMath.MEDIAN | afwMath.STDEVCLIP, statsControl)
2428  skyLevel = stats.getValue(afwMath.MEDIAN)
2429  skySigma = stats.getValue(afwMath.STDEVCLIP)
2430  self.log.info("Flattened sky level: %f +/- %f.", skyLevel, skySigma)
2431  metadata = exposure.getMetadata()
2432  metadata.set('SKYLEVEL', skyLevel)
2433  metadata.set('SKYSIGMA', skySigma)
2434 
2435  # calcluating flatlevel over the subgrids
2436  stat = afwMath.MEANCLIP if IsrQaConfig.flatness.doClip else afwMath.MEAN
2437  meshXHalf = int(IsrQaConfig.flatness.meshX/2.)
2438  meshYHalf = int(IsrQaConfig.flatness.meshY/2.)
2439  nX = int((exposure.getWidth() + meshXHalf) / IsrQaConfig.flatness.meshX)
2440  nY = int((exposure.getHeight() + meshYHalf) / IsrQaConfig.flatness.meshY)
2441  skyLevels = numpy.zeros((nX, nY))
2442 
2443  for j in range(nY):
2444  yc = meshYHalf + j * IsrQaConfig.flatness.meshY
2445  for i in range(nX):
2446  xc = meshXHalf + i * IsrQaConfig.flatness.meshX
2447 
2448  xLLC = xc - meshXHalf
2449  yLLC = yc - meshYHalf
2450  xURC = xc + meshXHalf - 1
2451  yURC = yc + meshYHalf - 1
2452 
2453  bbox = lsst.geom.Box2I(lsst.geom.Point2I(xLLC, yLLC), lsst.geom.Point2I(xURC, yURC))
2454  miMesh = maskedImage.Factory(exposure.getMaskedImage(), bbox, afwImage.LOCAL)
2455 
2456  skyLevels[i, j] = afwMath.makeStatistics(miMesh, stat, statsControl).getValue()
2457 
2458  good = numpy.where(numpy.isfinite(skyLevels))
2459  skyMedian = numpy.median(skyLevels[good])
2460  flatness = (skyLevels[good] - skyMedian) / skyMedian
2461  flatness_rms = numpy.std(flatness)
2462  flatness_pp = flatness.max() - flatness.min() if len(flatness) > 0 else numpy.nan
2463 
2464  self.log.info("Measuring sky levels in %dx%d grids: %f.", nX, nY, skyMedian)
2465  self.log.info("Sky flatness in %dx%d grids - pp: %f rms: %f.",
2466  nX, nY, flatness_pp, flatness_rms)
2467 
2468  metadata.set('FLATNESS_PP', float(flatness_pp))
2469  metadata.set('FLATNESS_RMS', float(flatness_rms))
2470  metadata.set('FLATNESS_NGRIDS', '%dx%d' % (nX, nY))
2471  metadata.set('FLATNESS_MESHX', IsrQaConfig.flatness.meshX)
2472  metadata.set('FLATNESS_MESHY', IsrQaConfig.flatness.meshY)
2473 
2474  def roughZeroPoint(self, exposure):
2475  """Set an approximate magnitude zero point for the exposure.
2476 
2477  Parameters
2478  ----------
2479  exposure : `lsst.afw.image.Exposure`
2480  Exposure to process.
2481  """
2482  # TODO DM-28093: change this to: exposure.getFilterLabel().physicalLabel
2483  filterName = afwImage.Filter(exposure.getFilter().getId()).getName() # Canonical name for filter
2484  if filterName in self.config.fluxMag0T1:
2485  fluxMag0 = self.config.fluxMag0T1[filterName]
2486  else:
2487  self.log.warn("No rough magnitude zero point set for filter %s.", filterName)
2488  fluxMag0 = self.config.defaultFluxMag0T1
2489 
2490  expTime = exposure.getInfo().getVisitInfo().getExposureTime()
2491  if not expTime > 0: # handle NaN as well as <= 0
2492  self.log.warn("Non-positive exposure time; skipping rough zero point.")
2493  return
2494 
2495  self.log.info("Setting rough magnitude zero point: %f", 2.5*math.log10(fluxMag0*expTime))
2496  exposure.setPhotoCalib(afwImage.makePhotoCalibFromCalibZeroPoint(fluxMag0*expTime, 0.0))
2497 
2498  def setValidPolygonIntersect(self, ccdExposure, fpPolygon):
2499  """Set the valid polygon as the intersection of fpPolygon and the ccd corners.
2500 
2501  Parameters
2502  ----------
2503  ccdExposure : `lsst.afw.image.Exposure`
2504  Exposure to process.
2505  fpPolygon : `lsst.afw.geom.Polygon`
2506  Polygon in focal plane coordinates.
2507  """
2508  # Get ccd corners in focal plane coordinates
2509  ccd = ccdExposure.getDetector()
2510  fpCorners = ccd.getCorners(FOCAL_PLANE)
2511  ccdPolygon = Polygon(fpCorners)
2512 
2513  # Get intersection of ccd corners with fpPolygon
2514  intersect = ccdPolygon.intersectionSingle(fpPolygon)
2515 
2516  # Transform back to pixel positions and build new polygon
2517  ccdPoints = ccd.transform(intersect, FOCAL_PLANE, PIXELS)
2518  validPolygon = Polygon(ccdPoints)
2519  ccdExposure.getInfo().setValidPolygon(validPolygon)
2520 
2521  @contextmanager
2522  def flatContext(self, exp, flat, dark=None):
2523  """Context manager that applies and removes flats and darks,
2524  if the task is configured to apply them.
2525 
2526  Parameters
2527  ----------
2528  exp : `lsst.afw.image.Exposure`
2529  Exposure to process.
2530  flat : `lsst.afw.image.Exposure`
2531  Flat exposure the same size as ``exp``.
2532  dark : `lsst.afw.image.Exposure`, optional
2533  Dark exposure the same size as ``exp``.
2534 
2535  Yields
2536  ------
2537  exp : `lsst.afw.image.Exposure`
2538  The flat and dark corrected exposure.
2539  """
2540  if self.config.doDark and dark is not None:
2541  self.darkCorrectiondarkCorrection(exp, dark)
2542  if self.config.doFlat:
2543  self.flatCorrectionflatCorrection(exp, flat)
2544  try:
2545  yield exp
2546  finally:
2547  if self.config.doFlat:
2548  self.flatCorrectionflatCorrection(exp, flat, invert=True)
2549  if self.config.doDark and dark is not None:
2550  self.darkCorrectiondarkCorrection(exp, dark, invert=True)
2551 
2552  def debugView(self, exposure, stepname):
2553  """Utility function to examine ISR exposure at different stages.
2554 
2555  Parameters
2556  ----------
2557  exposure : `lsst.afw.image.Exposure`
2558  Exposure to view.
2559  stepname : `str`
2560  State of processing to view.
2561  """
2562  frame = getDebugFrame(self._display, stepname)
2563  if frame:
2564  display = getDisplay(frame)
2565  display.scale('asinh', 'zscale')
2566  display.mtv(exposure)
2567  prompt = "Press Enter to continue [c]... "
2568  while True:
2569  ans = input(prompt).lower()
2570  if ans in ("", "c",):
2571  break
2572 
2573 
2574 class FakeAmp(object):
2575  """A Detector-like object that supports returning gain and saturation level
2576 
2577  This is used when the input exposure does not have a detector.
2578 
2579  Parameters
2580  ----------
2581  exposure : `lsst.afw.image.Exposure`
2582  Exposure to generate a fake amplifier for.
2583  config : `lsst.ip.isr.isrTaskConfig`
2584  Configuration to apply to the fake amplifier.
2585  """
2586 
2587  def __init__(self, exposure, config):
2588  self._bbox_bbox = exposure.getBBox(afwImage.LOCAL)
2589  self._RawHorizontalOverscanBBox_RawHorizontalOverscanBBox = lsst.geom.Box2I()
2590  self._gain_gain = config.gain
2591  self._readNoise_readNoise = config.readNoise
2592  self._saturation_saturation = config.saturation
2593 
2594  def getBBox(self):
2595  return self._bbox_bbox
2596 
2597  def getRawBBox(self):
2598  return self._bbox_bbox
2599 
2601  return self._RawHorizontalOverscanBBox_RawHorizontalOverscanBBox
2602 
2603  def getGain(self):
2604  return self._gain_gain
2605 
2606  def getReadNoise(self):
2607  return self._readNoise_readNoise
2608 
2609  def getSaturation(self):
2610  return self._saturation_saturation
2611 
2612  def getSuspectLevel(self):
2613  return float("NaN")
2614 
2615 
2616 class RunIsrConfig(pexConfig.Config):
2617  isr = pexConfig.ConfigurableField(target=IsrTask, doc="Instrument signature removal")
2618 
2619 
2620 class RunIsrTask(pipeBase.CmdLineTask):
2621  """Task to wrap the default IsrTask to allow it to be retargeted.
2622 
2623  The standard IsrTask can be called directly from a command line
2624  program, but doing so removes the ability of the task to be
2625  retargeted. As most cameras override some set of the IsrTask
2626  methods, this would remove those data-specific methods in the
2627  output post-ISR images. This wrapping class fixes the issue,
2628  allowing identical post-ISR images to be generated by both the
2629  processCcd and isrTask code.
2630  """
2631  ConfigClass = RunIsrConfig
2632  _DefaultName = "runIsr"
2633 
2634  def __init__(self, *args, **kwargs):
2635  super().__init__(*args, **kwargs)
2636  self.makeSubtask("isr")
2637 
2638  def runDataRef(self, dataRef):
2639  """
2640  Parameters
2641  ----------
2642  dataRef : `lsst.daf.persistence.ButlerDataRef`
2643  data reference of the detector data to be processed
2644 
2645  Returns
2646  -------
2647  result : `pipeBase.Struct`
2648  Result struct with component:
2649 
2650  - exposure : `lsst.afw.image.Exposure`
2651  Post-ISR processed exposure.
2652  """
2653  return self.isr.runDataRef(dataRef)
def getRawHorizontalOverscanBBox(self)
Definition: isrTask.py:2600
def __init__(self, exposure, config)
Definition: isrTask.py:2587
def __init__(self, *config=None)
Definition: isrTask.py:258
def flatCorrection(self, exposure, flatExposure, invert=False)
Definition: isrTask.py:2174
def maskAndInterpolateNan(self, exposure)
Definition: isrTask.py:2390
def saturationInterpolation(self, exposure)
Definition: isrTask.py:2223
def runDataRef(self, sensorRef)
Definition: isrTask.py:1679
def maskNan(self, exposure)
Definition: isrTask.py:2364
def maskAmplifier(self, ccdExposure, amp, defects)
Definition: isrTask.py:1863
def debugView(self, exposure, stepname)
Definition: isrTask.py:2552
def getIsrExposure(self, dataRef, datasetType, dateObs=None, immediate=True)
Definition: isrTask.py:1730
def saturationDetection(self, exposure, amp)
Definition: isrTask.py:2199
def maskDefect(self, exposure, defectBaseList)
Definition: isrTask.py:2283
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:2306
def overscanCorrection(self, ccdExposure, amp)
Definition: isrTask.py:1934
def measureBackground(self, exposure, IsrQaConfig=None)
Definition: isrTask.py:2410
def roughZeroPoint(self, exposure)
Definition: isrTask.py:2474
def maskAndInterpolateDefects(self, exposure, defectBaseList)
Definition: isrTask.py:2339
def setValidPolygonIntersect(self, ccdExposure, fpPolygon)
Definition: isrTask.py:2498
def readIsrData(self, dataRef, rawExposure)
Definition: isrTask.py:1016
def ensureExposure(self, inputExp, camera, detectorNum)
Definition: isrTask.py:1778
def updateVariance(self, ampExposure, amp, overscanImage=None)
Definition: isrTask.py:2062
def doLinearize(self, detector)
Definition: isrTask.py:2155
def flatContext(self, exp, flat, dark=None)
Definition: isrTask.py:2522
def convertIntToFloat(self, exposure)
Definition: isrTask.py:1826
def suspectDetection(self, exposure, amp)
Definition: isrTask.py:2248
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:1199
def darkCorrection(self, exposure, darkExposure, invert=False)
Definition: isrTask.py:2112
def __init__(self, *args, **kwargs)
Definition: isrTask.py:2634
def runDataRef(self, dataRef)
Definition: isrTask.py:2638
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