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