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