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