LSSTApplications  20.0.0
LSSTDataManagementBasePackage
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  # Initial CCD-level background statistics options.
639  doMeasureBackground = pexConfig.Field(
640  dtype=bool,
641  doc="Measure the background level on the reduced image?",
642  default=False,
643  )
644 
645  # Camera-specific masking configuration.
646  doCameraSpecificMasking = pexConfig.Field(
647  dtype=bool,
648  doc="Mask camera-specific bad regions?",
649  default=False,
650  )
651  masking = pexConfig.ConfigurableField(
652  target=MaskingTask,
653  doc="Masking task."
654  )
655 
656  # Interpolation options.
657 
658  doInterpolate = pexConfig.Field(
659  dtype=bool,
660  doc="Interpolate masked pixels?",
661  default=True,
662  )
663  doSaturationInterpolation = pexConfig.Field(
664  dtype=bool,
665  doc="Perform interpolation over pixels masked as saturated?"
666  " NB: This is independent of doSaturation; if that is False this plane"
667  " will likely be blank, resulting in a no-op here.",
668  default=True,
669  )
670  doNanInterpolation = pexConfig.Field(
671  dtype=bool,
672  doc="Perform interpolation over pixels masked as NaN?"
673  " NB: This is independent of doNanMasking; if that is False this plane"
674  " will likely be blank, resulting in a no-op here.",
675  default=True,
676  )
677  doNanInterpAfterFlat = pexConfig.Field(
678  dtype=bool,
679  doc=("If True, ensure we interpolate NaNs after flat-fielding, even if we "
680  "also have to interpolate them before flat-fielding."),
681  default=False,
682  )
683  maskListToInterpolate = pexConfig.ListField(
684  dtype=str,
685  doc="List of mask planes that should be interpolated.",
686  default=['SAT', 'BAD', 'UNMASKEDNAN'],
687  )
688  doSaveInterpPixels = pexConfig.Field(
689  dtype=bool,
690  doc="Save a copy of the pre-interpolated pixel values?",
691  default=False,
692  )
693 
694  # Default photometric calibration options.
695  fluxMag0T1 = pexConfig.DictField(
696  keytype=str,
697  itemtype=float,
698  doc="The approximate flux of a zero-magnitude object in a one-second exposure, per filter.",
699  default=dict((f, pow(10.0, 0.4*m)) for f, m in (("Unknown", 28.0),
700  ))
701  )
702  defaultFluxMag0T1 = pexConfig.Field(
703  dtype=float,
704  doc="Default value for fluxMag0T1 (for an unrecognized filter).",
705  default=pow(10.0, 0.4*28.0)
706  )
707 
708  # Vignette correction configuration.
709  doVignette = pexConfig.Field(
710  dtype=bool,
711  doc="Apply vignetting parameters?",
712  default=False,
713  )
714  vignette = pexConfig.ConfigurableField(
715  target=VignetteTask,
716  doc="Vignetting task.",
717  )
718 
719  # Transmission curve configuration.
720  doAttachTransmissionCurve = pexConfig.Field(
721  dtype=bool,
722  default=False,
723  doc="Construct and attach a wavelength-dependent throughput curve for this CCD image?"
724  )
725  doUseOpticsTransmission = pexConfig.Field(
726  dtype=bool,
727  default=True,
728  doc="Load and use transmission_optics (if doAttachTransmissionCurve is True)?"
729  )
730  doUseFilterTransmission = pexConfig.Field(
731  dtype=bool,
732  default=True,
733  doc="Load and use transmission_filter (if doAttachTransmissionCurve is True)?"
734  )
735  doUseSensorTransmission = pexConfig.Field(
736  dtype=bool,
737  default=True,
738  doc="Load and use transmission_sensor (if doAttachTransmissionCurve is True)?"
739  )
740  doUseAtmosphereTransmission = pexConfig.Field(
741  dtype=bool,
742  default=True,
743  doc="Load and use transmission_atmosphere (if doAttachTransmissionCurve is True)?"
744  )
745 
746  # Illumination correction.
747  doIlluminationCorrection = pexConfig.Field(
748  dtype=bool,
749  default=False,
750  doc="Perform illumination correction?"
751  )
752  illuminationCorrectionDataProductName = pexConfig.Field(
753  dtype=str,
754  doc="Name of the illumination correction data product.",
755  default="illumcor",
756  )
757  illumScale = pexConfig.Field(
758  dtype=float,
759  doc="Scale factor for the illumination correction.",
760  default=1.0,
761  )
762  illumFilters = pexConfig.ListField(
763  dtype=str,
764  default=[],
765  doc="Only perform illumination correction for these filters."
766  )
767 
768  # Write the outputs to disk. If ISR is run as a subtask, this may not be needed.
769  doWrite = pexConfig.Field(
770  dtype=bool,
771  doc="Persist postISRCCD?",
772  default=True,
773  )
774 
775  def validate(self):
776  super().validate()
777  if self.doFlat and self.doApplyGains:
778  raise ValueError("You may not specify both doFlat and doApplyGains")
779  if self.doSaturationInterpolation and "SAT" not in self.maskListToInterpolate:
780  self.config.maskListToInterpolate.append("SAT")
781  if self.doNanInterpolation and "UNMASKEDNAN" not in self.maskListToInterpolate:
782  self.config.maskListToInterpolate.append("UNMASKEDNAN")
783 
784 
785 class IsrTask(pipeBase.PipelineTask, pipeBase.CmdLineTask):
786  """Apply common instrument signature correction algorithms to a raw frame.
787 
788  The process for correcting imaging data is very similar from
789  camera to camera. This task provides a vanilla implementation of
790  doing these corrections, including the ability to turn certain
791  corrections off if they are not needed. The inputs to the primary
792  method, `run()`, are a raw exposure to be corrected and the
793  calibration data products. The raw input is a single chip sized
794  mosaic of all amps including overscans and other non-science
795  pixels. The method `runDataRef()` identifies and defines the
796  calibration data products, and is intended for use by a
797  `lsst.pipe.base.cmdLineTask.CmdLineTask` and takes as input only a
798  `daf.persistence.butlerSubset.ButlerDataRef`. This task may be
799  subclassed for different camera, although the most camera specific
800  methods have been split into subtasks that can be redirected
801  appropriately.
802 
803  The __init__ method sets up the subtasks for ISR processing, using
804  the defaults from `lsst.ip.isr`.
805 
806  Parameters
807  ----------
808  args : `list`
809  Positional arguments passed to the Task constructor. None used at this time.
810  kwargs : `dict`, optional
811  Keyword arguments passed on to the Task constructor. None used at this time.
812  """
813  ConfigClass = IsrTaskConfig
814  _DefaultName = "isr"
815 
816  def __init__(self, **kwargs):
817  super().__init__(**kwargs)
818  self.makeSubtask("assembleCcd")
819  self.makeSubtask("crosstalk")
820  self.makeSubtask("strayLight")
821  self.makeSubtask("fringe")
822  self.makeSubtask("masking")
823  self.makeSubtask("overscan")
824  self.makeSubtask("vignette")
825 
826  def runQuantum(self, butlerQC, inputRefs, outputRefs):
827  inputs = butlerQC.get(inputRefs)
828 
829  try:
830  inputs['detectorNum'] = inputRefs.ccdExposure.dataId['detector']
831  except Exception as e:
832  raise ValueError("Failure to find valid detectorNum value for Dataset %s: %s." %
833  (inputRefs, e))
834 
835  inputs['isGen3'] = True
836 
837  detector = inputs['ccdExposure'].getDetector()
838 
839  if self.doLinearize(detector) is True:
840  if 'linearizer' in inputs and isinstance(inputs['linearizer'], dict):
841  linearizer = linearize.Linearizer(detector=detector, log=self.log)
842  linearizer.fromYaml(inputs['linearizer'])
843  else:
844  linearizer = linearize.Linearizer(table=inputs.get('linearizer', None), detector=detector,
845  log=self.log)
846  inputs['linearizer'] = linearizer
847 
848  if self.config.doDefect is True:
849  if "defects" in inputs and inputs['defects'] is not None:
850  # defects is loaded as a BaseCatalog with columns x0, y0, width, height.
851  # masking expects a list of defects defined by their bounding box
852  if not isinstance(inputs["defects"], Defects):
853  inputs["defects"] = Defects.fromTable(inputs["defects"])
854 
855  # Load the correct style of brighter fatter kernel, and repack
856  # the information as a numpy array.
857  if self.config.doBrighterFatter:
858  brighterFatterKernel = inputs.pop('newBFKernel', None)
859  if brighterFatterKernel is None:
860  brighterFatterKernel = inputs.get('bfKernel', None)
861 
862  if brighterFatterKernel is not None and not isinstance(brighterFatterKernel, numpy.ndarray):
863  detId = detector.getId()
864  inputs['bfGains'] = brighterFatterKernel.gain
865  # If the kernel is not an ndarray, it's the cp_pipe version
866  # so extract the kernel for this detector, or raise an error
867  if self.config.brighterFatterLevel == 'DETECTOR':
868  if brighterFatterKernel.detectorKernel:
869  inputs['bfKernel'] = brighterFatterKernel.detectorKernel[detId]
870  elif brighterFatterKernel.detectorKernelFromAmpKernels:
871  inputs['bfKernel'] = brighterFatterKernel.detectorKernelFromAmpKernels[detId]
872  else:
873  raise RuntimeError("Failed to extract kernel from new-style BF kernel.")
874  else:
875  # TODO DM-15631 for implementing this
876  raise NotImplementedError("Per-amplifier brighter-fatter correction not implemented")
877 
878  # Broken: DM-17169
879  # ci_hsc does not use crosstalkSources, as it's intra-CCD CT only. This needs to be
880  # fixed for non-HSC cameras in the future.
881  # inputs['crosstalkSources'] = (self.crosstalk.prepCrosstalk(inputsIds['ccdExposure'])
882  # if self.config.doCrosstalk else None)
883 
884  if self.config.doFringe is True and self.fringe.checkFilter(inputs['ccdExposure']):
885  expId = inputs['ccdExposure'].getInfo().getVisitInfo().getExposureId()
886  inputs['fringes'] = self.fringe.loadFringes(inputs['fringes'],
887  expId=expId,
888  assembler=self.assembleCcd
889  if self.config.doAssembleIsrExposures else None)
890  else:
891  inputs['fringes'] = pipeBase.Struct(fringes=None)
892 
893  if self.config.doStrayLight is True and self.strayLight.checkFilter(inputs['ccdExposure']):
894  if 'strayLightData' not in inputs:
895  inputs['strayLightData'] = None
896 
897  outputs = self.run(**inputs)
898  butlerQC.put(outputs, outputRefs)
899 
900  def readIsrData(self, dataRef, rawExposure):
901  """!Retrieve necessary frames for instrument signature removal.
902 
903  Pre-fetching all required ISR data products limits the IO
904  required by the ISR. Any conflict between the calibration data
905  available and that needed for ISR is also detected prior to
906  doing processing, allowing it to fail quickly.
907 
908  Parameters
909  ----------
910  dataRef : `daf.persistence.butlerSubset.ButlerDataRef`
911  Butler reference of the detector data to be processed
912  rawExposure : `afw.image.Exposure`
913  The raw exposure that will later be corrected with the
914  retrieved calibration data; should not be modified in this
915  method.
916 
917  Returns
918  -------
919  result : `lsst.pipe.base.Struct`
920  Result struct with components (which may be `None`):
921  - ``bias``: bias calibration frame (`afw.image.Exposure`)
922  - ``linearizer``: functor for linearization (`ip.isr.linearize.LinearizeBase`)
923  - ``crosstalkSources``: list of possible crosstalk sources (`list`)
924  - ``dark``: dark calibration frame (`afw.image.Exposure`)
925  - ``flat``: flat calibration frame (`afw.image.Exposure`)
926  - ``bfKernel``: Brighter-Fatter kernel (`numpy.ndarray`)
927  - ``defects``: list of defects (`lsst.meas.algorithms.Defects`)
928  - ``fringes``: `lsst.pipe.base.Struct` with components:
929  - ``fringes``: fringe calibration frame (`afw.image.Exposure`)
930  - ``seed``: random seed derived from the ccdExposureId for random
931  number generator (`uint32`).
932  - ``opticsTransmission``: `lsst.afw.image.TransmissionCurve`
933  A ``TransmissionCurve`` that represents the throughput of the optics,
934  to be evaluated in focal-plane coordinates.
935  - ``filterTransmission`` : `lsst.afw.image.TransmissionCurve`
936  A ``TransmissionCurve`` that represents the throughput of the filter
937  itself, to be evaluated in focal-plane coordinates.
938  - ``sensorTransmission`` : `lsst.afw.image.TransmissionCurve`
939  A ``TransmissionCurve`` that represents the throughput of the sensor
940  itself, to be evaluated in post-assembly trimmed detector coordinates.
941  - ``atmosphereTransmission`` : `lsst.afw.image.TransmissionCurve`
942  A ``TransmissionCurve`` that represents the throughput of the
943  atmosphere, assumed to be spatially constant.
944  - ``strayLightData`` : `object`
945  An opaque object containing calibration information for
946  stray-light correction. If `None`, no correction will be
947  performed.
948  - ``illumMaskedImage`` : illumination correction image (`lsst.afw.image.MaskedImage`)
949 
950  Raises
951  ------
952  NotImplementedError :
953  Raised if a per-amplifier brighter-fatter kernel is requested by the configuration.
954  """
955  try:
956  dateObs = rawExposure.getInfo().getVisitInfo().getDate()
957  dateObs = dateObs.toPython().isoformat()
958  except RuntimeError:
959  self.log.warn("Unable to identify dateObs for rawExposure.")
960  dateObs = None
961 
962  ccd = rawExposure.getDetector()
963  filterName = afwImage.Filter(rawExposure.getFilter().getId()).getName() # Canonical name for filter
964  rawExposure.mask.addMaskPlane("UNMASKEDNAN") # needed to match pre DM-15862 processing.
965  biasExposure = (self.getIsrExposure(dataRef, self.config.biasDataProductName)
966  if self.config.doBias else None)
967  # immediate=True required for functors and linearizers are functors; see ticket DM-6515
968  linearizer = (dataRef.get("linearizer", immediate=True)
969  if self.doLinearize(ccd) else None)
970  if linearizer is not None and not isinstance(linearizer, numpy.ndarray):
971  linearizer.log = self.log
972  if isinstance(linearizer, numpy.ndarray):
973  linearizer = linearize.Linearizer(table=linearizer, detector=ccd)
974  crosstalkSources = (self.crosstalk.prepCrosstalk(dataRef)
975  if self.config.doCrosstalk else None)
976  darkExposure = (self.getIsrExposure(dataRef, self.config.darkDataProductName)
977  if self.config.doDark else None)
978  flatExposure = (self.getIsrExposure(dataRef, self.config.flatDataProductName,
979  dateObs=dateObs)
980  if self.config.doFlat else None)
981 
982  brighterFatterKernel = None
983  brighterFatterGains = None
984  if self.config.doBrighterFatter is True:
985  try:
986  # Use the new-style cp_pipe version of the kernel if it exists
987  # If using a new-style kernel, always use the self-consistent
988  # gains, i.e. the ones inside the kernel object itself
989  brighterFatterKernel = dataRef.get("brighterFatterKernel")
990  brighterFatterGains = brighterFatterKernel.gain
991  self.log.info("New style bright-fatter kernel (brighterFatterKernel) loaded")
992  except NoResults:
993  try: # Fall back to the old-style numpy-ndarray style kernel if necessary.
994  brighterFatterKernel = dataRef.get("bfKernel")
995  self.log.info("Old style bright-fatter kernel (np.array) loaded")
996  except NoResults:
997  brighterFatterKernel = None
998  if brighterFatterKernel is not None and not isinstance(brighterFatterKernel, numpy.ndarray):
999  # If the kernel is not an ndarray, it's the cp_pipe version
1000  # so extract the kernel for this detector, or raise an error
1001  if self.config.brighterFatterLevel == 'DETECTOR':
1002  if brighterFatterKernel.detectorKernel:
1003  brighterFatterKernel = brighterFatterKernel.detectorKernel[ccd.getId()]
1004  elif brighterFatterKernel.detectorKernelFromAmpKernels:
1005  brighterFatterKernel = brighterFatterKernel.detectorKernelFromAmpKernels[ccd.getId()]
1006  else:
1007  raise RuntimeError("Failed to extract kernel from new-style BF kernel.")
1008  else:
1009  # TODO DM-15631 for implementing this
1010  raise NotImplementedError("Per-amplifier brighter-fatter correction not implemented")
1011 
1012  defectList = (dataRef.get("defects")
1013  if self.config.doDefect else None)
1014  fringeStruct = (self.fringe.readFringes(dataRef, assembler=self.assembleCcd
1015  if self.config.doAssembleIsrExposures else None)
1016  if self.config.doFringe and self.fringe.checkFilter(rawExposure)
1017  else pipeBase.Struct(fringes=None))
1018 
1019  if self.config.doAttachTransmissionCurve:
1020  opticsTransmission = (dataRef.get("transmission_optics")
1021  if self.config.doUseOpticsTransmission else None)
1022  filterTransmission = (dataRef.get("transmission_filter")
1023  if self.config.doUseFilterTransmission else None)
1024  sensorTransmission = (dataRef.get("transmission_sensor")
1025  if self.config.doUseSensorTransmission else None)
1026  atmosphereTransmission = (dataRef.get("transmission_atmosphere")
1027  if self.config.doUseAtmosphereTransmission else None)
1028  else:
1029  opticsTransmission = None
1030  filterTransmission = None
1031  sensorTransmission = None
1032  atmosphereTransmission = None
1033 
1034  if self.config.doStrayLight:
1035  strayLightData = self.strayLight.readIsrData(dataRef, rawExposure)
1036  else:
1037  strayLightData = None
1038 
1039  illumMaskedImage = (self.getIsrExposure(dataRef,
1040  self.config.illuminationCorrectionDataProductName).getMaskedImage()
1041  if (self.config.doIlluminationCorrection and
1042  filterName in self.config.illumFilters)
1043  else None)
1044 
1045  # Struct should include only kwargs to run()
1046  return pipeBase.Struct(bias=biasExposure,
1047  linearizer=linearizer,
1048  crosstalkSources=crosstalkSources,
1049  dark=darkExposure,
1050  flat=flatExposure,
1051  bfKernel=brighterFatterKernel,
1052  bfGains=brighterFatterGains,
1053  defects=defectList,
1054  fringes=fringeStruct,
1055  opticsTransmission=opticsTransmission,
1056  filterTransmission=filterTransmission,
1057  sensorTransmission=sensorTransmission,
1058  atmosphereTransmission=atmosphereTransmission,
1059  strayLightData=strayLightData,
1060  illumMaskedImage=illumMaskedImage
1061  )
1062 
1063  @pipeBase.timeMethod
1064  def run(self, ccdExposure, camera=None, bias=None, linearizer=None, crosstalkSources=None,
1065  dark=None, flat=None, bfKernel=None, bfGains=None, defects=None,
1066  fringes=pipeBase.Struct(fringes=None), opticsTransmission=None, filterTransmission=None,
1067  sensorTransmission=None, atmosphereTransmission=None,
1068  detectorNum=None, strayLightData=None, illumMaskedImage=None,
1069  isGen3=False,
1070  ):
1071  """!Perform instrument signature removal on an exposure.
1072 
1073  Steps included in the ISR processing, in order performed, are:
1074  - saturation and suspect pixel masking
1075  - overscan subtraction
1076  - CCD assembly of individual amplifiers
1077  - bias subtraction
1078  - variance image construction
1079  - linearization of non-linear response
1080  - crosstalk masking
1081  - brighter-fatter correction
1082  - dark subtraction
1083  - fringe correction
1084  - stray light subtraction
1085  - flat correction
1086  - masking of known defects and camera specific features
1087  - vignette calculation
1088  - appending transmission curve and distortion model
1089 
1090  Parameters
1091  ----------
1092  ccdExposure : `lsst.afw.image.Exposure`
1093  The raw exposure that is to be run through ISR. The
1094  exposure is modified by this method.
1095  camera : `lsst.afw.cameraGeom.Camera`, optional
1096  The camera geometry for this exposure. Required if ``isGen3`` is
1097  `True` and one or more of ``ccdExposure``, ``bias``, ``dark``, or
1098  ``flat`` does not have an associated detector.
1099  bias : `lsst.afw.image.Exposure`, optional
1100  Bias calibration frame.
1101  linearizer : `lsst.ip.isr.linearize.LinearizeBase`, optional
1102  Functor for linearization.
1103  crosstalkSources : `list`, optional
1104  List of possible crosstalk sources.
1105  dark : `lsst.afw.image.Exposure`, optional
1106  Dark calibration frame.
1107  flat : `lsst.afw.image.Exposure`, optional
1108  Flat calibration frame.
1109  bfKernel : `numpy.ndarray`, optional
1110  Brighter-fatter kernel.
1111  bfGains : `dict` of `float`, optional
1112  Gains used to override the detector's nominal gains for the
1113  brighter-fatter correction. A dict keyed by amplifier name for
1114  the detector in question.
1115  defects : `lsst.meas.algorithms.Defects`, optional
1116  List of defects.
1117  fringes : `lsst.pipe.base.Struct`, optional
1118  Struct containing the fringe correction data, with
1119  elements:
1120  - ``fringes``: fringe calibration frame (`afw.image.Exposure`)
1121  - ``seed``: random seed derived from the ccdExposureId for random
1122  number generator (`uint32`)
1123  opticsTransmission: `lsst.afw.image.TransmissionCurve`, optional
1124  A ``TransmissionCurve`` that represents the throughput of the optics,
1125  to be evaluated in focal-plane coordinates.
1126  filterTransmission : `lsst.afw.image.TransmissionCurve`
1127  A ``TransmissionCurve`` that represents the throughput of the filter
1128  itself, to be evaluated in focal-plane coordinates.
1129  sensorTransmission : `lsst.afw.image.TransmissionCurve`
1130  A ``TransmissionCurve`` that represents the throughput of the sensor
1131  itself, to be evaluated in post-assembly trimmed detector coordinates.
1132  atmosphereTransmission : `lsst.afw.image.TransmissionCurve`
1133  A ``TransmissionCurve`` that represents the throughput of the
1134  atmosphere, assumed to be spatially constant.
1135  detectorNum : `int`, optional
1136  The integer number for the detector to process.
1137  isGen3 : bool, optional
1138  Flag this call to run() as using the Gen3 butler environment.
1139  strayLightData : `object`, optional
1140  Opaque object containing calibration information for stray-light
1141  correction. If `None`, no correction will be performed.
1142  illumMaskedImage : `lsst.afw.image.MaskedImage`, optional
1143  Illumination correction image.
1144 
1145  Returns
1146  -------
1147  result : `lsst.pipe.base.Struct`
1148  Result struct with component:
1149  - ``exposure`` : `afw.image.Exposure`
1150  The fully ISR corrected exposure.
1151  - ``outputExposure`` : `afw.image.Exposure`
1152  An alias for `exposure`
1153  - ``ossThumb`` : `numpy.ndarray`
1154  Thumbnail image of the exposure after overscan subtraction.
1155  - ``flattenedThumb`` : `numpy.ndarray`
1156  Thumbnail image of the exposure after flat-field correction.
1157 
1158  Raises
1159  ------
1160  RuntimeError
1161  Raised if a configuration option is set to True, but the
1162  required calibration data has not been specified.
1163 
1164  Notes
1165  -----
1166  The current processed exposure can be viewed by setting the
1167  appropriate lsstDebug entries in the `debug.display`
1168  dictionary. The names of these entries correspond to some of
1169  the IsrTaskConfig Boolean options, with the value denoting the
1170  frame to use. The exposure is shown inside the matching
1171  option check and after the processing of that step has
1172  finished. The steps with debug points are:
1173 
1174  doAssembleCcd
1175  doBias
1176  doCrosstalk
1177  doBrighterFatter
1178  doDark
1179  doFringe
1180  doStrayLight
1181  doFlat
1182 
1183  In addition, setting the "postISRCCD" entry displays the
1184  exposure after all ISR processing has finished.
1185 
1186  """
1187 
1188  if isGen3 is True:
1189  # Gen3 currently cannot automatically do configuration overrides.
1190  # DM-15257 looks to discuss this issue.
1191  # Configure input exposures;
1192  if detectorNum is None:
1193  raise RuntimeError("Must supply the detectorNum if running as Gen3.")
1194 
1195  ccdExposure = self.ensureExposure(ccdExposure, camera, detectorNum)
1196  bias = self.ensureExposure(bias, camera, detectorNum)
1197  dark = self.ensureExposure(dark, camera, detectorNum)
1198  flat = self.ensureExposure(flat, camera, detectorNum)
1199  else:
1200  if isinstance(ccdExposure, ButlerDataRef):
1201  return self.runDataRef(ccdExposure)
1202 
1203  ccd = ccdExposure.getDetector()
1204  filterName = afwImage.Filter(ccdExposure.getFilter().getId()).getName() # Canonical name for filter
1205 
1206  if not ccd:
1207  assert not self.config.doAssembleCcd, "You need a Detector to run assembleCcd."
1208  ccd = [FakeAmp(ccdExposure, self.config)]
1209 
1210  # Validate Input
1211  if self.config.doBias and bias is None:
1212  raise RuntimeError("Must supply a bias exposure if config.doBias=True.")
1213  if self.doLinearize(ccd) and linearizer is None:
1214  raise RuntimeError("Must supply a linearizer if config.doLinearize=True for this detector.")
1215  if self.config.doBrighterFatter and bfKernel is None:
1216  raise RuntimeError("Must supply a kernel if config.doBrighterFatter=True.")
1217  if self.config.doDark and dark is None:
1218  raise RuntimeError("Must supply a dark exposure if config.doDark=True.")
1219  if self.config.doFlat and flat is None:
1220  raise RuntimeError("Must supply a flat exposure if config.doFlat=True.")
1221  if self.config.doDefect and defects is None:
1222  raise RuntimeError("Must supply defects if config.doDefect=True.")
1223  if (self.config.doFringe and filterName in self.fringe.config.filters and
1224  fringes.fringes is None):
1225  # The `fringes` object needs to be a pipeBase.Struct, as
1226  # we use it as a `dict` for the parameters of
1227  # `FringeTask.run()`. The `fringes.fringes` `list` may
1228  # not be `None` if `doFringe=True`. Otherwise, raise.
1229  raise RuntimeError("Must supply fringe exposure as a pipeBase.Struct.")
1230  if (self.config.doIlluminationCorrection and filterName in self.config.illumFilters and
1231  illumMaskedImage is None):
1232  raise RuntimeError("Must supply an illumcor if config.doIlluminationCorrection=True.")
1233 
1234  # Begin ISR processing.
1235  if self.config.doConvertIntToFloat:
1236  self.log.info("Converting exposure to floating point values.")
1237  ccdExposure = self.convertIntToFloat(ccdExposure)
1238 
1239  # Amplifier level processing.
1240  overscans = []
1241  for amp in ccd:
1242  # if ccdExposure is one amp, check for coverage to prevent performing ops multiple times
1243  if ccdExposure.getBBox().contains(amp.getBBox()):
1244  # Check for fully masked bad amplifiers, and generate masks for SUSPECT and SATURATED values.
1245  badAmp = self.maskAmplifier(ccdExposure, amp, defects)
1246 
1247  if self.config.doOverscan and not badAmp:
1248  # Overscan correction on amp-by-amp basis.
1249  overscanResults = self.overscanCorrection(ccdExposure, amp)
1250  self.log.debug("Corrected overscan for amplifier %s.", amp.getName())
1251  if overscanResults is not None and \
1252  self.config.qa is not None and self.config.qa.saveStats is True:
1253  if isinstance(overscanResults.overscanFit, float):
1254  qaMedian = overscanResults.overscanFit
1255  qaStdev = float("NaN")
1256  else:
1257  qaStats = afwMath.makeStatistics(overscanResults.overscanFit,
1258  afwMath.MEDIAN | afwMath.STDEVCLIP)
1259  qaMedian = qaStats.getValue(afwMath.MEDIAN)
1260  qaStdev = qaStats.getValue(afwMath.STDEVCLIP)
1261 
1262  self.metadata.set(f"ISR OSCAN {amp.getName()} MEDIAN", qaMedian)
1263  self.metadata.set(f"ISR OSCAN {amp.getName()} STDEV", qaStdev)
1264  self.log.debug(" Overscan stats for amplifer %s: %f +/- %f",
1265  amp.getName(), qaMedian, qaStdev)
1266  ccdExposure.getMetadata().set('OVERSCAN', "Overscan corrected")
1267  else:
1268  if badAmp:
1269  self.log.warn("Amplifier %s is bad.", amp.getName())
1270  overscanResults = None
1271 
1272  overscans.append(overscanResults if overscanResults is not None else None)
1273  else:
1274  self.log.info("Skipped OSCAN for %s.", amp.getName())
1275 
1276  if self.config.doCrosstalk and self.config.doCrosstalkBeforeAssemble:
1277  self.log.info("Applying crosstalk correction.")
1278  self.crosstalk.run(ccdExposure, crosstalkSources=crosstalkSources)
1279  self.debugView(ccdExposure, "doCrosstalk")
1280 
1281  if self.config.doAssembleCcd:
1282  self.log.info("Assembling CCD from amplifiers.")
1283  ccdExposure = self.assembleCcd.assembleCcd(ccdExposure)
1284 
1285  if self.config.expectWcs and not ccdExposure.getWcs():
1286  self.log.warn("No WCS found in input exposure.")
1287  self.debugView(ccdExposure, "doAssembleCcd")
1288 
1289  ossThumb = None
1290  if self.config.qa.doThumbnailOss:
1291  ossThumb = isrQa.makeThumbnail(ccdExposure, isrQaConfig=self.config.qa)
1292 
1293  if self.config.doBias:
1294  self.log.info("Applying bias correction.")
1295  isrFunctions.biasCorrection(ccdExposure.getMaskedImage(), bias.getMaskedImage(),
1296  trimToFit=self.config.doTrimToMatchCalib)
1297  self.debugView(ccdExposure, "doBias")
1298 
1299  if self.config.doVariance:
1300  for amp, overscanResults in zip(ccd, overscans):
1301  if ccdExposure.getBBox().contains(amp.getBBox()):
1302  self.log.debug("Constructing variance map for amplifer %s.", amp.getName())
1303  ampExposure = ccdExposure.Factory(ccdExposure, amp.getBBox())
1304  if overscanResults is not None:
1305  self.updateVariance(ampExposure, amp,
1306  overscanImage=overscanResults.overscanImage)
1307  else:
1308  self.updateVariance(ampExposure, amp,
1309  overscanImage=None)
1310  if self.config.qa is not None and self.config.qa.saveStats is True:
1311  qaStats = afwMath.makeStatistics(ampExposure.getVariance(),
1312  afwMath.MEDIAN | afwMath.STDEVCLIP)
1313  self.metadata.set(f"ISR VARIANCE {amp.getName()} MEDIAN",
1314  qaStats.getValue(afwMath.MEDIAN))
1315  self.metadata.set(f"ISR VARIANCE {amp.getName()} STDEV",
1316  qaStats.getValue(afwMath.STDEVCLIP))
1317  self.log.debug(" Variance stats for amplifer %s: %f +/- %f.",
1318  amp.getName(), qaStats.getValue(afwMath.MEDIAN),
1319  qaStats.getValue(afwMath.STDEVCLIP))
1320 
1321  if self.doLinearize(ccd):
1322  self.log.info("Applying linearizer.")
1323  linearizer.applyLinearity(image=ccdExposure.getMaskedImage().getImage(),
1324  detector=ccd, log=self.log)
1325 
1326  if self.config.doCrosstalk and not self.config.doCrosstalkBeforeAssemble:
1327  self.log.info("Applying crosstalk correction.")
1328  self.crosstalk.run(ccdExposure, crosstalkSources=crosstalkSources, isTrimmed=True)
1329  self.debugView(ccdExposure, "doCrosstalk")
1330 
1331  # Masking block. Optionally mask known defects, NAN pixels, widen trails, and do
1332  # anything else the camera needs. Saturated and suspect pixels have already been masked.
1333  if self.config.doDefect:
1334  self.log.info("Masking defects.")
1335  self.maskDefect(ccdExposure, defects)
1336 
1337  if self.config.numEdgeSuspect > 0:
1338  self.log.info("Masking edges as SUSPECT.")
1339  self.maskEdges(ccdExposure, numEdgePixels=self.config.numEdgeSuspect,
1340  maskPlane="SUSPECT")
1341 
1342  if self.config.doNanMasking:
1343  self.log.info("Masking NAN value pixels.")
1344  self.maskNan(ccdExposure)
1345 
1346  if self.config.doWidenSaturationTrails:
1347  self.log.info("Widening saturation trails.")
1348  isrFunctions.widenSaturationTrails(ccdExposure.getMaskedImage().getMask())
1349 
1350  if self.config.doCameraSpecificMasking:
1351  self.log.info("Masking regions for camera specific reasons.")
1352  self.masking.run(ccdExposure)
1353 
1354  if self.config.doBrighterFatter:
1355  # We need to apply flats and darks before we can interpolate, and we
1356  # need to interpolate before we do B-F, but we do B-F without the
1357  # flats and darks applied so we can work in units of electrons or holes.
1358  # This context manager applies and then removes the darks and flats.
1359  #
1360  # We also do not want to interpolate values here, so operate on temporary
1361  # images so we can apply only the BF-correction and roll back the
1362  # interpolation.
1363  interpExp = ccdExposure.clone()
1364  with self.flatContext(interpExp, flat, dark):
1365  isrFunctions.interpolateFromMask(
1366  maskedImage=interpExp.getMaskedImage(),
1367  fwhm=self.config.fwhm,
1368  growSaturatedFootprints=self.config.growSaturationFootprintSize,
1369  maskNameList=self.config.maskListToInterpolate
1370  )
1371  bfExp = interpExp.clone()
1372 
1373  self.log.info("Applying brighter fatter correction using kernel type %s / gains %s.",
1374  type(bfKernel), type(bfGains))
1375  bfResults = isrFunctions.brighterFatterCorrection(bfExp, bfKernel,
1376  self.config.brighterFatterMaxIter,
1377  self.config.brighterFatterThreshold,
1378  self.config.brighterFatterApplyGain,
1379  bfGains)
1380  if bfResults[1] == self.config.brighterFatterMaxIter:
1381  self.log.warn("Brighter fatter correction did not converge, final difference %f.",
1382  bfResults[0])
1383  else:
1384  self.log.info("Finished brighter fatter correction in %d iterations.",
1385  bfResults[1])
1386  image = ccdExposure.getMaskedImage().getImage()
1387  bfCorr = bfExp.getMaskedImage().getImage()
1388  bfCorr -= interpExp.getMaskedImage().getImage()
1389  image += bfCorr
1390 
1391  # Applying the brighter-fatter correction applies a
1392  # convolution to the science image. At the edges this
1393  # convolution may not have sufficient valid pixels to
1394  # produce a valid correction. Mark pixels within the size
1395  # of the brighter-fatter kernel as EDGE to warn of this
1396  # fact.
1397  self.log.info("Ensuring image edges are masked as SUSPECT to the brighter-fatter kernel size.")
1398  self.maskEdges(ccdExposure, numEdgePixels=numpy.max(bfKernel.shape) // 2,
1399  maskPlane="EDGE")
1400 
1401  if self.config.brighterFatterMaskGrowSize > 0:
1402  self.log.info("Growing masks to account for brighter-fatter kernel convolution.")
1403  for maskPlane in self.config.maskListToInterpolate:
1404  isrFunctions.growMasks(ccdExposure.getMask(),
1405  radius=self.config.brighterFatterMaskGrowSize,
1406  maskNameList=maskPlane,
1407  maskValue=maskPlane)
1408 
1409  self.debugView(ccdExposure, "doBrighterFatter")
1410 
1411  if self.config.doDark:
1412  self.log.info("Applying dark correction.")
1413  self.darkCorrection(ccdExposure, dark)
1414  self.debugView(ccdExposure, "doDark")
1415 
1416  if self.config.doFringe and not self.config.fringeAfterFlat:
1417  self.log.info("Applying fringe correction before flat.")
1418  self.fringe.run(ccdExposure, **fringes.getDict())
1419  self.debugView(ccdExposure, "doFringe")
1420 
1421  if self.config.doStrayLight and self.strayLight.check(ccdExposure):
1422  self.log.info("Checking strayLight correction.")
1423  self.strayLight.run(ccdExposure, strayLightData)
1424  self.debugView(ccdExposure, "doStrayLight")
1425 
1426  if self.config.doFlat:
1427  self.log.info("Applying flat correction.")
1428  self.flatCorrection(ccdExposure, flat)
1429  self.debugView(ccdExposure, "doFlat")
1430 
1431  if self.config.doApplyGains:
1432  self.log.info("Applying gain correction instead of flat.")
1433  isrFunctions.applyGains(ccdExposure, self.config.normalizeGains)
1434 
1435  if self.config.doFringe and self.config.fringeAfterFlat:
1436  self.log.info("Applying fringe correction after flat.")
1437  self.fringe.run(ccdExposure, **fringes.getDict())
1438 
1439  if self.config.doVignette:
1440  self.log.info("Constructing Vignette polygon.")
1441  self.vignettePolygon = self.vignette.run(ccdExposure)
1442 
1443  if self.config.vignette.doWriteVignettePolygon:
1444  self.setValidPolygonIntersect(ccdExposure, self.vignettePolygon)
1445 
1446  if self.config.doAttachTransmissionCurve:
1447  self.log.info("Adding transmission curves.")
1448  isrFunctions.attachTransmissionCurve(ccdExposure, opticsTransmission=opticsTransmission,
1449  filterTransmission=filterTransmission,
1450  sensorTransmission=sensorTransmission,
1451  atmosphereTransmission=atmosphereTransmission)
1452 
1453  flattenedThumb = None
1454  if self.config.qa.doThumbnailFlattened:
1455  flattenedThumb = isrQa.makeThumbnail(ccdExposure, isrQaConfig=self.config.qa)
1456 
1457  if self.config.doIlluminationCorrection and filterName in self.config.illumFilters:
1458  self.log.info("Performing illumination correction.")
1459  isrFunctions.illuminationCorrection(ccdExposure.getMaskedImage(),
1460  illumMaskedImage, illumScale=self.config.illumScale,
1461  trimToFit=self.config.doTrimToMatchCalib)
1462 
1463  preInterpExp = None
1464  if self.config.doSaveInterpPixels:
1465  preInterpExp = ccdExposure.clone()
1466 
1467  # Reset and interpolate bad pixels.
1468  #
1469  # Large contiguous bad regions (which should have the BAD mask
1470  # bit set) should have their values set to the image median.
1471  # This group should include defects and bad amplifiers. As the
1472  # area covered by these defects are large, there's little
1473  # reason to expect that interpolation would provide a more
1474  # useful value.
1475  #
1476  # Smaller defects can be safely interpolated after the larger
1477  # regions have had their pixel values reset. This ensures
1478  # that the remaining defects adjacent to bad amplifiers (as an
1479  # example) do not attempt to interpolate extreme values.
1480  if self.config.doSetBadRegions:
1481  badPixelCount, badPixelValue = isrFunctions.setBadRegions(ccdExposure)
1482  if badPixelCount > 0:
1483  self.log.info("Set %d BAD pixels to %f.", badPixelCount, badPixelValue)
1484 
1485  if self.config.doInterpolate:
1486  self.log.info("Interpolating masked pixels.")
1487  isrFunctions.interpolateFromMask(
1488  maskedImage=ccdExposure.getMaskedImage(),
1489  fwhm=self.config.fwhm,
1490  growSaturatedFootprints=self.config.growSaturationFootprintSize,
1491  maskNameList=list(self.config.maskListToInterpolate)
1492  )
1493 
1494  self.roughZeroPoint(ccdExposure)
1495 
1496  if self.config.doMeasureBackground:
1497  self.log.info("Measuring background level.")
1498  self.measureBackground(ccdExposure, self.config.qa)
1499 
1500  if self.config.qa is not None and self.config.qa.saveStats is True:
1501  for amp in ccd:
1502  ampExposure = ccdExposure.Factory(ccdExposure, amp.getBBox())
1503  qaStats = afwMath.makeStatistics(ampExposure.getImage(),
1504  afwMath.MEDIAN | afwMath.STDEVCLIP)
1505  self.metadata.set("ISR BACKGROUND {} MEDIAN".format(amp.getName()),
1506  qaStats.getValue(afwMath.MEDIAN))
1507  self.metadata.set("ISR BACKGROUND {} STDEV".format(amp.getName()),
1508  qaStats.getValue(afwMath.STDEVCLIP))
1509  self.log.debug(" Background stats for amplifer %s: %f +/- %f",
1510  amp.getName(), qaStats.getValue(afwMath.MEDIAN),
1511  qaStats.getValue(afwMath.STDEVCLIP))
1512 
1513  self.debugView(ccdExposure, "postISRCCD")
1514 
1515  return pipeBase.Struct(
1516  exposure=ccdExposure,
1517  ossThumb=ossThumb,
1518  flattenedThumb=flattenedThumb,
1519 
1520  preInterpolatedExposure=preInterpExp,
1521  outputExposure=ccdExposure,
1522  outputOssThumbnail=ossThumb,
1523  outputFlattenedThumbnail=flattenedThumb,
1524  )
1525 
1526  @pipeBase.timeMethod
1527  def runDataRef(self, sensorRef):
1528  """Perform instrument signature removal on a ButlerDataRef of a Sensor.
1529 
1530  This method contains the `CmdLineTask` interface to the ISR
1531  processing. All IO is handled here, freeing the `run()` method
1532  to manage only pixel-level calculations. The steps performed
1533  are:
1534  - Read in necessary detrending/isr/calibration data.
1535  - Process raw exposure in `run()`.
1536  - Persist the ISR-corrected exposure as "postISRCCD" if
1537  config.doWrite=True.
1538 
1539  Parameters
1540  ----------
1541  sensorRef : `daf.persistence.butlerSubset.ButlerDataRef`
1542  DataRef of the detector data to be processed
1543 
1544  Returns
1545  -------
1546  result : `lsst.pipe.base.Struct`
1547  Result struct with component:
1548  - ``exposure`` : `afw.image.Exposure`
1549  The fully ISR corrected exposure.
1550 
1551  Raises
1552  ------
1553  RuntimeError
1554  Raised if a configuration option is set to True, but the
1555  required calibration data does not exist.
1556 
1557  """
1558  self.log.info("Performing ISR on sensor %s.", sensorRef.dataId)
1559 
1560  ccdExposure = sensorRef.get(self.config.datasetType)
1561 
1562  camera = sensorRef.get("camera")
1563  isrData = self.readIsrData(sensorRef, ccdExposure)
1564 
1565  result = self.run(ccdExposure, camera=camera, **isrData.getDict())
1566 
1567  if self.config.doWrite:
1568  sensorRef.put(result.exposure, "postISRCCD")
1569  if result.preInterpolatedExposure is not None:
1570  sensorRef.put(result.preInterpolatedExposure, "postISRCCD_uninterpolated")
1571  if result.ossThumb is not None:
1572  isrQa.writeThumbnail(sensorRef, result.ossThumb, "ossThumb")
1573  if result.flattenedThumb is not None:
1574  isrQa.writeThumbnail(sensorRef, result.flattenedThumb, "flattenedThumb")
1575 
1576  return result
1577 
1578  def getIsrExposure(self, dataRef, datasetType, dateObs=None, immediate=True):
1579  """!Retrieve a calibration dataset for removing instrument signature.
1580 
1581  Parameters
1582  ----------
1583 
1584  dataRef : `daf.persistence.butlerSubset.ButlerDataRef`
1585  DataRef of the detector data to find calibration datasets
1586  for.
1587  datasetType : `str`
1588  Type of dataset to retrieve (e.g. 'bias', 'flat', etc).
1589  dateObs : `str`, optional
1590  Date of the observation. Used to correct butler failures
1591  when using fallback filters.
1592  immediate : `Bool`
1593  If True, disable butler proxies to enable error handling
1594  within this routine.
1595 
1596  Returns
1597  -------
1598  exposure : `lsst.afw.image.Exposure`
1599  Requested calibration frame.
1600 
1601  Raises
1602  ------
1603  RuntimeError
1604  Raised if no matching calibration frame can be found.
1605  """
1606  try:
1607  exp = dataRef.get(datasetType, immediate=immediate)
1608  except Exception as exc1:
1609  if not self.config.fallbackFilterName:
1610  raise RuntimeError("Unable to retrieve %s for %s: %s." % (datasetType, dataRef.dataId, exc1))
1611  try:
1612  if self.config.useFallbackDate and dateObs:
1613  exp = dataRef.get(datasetType, filter=self.config.fallbackFilterName,
1614  dateObs=dateObs, immediate=immediate)
1615  else:
1616  exp = dataRef.get(datasetType, filter=self.config.fallbackFilterName, immediate=immediate)
1617  except Exception as exc2:
1618  raise RuntimeError("Unable to retrieve %s for %s, even with fallback filter %s: %s AND %s." %
1619  (datasetType, dataRef.dataId, self.config.fallbackFilterName, exc1, exc2))
1620  self.log.warn("Using fallback calibration from filter %s.", self.config.fallbackFilterName)
1621 
1622  if self.config.doAssembleIsrExposures:
1623  exp = self.assembleCcd.assembleCcd(exp)
1624  return exp
1625 
1626  def ensureExposure(self, inputExp, camera, detectorNum):
1627  """Ensure that the data returned by Butler is a fully constructed exposure.
1628 
1629  ISR requires exposure-level image data for historical reasons, so if we did
1630  not recieve that from Butler, construct it from what we have, modifying the
1631  input in place.
1632 
1633  Parameters
1634  ----------
1635  inputExp : `lsst.afw.image.Exposure`, `lsst.afw.image.DecoratedImageU`, or
1636  `lsst.afw.image.ImageF`
1637  The input data structure obtained from Butler.
1638  camera : `lsst.afw.cameraGeom.camera`
1639  The camera associated with the image. Used to find the appropriate
1640  detector.
1641  detectorNum : `int`
1642  The detector this exposure should match.
1643 
1644  Returns
1645  -------
1646  inputExp : `lsst.afw.image.Exposure`
1647  The re-constructed exposure, with appropriate detector parameters.
1648 
1649  Raises
1650  ------
1651  TypeError
1652  Raised if the input data cannot be used to construct an exposure.
1653  """
1654  if isinstance(inputExp, afwImage.DecoratedImageU):
1655  inputExp = afwImage.makeExposure(afwImage.makeMaskedImage(inputExp))
1656  elif isinstance(inputExp, afwImage.ImageF):
1657  inputExp = afwImage.makeExposure(afwImage.makeMaskedImage(inputExp))
1658  elif isinstance(inputExp, afwImage.MaskedImageF):
1659  inputExp = afwImage.makeExposure(inputExp)
1660  elif isinstance(inputExp, afwImage.Exposure):
1661  pass
1662  elif inputExp is None:
1663  # Assume this will be caught by the setup if it is a problem.
1664  return inputExp
1665  else:
1666  raise TypeError("Input Exposure is not known type in isrTask.ensureExposure: %s." %
1667  (type(inputExp), ))
1668 
1669  if inputExp.getDetector() is None:
1670  inputExp.setDetector(camera[detectorNum])
1671 
1672  return inputExp
1673 
1674  def convertIntToFloat(self, exposure):
1675  """Convert exposure image from uint16 to float.
1676 
1677  If the exposure does not need to be converted, the input is
1678  immediately returned. For exposures that are converted to use
1679  floating point pixels, the variance is set to unity and the
1680  mask to zero.
1681 
1682  Parameters
1683  ----------
1684  exposure : `lsst.afw.image.Exposure`
1685  The raw exposure to be converted.
1686 
1687  Returns
1688  -------
1689  newexposure : `lsst.afw.image.Exposure`
1690  The input ``exposure``, converted to floating point pixels.
1691 
1692  Raises
1693  ------
1694  RuntimeError
1695  Raised if the exposure type cannot be converted to float.
1696 
1697  """
1698  if isinstance(exposure, afwImage.ExposureF):
1699  # Nothing to be done
1700  self.log.debug("Exposure already of type float.")
1701  return exposure
1702  if not hasattr(exposure, "convertF"):
1703  raise RuntimeError("Unable to convert exposure (%s) to float." % type(exposure))
1704 
1705  newexposure = exposure.convertF()
1706  newexposure.variance[:] = 1
1707  newexposure.mask[:] = 0x0
1708 
1709  return newexposure
1710 
1711  def maskAmplifier(self, ccdExposure, amp, defects):
1712  """Identify bad amplifiers, saturated and suspect pixels.
1713 
1714  Parameters
1715  ----------
1716  ccdExposure : `lsst.afw.image.Exposure`
1717  Input exposure to be masked.
1718  amp : `lsst.afw.table.AmpInfoCatalog`
1719  Catalog of parameters defining the amplifier on this
1720  exposure to mask.
1721  defects : `lsst.meas.algorithms.Defects`
1722  List of defects. Used to determine if the entire
1723  amplifier is bad.
1724 
1725  Returns
1726  -------
1727  badAmp : `Bool`
1728  If this is true, the entire amplifier area is covered by
1729  defects and unusable.
1730 
1731  """
1732  maskedImage = ccdExposure.getMaskedImage()
1733 
1734  badAmp = False
1735 
1736  # Check if entire amp region is defined as a defect (need to use amp.getBBox() for correct
1737  # comparison with current defects definition.
1738  if defects is not None:
1739  badAmp = bool(sum([v.getBBox().contains(amp.getBBox()) for v in defects]))
1740 
1741  # In the case of a bad amp, we will set mask to "BAD" (here use amp.getRawBBox() for correct
1742  # association with pixels in current ccdExposure).
1743  if badAmp:
1744  dataView = afwImage.MaskedImageF(maskedImage, amp.getRawBBox(),
1745  afwImage.PARENT)
1746  maskView = dataView.getMask()
1747  maskView |= maskView.getPlaneBitMask("BAD")
1748  del maskView
1749  return badAmp
1750 
1751  # Mask remaining defects after assembleCcd() to allow for defects that cross amplifier boundaries.
1752  # Saturation and suspect pixels can be masked now, though.
1753  limits = dict()
1754  if self.config.doSaturation and not badAmp:
1755  limits.update({self.config.saturatedMaskName: amp.getSaturation()})
1756  if self.config.doSuspect and not badAmp:
1757  limits.update({self.config.suspectMaskName: amp.getSuspectLevel()})
1758  if math.isfinite(self.config.saturation):
1759  limits.update({self.config.saturatedMaskName: self.config.saturation})
1760 
1761  for maskName, maskThreshold in limits.items():
1762  if not math.isnan(maskThreshold):
1763  dataView = maskedImage.Factory(maskedImage, amp.getRawBBox())
1764  isrFunctions.makeThresholdMask(
1765  maskedImage=dataView,
1766  threshold=maskThreshold,
1767  growFootprints=0,
1768  maskName=maskName
1769  )
1770 
1771  # Determine if we've fully masked this amplifier with SUSPECT and SAT pixels.
1772  maskView = afwImage.Mask(maskedImage.getMask(), amp.getRawDataBBox(),
1773  afwImage.PARENT)
1774  maskVal = maskView.getPlaneBitMask([self.config.saturatedMaskName,
1775  self.config.suspectMaskName])
1776  if numpy.all(maskView.getArray() & maskVal > 0):
1777  badAmp = True
1778  maskView |= maskView.getPlaneBitMask("BAD")
1779 
1780  return badAmp
1781 
1782  def overscanCorrection(self, ccdExposure, amp):
1783  """Apply overscan correction in place.
1784 
1785  This method does initial pixel rejection of the overscan
1786  region. The overscan can also be optionally segmented to
1787  allow for discontinuous overscan responses to be fit
1788  separately. The actual overscan subtraction is performed by
1789  the `lsst.ip.isr.isrFunctions.overscanCorrection` function,
1790  which is called here after the amplifier is preprocessed.
1791 
1792  Parameters
1793  ----------
1794  ccdExposure : `lsst.afw.image.Exposure`
1795  Exposure to have overscan correction performed.
1796  amp : `lsst.afw.table.AmpInfoCatalog`
1797  The amplifier to consider while correcting the overscan.
1798 
1799  Returns
1800  -------
1801  overscanResults : `lsst.pipe.base.Struct`
1802  Result struct with components:
1803  - ``imageFit`` : scalar or `lsst.afw.image.Image`
1804  Value or fit subtracted from the amplifier image data.
1805  - ``overscanFit`` : scalar or `lsst.afw.image.Image`
1806  Value or fit subtracted from the overscan image data.
1807  - ``overscanImage`` : `lsst.afw.image.Image`
1808  Image of the overscan region with the overscan
1809  correction applied. This quantity is used to estimate
1810  the amplifier read noise empirically.
1811 
1812  Raises
1813  ------
1814  RuntimeError
1815  Raised if the ``amp`` does not contain raw pixel information.
1816 
1817  See Also
1818  --------
1819  lsst.ip.isr.isrFunctions.overscanCorrection
1820  """
1821  if amp.getRawHorizontalOverscanBBox().isEmpty():
1822  self.log.info("ISR_OSCAN: No overscan region. Not performing overscan correction.")
1823  return None
1824 
1825  statControl = afwMath.StatisticsControl()
1826  statControl.setAndMask(ccdExposure.mask.getPlaneBitMask("SAT"))
1827 
1828  # Determine the bounding boxes
1829  dataBBox = amp.getRawDataBBox()
1830  oscanBBox = amp.getRawHorizontalOverscanBBox()
1831  dx0 = 0
1832  dx1 = 0
1833 
1834  prescanBBox = amp.getRawPrescanBBox()
1835  if (oscanBBox.getBeginX() > prescanBBox.getBeginX()): # amp is at the right
1836  dx0 += self.config.overscanNumLeadingColumnsToSkip
1837  dx1 -= self.config.overscanNumTrailingColumnsToSkip
1838  else:
1839  dx0 += self.config.overscanNumTrailingColumnsToSkip
1840  dx1 -= self.config.overscanNumLeadingColumnsToSkip
1841 
1842  # Determine if we need to work on subregions of the amplifier and overscan.
1843  imageBBoxes = []
1844  overscanBBoxes = []
1845 
1846  if ((self.config.overscanBiasJump and
1847  self.config.overscanBiasJumpLocation) and
1848  (ccdExposure.getMetadata().exists(self.config.overscanBiasJumpKeyword) and
1849  ccdExposure.getMetadata().getScalar(self.config.overscanBiasJumpKeyword) in
1850  self.config.overscanBiasJumpDevices)):
1851  if amp.getReadoutCorner() in (ReadoutCorner.LL, ReadoutCorner.LR):
1852  yLower = self.config.overscanBiasJumpLocation
1853  yUpper = dataBBox.getHeight() - yLower
1854  else:
1855  yUpper = self.config.overscanBiasJumpLocation
1856  yLower = dataBBox.getHeight() - yUpper
1857 
1858  imageBBoxes.append(lsst.geom.Box2I(dataBBox.getBegin(),
1859  lsst.geom.Extent2I(dataBBox.getWidth(), yLower)))
1860  overscanBBoxes.append(lsst.geom.Box2I(oscanBBox.getBegin() +
1861  lsst.geom.Extent2I(dx0, 0),
1862  lsst.geom.Extent2I(oscanBBox.getWidth() - dx0 + dx1,
1863  yLower)))
1864 
1865  imageBBoxes.append(lsst.geom.Box2I(dataBBox.getBegin() + lsst.geom.Extent2I(0, yLower),
1866  lsst.geom.Extent2I(dataBBox.getWidth(), yUpper)))
1867  overscanBBoxes.append(lsst.geom.Box2I(oscanBBox.getBegin() + lsst.geom.Extent2I(dx0, yLower),
1868  lsst.geom.Extent2I(oscanBBox.getWidth() - dx0 + dx1,
1869  yUpper)))
1870  else:
1871  imageBBoxes.append(lsst.geom.Box2I(dataBBox.getBegin(),
1872  lsst.geom.Extent2I(dataBBox.getWidth(), dataBBox.getHeight())))
1873  overscanBBoxes.append(lsst.geom.Box2I(oscanBBox.getBegin() + lsst.geom.Extent2I(dx0, 0),
1874  lsst.geom.Extent2I(oscanBBox.getWidth() - dx0 + dx1,
1875  oscanBBox.getHeight())))
1876 
1877  # Perform overscan correction on subregions, ensuring saturated pixels are masked.
1878  for imageBBox, overscanBBox in zip(imageBBoxes, overscanBBoxes):
1879  ampImage = ccdExposure.maskedImage[imageBBox]
1880  overscanImage = ccdExposure.maskedImage[overscanBBox]
1881 
1882  overscanArray = overscanImage.image.array
1883  median = numpy.ma.median(numpy.ma.masked_where(overscanImage.mask.array, overscanArray))
1884  bad = numpy.where(numpy.abs(overscanArray - median) > self.config.overscanMaxDev)
1885  overscanImage.mask.array[bad] = overscanImage.mask.getPlaneBitMask("SAT")
1886 
1887  statControl = afwMath.StatisticsControl()
1888  statControl.setAndMask(ccdExposure.mask.getPlaneBitMask("SAT"))
1889 
1890  overscanResults = self.overscan.run(ampImage.getImage(), overscanImage)
1891 
1892  # Measure average overscan levels and record them in the metadata.
1893  levelStat = afwMath.MEDIAN
1894  sigmaStat = afwMath.STDEVCLIP
1895 
1896  sctrl = afwMath.StatisticsControl(self.config.qa.flatness.clipSigma,
1897  self.config.qa.flatness.nIter)
1898  metadata = ccdExposure.getMetadata()
1899  ampNum = amp.getName()
1900  # if self.config.overscanFitType in ("MEDIAN", "MEAN", "MEANCLIP"):
1901  if isinstance(overscanResults.overscanFit, float):
1902  metadata.set("ISR_OSCAN_LEVEL%s" % ampNum, overscanResults.overscanFit)
1903  metadata.set("ISR_OSCAN_SIGMA%s" % ampNum, 0.0)
1904  else:
1905  stats = afwMath.makeStatistics(overscanResults.overscanFit, levelStat | sigmaStat, sctrl)
1906  metadata.set("ISR_OSCAN_LEVEL%s" % ampNum, stats.getValue(levelStat))
1907  metadata.set("ISR_OSCAN_SIGMA%s" % ampNum, stats.getValue(sigmaStat))
1908 
1909  return overscanResults
1910 
1911  def updateVariance(self, ampExposure, amp, overscanImage=None):
1912  """Set the variance plane using the amplifier gain and read noise
1913 
1914  The read noise is calculated from the ``overscanImage`` if the
1915  ``doEmpiricalReadNoise`` option is set in the configuration; otherwise
1916  the value from the amplifier data is used.
1917 
1918  Parameters
1919  ----------
1920  ampExposure : `lsst.afw.image.Exposure`
1921  Exposure to process.
1922  amp : `lsst.afw.table.AmpInfoRecord` or `FakeAmp`
1923  Amplifier detector data.
1924  overscanImage : `lsst.afw.image.MaskedImage`, optional.
1925  Image of overscan, required only for empirical read noise.
1926 
1927  See also
1928  --------
1929  lsst.ip.isr.isrFunctions.updateVariance
1930  """
1931  maskPlanes = [self.config.saturatedMaskName, self.config.suspectMaskName]
1932  gain = amp.getGain()
1933 
1934  if math.isnan(gain):
1935  gain = 1.0
1936  self.log.warn("Gain set to NAN! Updating to 1.0 to generate Poisson variance.")
1937  elif gain <= 0:
1938  patchedGain = 1.0
1939  self.log.warn("Gain for amp %s == %g <= 0; setting to %f.",
1940  amp.getName(), gain, patchedGain)
1941  gain = patchedGain
1942 
1943  if self.config.doEmpiricalReadNoise and overscanImage is None:
1944  self.log.info("Overscan is none for EmpiricalReadNoise.")
1945 
1946  if self.config.doEmpiricalReadNoise and overscanImage is not None:
1947  stats = afwMath.StatisticsControl()
1948  stats.setAndMask(overscanImage.mask.getPlaneBitMask(maskPlanes))
1949  readNoise = afwMath.makeStatistics(overscanImage, afwMath.STDEVCLIP, stats).getValue()
1950  self.log.info("Calculated empirical read noise for amp %s: %f.",
1951  amp.getName(), readNoise)
1952  else:
1953  readNoise = amp.getReadNoise()
1954 
1955  isrFunctions.updateVariance(
1956  maskedImage=ampExposure.getMaskedImage(),
1957  gain=gain,
1958  readNoise=readNoise,
1959  )
1960 
1961  def darkCorrection(self, exposure, darkExposure, invert=False):
1962  """!Apply dark correction in place.
1963 
1964  Parameters
1965  ----------
1966  exposure : `lsst.afw.image.Exposure`
1967  Exposure to process.
1968  darkExposure : `lsst.afw.image.Exposure`
1969  Dark exposure of the same size as ``exposure``.
1970  invert : `Bool`, optional
1971  If True, re-add the dark to an already corrected image.
1972 
1973  Raises
1974  ------
1975  RuntimeError
1976  Raised if either ``exposure`` or ``darkExposure`` do not
1977  have their dark time defined.
1978 
1979  See Also
1980  --------
1981  lsst.ip.isr.isrFunctions.darkCorrection
1982  """
1983  expScale = exposure.getInfo().getVisitInfo().getDarkTime()
1984  if math.isnan(expScale):
1985  raise RuntimeError("Exposure darktime is NAN.")
1986  if darkExposure.getInfo().getVisitInfo() is not None \
1987  and not math.isnan(darkExposure.getInfo().getVisitInfo().getDarkTime()):
1988  darkScale = darkExposure.getInfo().getVisitInfo().getDarkTime()
1989  else:
1990  # DM-17444: darkExposure.getInfo.getVisitInfo() is None
1991  # so getDarkTime() does not exist.
1992  self.log.warn("darkExposure.getInfo().getVisitInfo() does not exist. Using darkScale = 1.0.")
1993  darkScale = 1.0
1994 
1995  isrFunctions.darkCorrection(
1996  maskedImage=exposure.getMaskedImage(),
1997  darkMaskedImage=darkExposure.getMaskedImage(),
1998  expScale=expScale,
1999  darkScale=darkScale,
2000  invert=invert,
2001  trimToFit=self.config.doTrimToMatchCalib
2002  )
2003 
2004  def doLinearize(self, detector):
2005  """!Check if linearization is needed for the detector cameraGeom.
2006 
2007  Checks config.doLinearize and the linearity type of the first
2008  amplifier.
2009 
2010  Parameters
2011  ----------
2012  detector : `lsst.afw.cameraGeom.Detector`
2013  Detector to get linearity type from.
2014 
2015  Returns
2016  -------
2017  doLinearize : `Bool`
2018  If True, linearization should be performed.
2019  """
2020  return self.config.doLinearize and \
2021  detector.getAmplifiers()[0].getLinearityType() != NullLinearityType
2022 
2023  def flatCorrection(self, exposure, flatExposure, invert=False):
2024  """!Apply flat correction in place.
2025 
2026  Parameters
2027  ----------
2028  exposure : `lsst.afw.image.Exposure`
2029  Exposure to process.
2030  flatExposure : `lsst.afw.image.Exposure`
2031  Flat exposure of the same size as ``exposure``.
2032  invert : `Bool`, optional
2033  If True, unflatten an already flattened image.
2034 
2035  See Also
2036  --------
2037  lsst.ip.isr.isrFunctions.flatCorrection
2038  """
2039  isrFunctions.flatCorrection(
2040  maskedImage=exposure.getMaskedImage(),
2041  flatMaskedImage=flatExposure.getMaskedImage(),
2042  scalingType=self.config.flatScalingType,
2043  userScale=self.config.flatUserScale,
2044  invert=invert,
2045  trimToFit=self.config.doTrimToMatchCalib
2046  )
2047 
2048  def saturationDetection(self, exposure, amp):
2049  """!Detect saturated pixels and mask them using mask plane config.saturatedMaskName, in place.
2050 
2051  Parameters
2052  ----------
2053  exposure : `lsst.afw.image.Exposure`
2054  Exposure to process. Only the amplifier DataSec is processed.
2055  amp : `lsst.afw.table.AmpInfoCatalog`
2056  Amplifier detector data.
2057 
2058  See Also
2059  --------
2060  lsst.ip.isr.isrFunctions.makeThresholdMask
2061  """
2062  if not math.isnan(amp.getSaturation()):
2063  maskedImage = exposure.getMaskedImage()
2064  dataView = maskedImage.Factory(maskedImage, amp.getRawBBox())
2065  isrFunctions.makeThresholdMask(
2066  maskedImage=dataView,
2067  threshold=amp.getSaturation(),
2068  growFootprints=0,
2069  maskName=self.config.saturatedMaskName,
2070  )
2071 
2072  def saturationInterpolation(self, exposure):
2073  """!Interpolate over saturated pixels, in place.
2074 
2075  This method should be called after `saturationDetection`, to
2076  ensure that the saturated pixels have been identified in the
2077  SAT mask. It should also be called after `assembleCcd`, since
2078  saturated regions may cross amplifier boundaries.
2079 
2080  Parameters
2081  ----------
2082  exposure : `lsst.afw.image.Exposure`
2083  Exposure to process.
2084 
2085  See Also
2086  --------
2087  lsst.ip.isr.isrTask.saturationDetection
2088  lsst.ip.isr.isrFunctions.interpolateFromMask
2089  """
2090  isrFunctions.interpolateFromMask(
2091  maskedImage=exposure.getMaskedImage(),
2092  fwhm=self.config.fwhm,
2093  growSaturatedFootprints=self.config.growSaturationFootprintSize,
2094  maskNameList=list(self.config.saturatedMaskName),
2095  )
2096 
2097  def suspectDetection(self, exposure, amp):
2098  """!Detect suspect pixels and mask them using mask plane config.suspectMaskName, in place.
2099 
2100  Parameters
2101  ----------
2102  exposure : `lsst.afw.image.Exposure`
2103  Exposure to process. Only the amplifier DataSec is processed.
2104  amp : `lsst.afw.table.AmpInfoCatalog`
2105  Amplifier detector data.
2106 
2107  See Also
2108  --------
2109  lsst.ip.isr.isrFunctions.makeThresholdMask
2110 
2111  Notes
2112  -----
2113  Suspect pixels are pixels whose value is greater than amp.getSuspectLevel().
2114  This is intended to indicate pixels that may be affected by unknown systematics;
2115  for example if non-linearity corrections above a certain level are unstable
2116  then that would be a useful value for suspectLevel. A value of `nan` indicates
2117  that no such level exists and no pixels are to be masked as suspicious.
2118  """
2119  suspectLevel = amp.getSuspectLevel()
2120  if math.isnan(suspectLevel):
2121  return
2122 
2123  maskedImage = exposure.getMaskedImage()
2124  dataView = maskedImage.Factory(maskedImage, amp.getRawBBox())
2125  isrFunctions.makeThresholdMask(
2126  maskedImage=dataView,
2127  threshold=suspectLevel,
2128  growFootprints=0,
2129  maskName=self.config.suspectMaskName,
2130  )
2131 
2132  def maskDefect(self, exposure, defectBaseList):
2133  """!Mask defects using mask plane "BAD", in place.
2134 
2135  Parameters
2136  ----------
2137  exposure : `lsst.afw.image.Exposure`
2138  Exposure to process.
2139  defectBaseList : `lsst.meas.algorithms.Defects` or `list` of
2140  `lsst.afw.image.DefectBase`.
2141  List of defects to mask.
2142 
2143  Notes
2144  -----
2145  Call this after CCD assembly, since defects may cross amplifier boundaries.
2146  """
2147  maskedImage = exposure.getMaskedImage()
2148  if not isinstance(defectBaseList, Defects):
2149  # Promotes DefectBase to Defect
2150  defectList = Defects(defectBaseList)
2151  else:
2152  defectList = defectBaseList
2153  defectList.maskPixels(maskedImage, maskName="BAD")
2154 
2155  def maskEdges(self, exposure, numEdgePixels=0, maskPlane="SUSPECT"):
2156  """!Mask edge pixels with applicable mask plane.
2157 
2158  Parameters
2159  ----------
2160  exposure : `lsst.afw.image.Exposure`
2161  Exposure to process.
2162  numEdgePixels : `int`, optional
2163  Number of edge pixels to mask.
2164  maskPlane : `str`, optional
2165  Mask plane name to use.
2166  """
2167  maskedImage = exposure.getMaskedImage()
2168  maskBitMask = maskedImage.getMask().getPlaneBitMask(maskPlane)
2169 
2170  if numEdgePixels > 0:
2171  goodBBox = maskedImage.getBBox()
2172  # This makes a bbox numEdgeSuspect pixels smaller than the image on each side
2173  goodBBox.grow(-numEdgePixels)
2174  # Mask pixels outside goodBBox
2175  SourceDetectionTask.setEdgeBits(
2176  maskedImage,
2177  goodBBox,
2178  maskBitMask
2179  )
2180 
2181  def maskAndInterpolateDefects(self, exposure, defectBaseList):
2182  """Mask and interpolate defects using mask plane "BAD", in place.
2183 
2184  Parameters
2185  ----------
2186  exposure : `lsst.afw.image.Exposure`
2187  Exposure to process.
2188  defectBaseList : `lsst.meas.algorithms.Defects` or `list` of
2189  `lsst.afw.image.DefectBase`.
2190  List of defects to mask and interpolate.
2191 
2192  See Also
2193  --------
2194  lsst.ip.isr.isrTask.maskDefect()
2195  """
2196  self.maskDefect(exposure, defectBaseList)
2197  self.maskEdges(exposure, numEdgePixels=self.config.numEdgeSuspect,
2198  maskPlane="SUSPECT")
2199  isrFunctions.interpolateFromMask(
2200  maskedImage=exposure.getMaskedImage(),
2201  fwhm=self.config.fwhm,
2202  growSaturatedFootprints=0,
2203  maskNameList=["BAD"],
2204  )
2205 
2206  def maskNan(self, exposure):
2207  """Mask NaNs using mask plane "UNMASKEDNAN", in place.
2208 
2209  Parameters
2210  ----------
2211  exposure : `lsst.afw.image.Exposure`
2212  Exposure to process.
2213 
2214  Notes
2215  -----
2216  We mask over all NaNs, including those that are masked with
2217  other bits (because those may or may not be interpolated over
2218  later, and we want to remove all NaNs). Despite this
2219  behaviour, the "UNMASKEDNAN" mask plane is used to preserve
2220  the historical name.
2221  """
2222  maskedImage = exposure.getMaskedImage()
2223 
2224  # Find and mask NaNs
2225  maskedImage.getMask().addMaskPlane("UNMASKEDNAN")
2226  maskVal = maskedImage.getMask().getPlaneBitMask("UNMASKEDNAN")
2227  numNans = maskNans(maskedImage, maskVal)
2228  self.metadata.set("NUMNANS", numNans)
2229  if numNans > 0:
2230  self.log.warn("There were %d unmasked NaNs.", numNans)
2231 
2232  def maskAndInterpolateNan(self, exposure):
2233  """"Mask and interpolate NaNs using mask plane "UNMASKEDNAN", in place.
2234 
2235  Parameters
2236  ----------
2237  exposure : `lsst.afw.image.Exposure`
2238  Exposure to process.
2239 
2240  See Also
2241  --------
2242  lsst.ip.isr.isrTask.maskNan()
2243  """
2244  self.maskNan(exposure)
2245  isrFunctions.interpolateFromMask(
2246  maskedImage=exposure.getMaskedImage(),
2247  fwhm=self.config.fwhm,
2248  growSaturatedFootprints=0,
2249  maskNameList=["UNMASKEDNAN"],
2250  )
2251 
2252  def measureBackground(self, exposure, IsrQaConfig=None):
2253  """Measure the image background in subgrids, for quality control purposes.
2254 
2255  Parameters
2256  ----------
2257  exposure : `lsst.afw.image.Exposure`
2258  Exposure to process.
2259  IsrQaConfig : `lsst.ip.isr.isrQa.IsrQaConfig`
2260  Configuration object containing parameters on which background
2261  statistics and subgrids to use.
2262  """
2263  if IsrQaConfig is not None:
2264  statsControl = afwMath.StatisticsControl(IsrQaConfig.flatness.clipSigma,
2265  IsrQaConfig.flatness.nIter)
2266  maskVal = exposure.getMaskedImage().getMask().getPlaneBitMask(["BAD", "SAT", "DETECTED"])
2267  statsControl.setAndMask(maskVal)
2268  maskedImage = exposure.getMaskedImage()
2269  stats = afwMath.makeStatistics(maskedImage, afwMath.MEDIAN | afwMath.STDEVCLIP, statsControl)
2270  skyLevel = stats.getValue(afwMath.MEDIAN)
2271  skySigma = stats.getValue(afwMath.STDEVCLIP)
2272  self.log.info("Flattened sky level: %f +/- %f.", skyLevel, skySigma)
2273  metadata = exposure.getMetadata()
2274  metadata.set('SKYLEVEL', skyLevel)
2275  metadata.set('SKYSIGMA', skySigma)
2276 
2277  # calcluating flatlevel over the subgrids
2278  stat = afwMath.MEANCLIP if IsrQaConfig.flatness.doClip else afwMath.MEAN
2279  meshXHalf = int(IsrQaConfig.flatness.meshX/2.)
2280  meshYHalf = int(IsrQaConfig.flatness.meshY/2.)
2281  nX = int((exposure.getWidth() + meshXHalf) / IsrQaConfig.flatness.meshX)
2282  nY = int((exposure.getHeight() + meshYHalf) / IsrQaConfig.flatness.meshY)
2283  skyLevels = numpy.zeros((nX, nY))
2284 
2285  for j in range(nY):
2286  yc = meshYHalf + j * IsrQaConfig.flatness.meshY
2287  for i in range(nX):
2288  xc = meshXHalf + i * IsrQaConfig.flatness.meshX
2289 
2290  xLLC = xc - meshXHalf
2291  yLLC = yc - meshYHalf
2292  xURC = xc + meshXHalf - 1
2293  yURC = yc + meshYHalf - 1
2294 
2295  bbox = lsst.geom.Box2I(lsst.geom.Point2I(xLLC, yLLC), lsst.geom.Point2I(xURC, yURC))
2296  miMesh = maskedImage.Factory(exposure.getMaskedImage(), bbox, afwImage.LOCAL)
2297 
2298  skyLevels[i, j] = afwMath.makeStatistics(miMesh, stat, statsControl).getValue()
2299 
2300  good = numpy.where(numpy.isfinite(skyLevels))
2301  skyMedian = numpy.median(skyLevels[good])
2302  flatness = (skyLevels[good] - skyMedian) / skyMedian
2303  flatness_rms = numpy.std(flatness)
2304  flatness_pp = flatness.max() - flatness.min() if len(flatness) > 0 else numpy.nan
2305 
2306  self.log.info("Measuring sky levels in %dx%d grids: %f.", nX, nY, skyMedian)
2307  self.log.info("Sky flatness in %dx%d grids - pp: %f rms: %f.",
2308  nX, nY, flatness_pp, flatness_rms)
2309 
2310  metadata.set('FLATNESS_PP', float(flatness_pp))
2311  metadata.set('FLATNESS_RMS', float(flatness_rms))
2312  metadata.set('FLATNESS_NGRIDS', '%dx%d' % (nX, nY))
2313  metadata.set('FLATNESS_MESHX', IsrQaConfig.flatness.meshX)
2314  metadata.set('FLATNESS_MESHY', IsrQaConfig.flatness.meshY)
2315 
2316  def roughZeroPoint(self, exposure):
2317  """Set an approximate magnitude zero point for the exposure.
2318 
2319  Parameters
2320  ----------
2321  exposure : `lsst.afw.image.Exposure`
2322  Exposure to process.
2323  """
2324  filterName = afwImage.Filter(exposure.getFilter().getId()).getName() # Canonical name for filter
2325  if filterName in self.config.fluxMag0T1:
2326  fluxMag0 = self.config.fluxMag0T1[filterName]
2327  else:
2328  self.log.warn("No rough magnitude zero point set for filter %s.", filterName)
2329  fluxMag0 = self.config.defaultFluxMag0T1
2330 
2331  expTime = exposure.getInfo().getVisitInfo().getExposureTime()
2332  if not expTime > 0: # handle NaN as well as <= 0
2333  self.log.warn("Non-positive exposure time; skipping rough zero point.")
2334  return
2335 
2336  self.log.info("Setting rough magnitude zero point: %f", 2.5*math.log10(fluxMag0*expTime))
2337  exposure.setPhotoCalib(afwImage.makePhotoCalibFromCalibZeroPoint(fluxMag0*expTime, 0.0))
2338 
2339  def setValidPolygonIntersect(self, ccdExposure, fpPolygon):
2340  """!Set the valid polygon as the intersection of fpPolygon and the ccd corners.
2341 
2342  Parameters
2343  ----------
2344  ccdExposure : `lsst.afw.image.Exposure`
2345  Exposure to process.
2346  fpPolygon : `lsst.afw.geom.Polygon`
2347  Polygon in focal plane coordinates.
2348  """
2349  # Get ccd corners in focal plane coordinates
2350  ccd = ccdExposure.getDetector()
2351  fpCorners = ccd.getCorners(FOCAL_PLANE)
2352  ccdPolygon = Polygon(fpCorners)
2353 
2354  # Get intersection of ccd corners with fpPolygon
2355  intersect = ccdPolygon.intersectionSingle(fpPolygon)
2356 
2357  # Transform back to pixel positions and build new polygon
2358  ccdPoints = ccd.transform(intersect, FOCAL_PLANE, PIXELS)
2359  validPolygon = Polygon(ccdPoints)
2360  ccdExposure.getInfo().setValidPolygon(validPolygon)
2361 
2362  @contextmanager
2363  def flatContext(self, exp, flat, dark=None):
2364  """Context manager that applies and removes flats and darks,
2365  if the task is configured to apply them.
2366 
2367  Parameters
2368  ----------
2369  exp : `lsst.afw.image.Exposure`
2370  Exposure to process.
2371  flat : `lsst.afw.image.Exposure`
2372  Flat exposure the same size as ``exp``.
2373  dark : `lsst.afw.image.Exposure`, optional
2374  Dark exposure the same size as ``exp``.
2375 
2376  Yields
2377  ------
2378  exp : `lsst.afw.image.Exposure`
2379  The flat and dark corrected exposure.
2380  """
2381  if self.config.doDark and dark is not None:
2382  self.darkCorrection(exp, dark)
2383  if self.config.doFlat:
2384  self.flatCorrection(exp, flat)
2385  try:
2386  yield exp
2387  finally:
2388  if self.config.doFlat:
2389  self.flatCorrection(exp, flat, invert=True)
2390  if self.config.doDark and dark is not None:
2391  self.darkCorrection(exp, dark, invert=True)
2392 
2393  def debugView(self, exposure, stepname):
2394  """Utility function to examine ISR exposure at different stages.
2395 
2396  Parameters
2397  ----------
2398  exposure : `lsst.afw.image.Exposure`
2399  Exposure to view.
2400  stepname : `str`
2401  State of processing to view.
2402  """
2403  frame = getDebugFrame(self._display, stepname)
2404  if frame:
2405  display = getDisplay(frame)
2406  display.scale('asinh', 'zscale')
2407  display.mtv(exposure)
2408  prompt = "Press Enter to continue [c]... "
2409  while True:
2410  ans = input(prompt).lower()
2411  if ans in ("", "c",):
2412  break
2413 
2414 
2416  """A Detector-like object that supports returning gain and saturation level
2417 
2418  This is used when the input exposure does not have a detector.
2419 
2420  Parameters
2421  ----------
2422  exposure : `lsst.afw.image.Exposure`
2423  Exposure to generate a fake amplifier for.
2424  config : `lsst.ip.isr.isrTaskConfig`
2425  Configuration to apply to the fake amplifier.
2426  """
2427 
2428  def __init__(self, exposure, config):
2429  self._bbox = exposure.getBBox(afwImage.LOCAL)
2431  self._gain = config.gain
2432  self._readNoise = config.readNoise
2433  self._saturation = config.saturation
2434 
2435  def getBBox(self):
2436  return self._bbox
2437 
2438  def getRawBBox(self):
2439  return self._bbox
2440 
2442  return self._RawHorizontalOverscanBBox
2443 
2444  def getGain(self):
2445  return self._gain
2446 
2447  def getReadNoise(self):
2448  return self._readNoise
2449 
2450  def getSaturation(self):
2451  return self._saturation
2452 
2453  def getSuspectLevel(self):
2454  return float("NaN")
2455 
2456 
2457 class RunIsrConfig(pexConfig.Config):
2458  isr = pexConfig.ConfigurableField(target=IsrTask, doc="Instrument signature removal")
2459 
2460 
2461 class RunIsrTask(pipeBase.CmdLineTask):
2462  """Task to wrap the default IsrTask to allow it to be retargeted.
2463 
2464  The standard IsrTask can be called directly from a command line
2465  program, but doing so removes the ability of the task to be
2466  retargeted. As most cameras override some set of the IsrTask
2467  methods, this would remove those data-specific methods in the
2468  output post-ISR images. This wrapping class fixes the issue,
2469  allowing identical post-ISR images to be generated by both the
2470  processCcd and isrTask code.
2471  """
2472  ConfigClass = RunIsrConfig
2473  _DefaultName = "runIsr"
2474 
2475  def __init__(self, *args, **kwargs):
2476  super().__init__(*args, **kwargs)
2477  self.makeSubtask("isr")
2478 
2479  def runDataRef(self, dataRef):
2480  """
2481  Parameters
2482  ----------
2483  dataRef : `lsst.daf.persistence.ButlerDataRef`
2484  data reference of the detector data to be processed
2485 
2486  Returns
2487  -------
2488  result : `pipeBase.Struct`
2489  Result struct with component:
2490 
2491  - exposure : `lsst.afw.image.Exposure`
2492  Post-ISR processed exposure.
2493  """
2494  return self.isr.runDataRef(dataRef)
lsst::ip::isr.isrTask.IsrTaskConfig.maskListToInterpolate
maskListToInterpolate
Definition: isrTask.py:683
lsst::ip::isr.isrTask.FakeAmp.getRawBBox
def getRawBBox(self)
Definition: isrTask.py:2438
lsst::afw::image
Backwards-compatibility support for depersisting the old Calib (FluxMag0/FluxMag0Err) objects.
Definition: imageAlgorithm.dox:1
lsst::ip::isr.linearize.Linearizer
Definition: linearize.py:40
lsst::log.log.logContinued.warn
def warn(fmt, *args)
Definition: logContinued.py:202
lsst::meas::algorithms.detection
Definition: detection.py:1
lsst::ip::isr.isrTask.IsrTask.debugView
def debugView(self, exposure, stepname)
Definition: isrTask.py:2393
lsst::ip::isr.isrTask.FakeAmp.getReadNoise
def getReadNoise(self)
Definition: isrTask.py:2447
lsst::log.log.logContinued.info
def info(fmt, *args)
Definition: logContinued.py:198
lsst::ip::isr.isrTask.IsrTask.maskAndInterpolateNan
def maskAndInterpolateNan(self, exposure)
Definition: isrTask.py:2232
lsst::ip::isr.isrTask.IsrTaskConfig.validate
def validate(self)
Definition: isrTask.py:775
lsst::afw::image::Mask
Represent a 2-dimensional array of bitmask pixels.
Definition: Mask.h:77
lsst::meas::algorithms.defects.Defects
Definition: defects.py:55
lsst::ip::isr.isrTask.IsrTaskConnections
Definition: isrTask.py:63
lsst::ip::isr.isrTask.FakeAmp.getRawHorizontalOverscanBBox
def getRawHorizontalOverscanBBox(self)
Definition: isrTask.py:2441
lsst::afw.display.ds9.getDisplay
getDisplay
Definition: ds9.py:34
lsst::ip::isr.isrTask.FakeAmp._RawHorizontalOverscanBBox
_RawHorizontalOverscanBBox
Definition: isrTask.py:2430
lsst::ip::isr.isrTask.IsrTask.measureBackground
def measureBackground(self, exposure, IsrQaConfig=None)
Definition: isrTask.py:2252
lsst::afw::image::Exposure
A class to contain the data, WCS, and other information needed to describe an image of the sky.
Definition: Exposure.h:72
lsst::ip::isr.isrTask.IsrTask.maskAmplifier
def maskAmplifier(self, ccdExposure, amp, defects)
Definition: isrTask.py:1711
lsst::ip::isr.isrTask.IsrTask.roughZeroPoint
def roughZeroPoint(self, exposure)
Definition: isrTask.py:2316
lsst::afw::image::Filter
Holds an integer identifier for an LSST filter.
Definition: Filter.h:141
lsst::ip::isr.isrTask.IsrTask.saturationDetection
def saturationDetection(self, exposure, amp)
Detect saturated pixels and mask them using mask plane config.saturatedMaskName, in place.
Definition: isrTask.py:2048
lsst.gdb.afw.printers.debug
bool debug
Definition: printers.py:9
lsst::afw::image::makeExposure
std::shared_ptr< Exposure< ImagePixelT, MaskPixelT, VariancePixelT > > makeExposure(MaskedImage< ImagePixelT, MaskPixelT, VariancePixelT > &mimage, std::shared_ptr< geom::SkyWcs const > wcs=std::shared_ptr< geom::SkyWcs const >())
A function to return an Exposure of the correct type (cf.
Definition: Exposure.h:442
lsst::ip::isr.isrTask.IsrTask.maskNan
def maskNan(self, exposure)
Definition: isrTask.py:2206
lsstDebug.getDebugFrame
def getDebugFrame(debugDisplay, name)
Definition: lsstDebug.py:90
pex.config.history.format
def format(config, name=None, writeSourceLine=True, prefix="", verbose=False)
Definition: history.py:174
lsst::ip::isr.isrTask.IsrTaskConfig.doNanInterpolation
doNanInterpolation
Definition: isrTask.py:670
lsst::ip::isr.isrTask.IsrTask.runDataRef
def runDataRef(self, sensorRef)
Definition: isrTask.py:1527
lsst::afw.display
Definition: __init__.py:1
lsst::afw::math::makeStatistics
Statistics makeStatistics(lsst::afw::math::MaskedVector< EntryT > const &mv, std::vector< WeightPixel > const &vweights, int const flags, StatisticsControl const &sctrl=StatisticsControl())
The makeStatistics() overload to handle lsst::afw::math::MaskedVector<>
Definition: Statistics.h:520
lsst::ip::isr.isrTask.IsrTask.suspectDetection
def suspectDetection(self, exposure, amp)
Detect suspect pixels and mask them using mask plane config.suspectMaskName, in place.
Definition: isrTask.py:2097
lsst::afw::image::makePhotoCalibFromCalibZeroPoint
std::shared_ptr< PhotoCalib > makePhotoCalibFromCalibZeroPoint(double instFluxMag0, double instFluxMag0Err)
Construct a PhotoCalib from the deprecated Calib-style instFluxMag0/instFluxMag0Err values.
Definition: PhotoCalib.cc:614
lsst::ip::isr.isrTask.IsrTask.maskEdges
def maskEdges(self, exposure, numEdgePixels=0, maskPlane="SUSPECT")
Mask edge pixels with applicable mask plane.
Definition: isrTask.py:2155
lsst::ip::isr.isrTask.FakeAmp.getBBox
def getBBox(self)
Definition: isrTask.py:2435
lsst::ip::isr.isrTask.FakeAmp
Definition: isrTask.py:2415
lsst::ip::isr.isrTask.IsrTask.maskDefect
def maskDefect(self, exposure, defectBaseList)
Mask defects using mask plane "BAD", in place.
Definition: isrTask.py:2132
lsst::ip::isr.isrTask.FakeAmp._readNoise
_readNoise
Definition: isrTask.py:2432
lsst::ip::isr.isrTask.IsrTaskConfig
Definition: isrTask.py:228
lsst::ip::isr.isrTask.IsrTask.flatCorrection
def flatCorrection(self, exposure, flatExposure, invert=False)
Apply flat correction in place.
Definition: isrTask.py:2023
lsst::ip::isr.isrQa.IsrQaConfig
Definition: isrQa.py:58
lsst::ip::isr.isrTask.FakeAmp.getSuspectLevel
def getSuspectLevel(self)
Definition: isrTask.py:2453
lsst::ip::isr.isrTask.IsrTask
Definition: isrTask.py:785
lsst::ip::isr.isrTask.IsrTask.overscanCorrection
def overscanCorrection(self, ccdExposure, amp)
Definition: isrTask.py:1782
lsst::ip::isr.isrTask.RunIsrTask
Definition: isrTask.py:2461
lsst::ip::isr.isrTask.IsrTask.doLinearize
def doLinearize(self, detector)
Check if linearization is needed for the detector cameraGeom.
Definition: isrTask.py:2004
lsst::afw::geom::polygon::Polygon
Cartesian polygons.
Definition: Polygon.h:59
lsst::ip::isr.isrTask.IsrTaskConfig.doFlat
doFlat
Definition: isrTask.py:579
lsst::ip::isr.isrTask.IsrTask.__init__
def __init__(self, **kwargs)
Definition: isrTask.py:816
astshim.fitsChanContinued.contains
def contains(self, name)
Definition: fitsChanContinued.py:127
lsst::ip::isr::maskNans
size_t maskNans(afw::image::MaskedImage< PixelT > const &mi, afw::image::MaskPixel maskVal, afw::image::MaskPixel allow=0)
Mask NANs in an image.
Definition: Isr.cc:35
lsstDebug.getInfo
getInfo
Definition: lsstDebug.py:87
lsst::ip::isr.isrTask.IsrTask.darkCorrection
def darkCorrection(self, exposure, darkExposure, invert=False)
Apply dark correction in place.
Definition: isrTask.py:1961
lsst::ip::isr.isrTask.IsrTask.readIsrData
def readIsrData(self, dataRef, rawExposure)
Retrieve necessary frames for instrument signature removal.
Definition: isrTask.py:900
lsst::ip::isr.isrTask.FakeAmp._gain
_gain
Definition: isrTask.py:2431
lsst::ip::isr.isrTask.IsrTask.convertIntToFloat
def convertIntToFloat(self, exposure)
Definition: isrTask.py:1674
lsst::ip::isr.isrTask.IsrTask.flatContext
def flatContext(self, exp, flat, dark=None)
Definition: isrTask.py:2363
lsst::ip::isr.isrTask.RunIsrTask.__init__
def __init__(self, *args, **kwargs)
Definition: isrTask.py:2475
lsst::ip::isr.isrTask.IsrTaskConfig.doSaturationInterpolation
doSaturationInterpolation
Definition: isrTask.py:663
lsst::ip::isr.isrTask.FakeAmp._bbox
_bbox
Definition: isrTask.py:2429
lsst::daf::persistence.butler
Definition: butler.py:1
lsst::ip::isr.isrTask.IsrTaskConfig.doApplyGains
doApplyGains
Definition: isrTask.py:611
lsst::ip::isr.isrTask.IsrTask.vignettePolygon
vignettePolygon
Definition: isrTask.py:1435
lsst::ip::isr.isrTask.IsrTask.updateVariance
def updateVariance(self, ampExposure, amp, overscanImage=None)
Definition: isrTask.py:1911
lsst::ip::isr.isrTask.IsrTask.maskAndInterpolateDefects
def maskAndInterpolateDefects(self, exposure, defectBaseList)
Definition: isrTask.py:2181
lsst::ip::isr.isrTask.IsrTask.run
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:1064
object
lsst::ip::isr.isrTask.FakeAmp._saturation
_saturation
Definition: isrTask.py:2433
lsst::ip::isr.isrTask.FakeAmp.__init__
def __init__(self, exposure, config)
Definition: isrTask.py:2428
lsst::ip::isr.isrTask.RunIsrTask.runDataRef
def runDataRef(self, dataRef)
Definition: isrTask.py:2479
lsst::ip::isr.isrTask.FakeAmp.getGain
def getGain(self)
Definition: isrTask.py:2444
lsst::afw::math::StatisticsControl
Pass parameters to a Statistics object.
Definition: Statistics.h:93
lsst::geom
Definition: geomOperators.dox:4
lsst::ip::isr.isrTask.IsrTask.ensureExposure
def ensureExposure(self, inputExp, camera, detectorNum)
Definition: isrTask.py:1626
lsst::afw::cameraGeom
Definition: Amplifier.h:33
lsst::ip::isr.isrTask.IsrTask.setValidPolygonIntersect
def setValidPolygonIntersect(self, ccdExposure, fpPolygon)
Set the valid polygon as the intersection of fpPolygon and the ccd corners.
Definition: isrTask.py:2339
list
daf::base::PropertyList * list
Definition: fits.cc:913
type
table::Key< int > type
Definition: Detector.cc:163
lsst::afw::math
Definition: statistics.dox:6
lsst::geom::Point< int, 2 >
lsst::geom::Box2I
An integer coordinate rectangle.
Definition: Box.h:55
lsst::ip::isr.isrTask.FakeAmp.getSaturation
def getSaturation(self)
Definition: isrTask.py:2450
lsst::ip::isr.isrTask.IsrTask.saturationInterpolation
def saturationInterpolation(self, exposure)
Interpolate over saturated pixels, in place.
Definition: isrTask.py:2072
lsst::afw::image::makeMaskedImage
MaskedImage< ImagePixelT, MaskPixelT, VariancePixelT > * makeMaskedImage(typename std::shared_ptr< Image< ImagePixelT >> image, typename std::shared_ptr< Mask< MaskPixelT >> mask=Mask< MaskPixelT >(), typename std::shared_ptr< Image< VariancePixelT >> variance=Image< VariancePixelT >())
A function to return a MaskedImage of the correct type (cf.
Definition: MaskedImage.h:1279
lsst::daf::persistence
Definition: Utils.h:50
lsst::ip::isr.isrTask.IsrTask.runQuantum
def runQuantum(self, butlerQC, inputRefs, outputRefs)
Definition: isrTask.py:826
lsst::ip::isr.isrTask.IsrTask.getIsrExposure
def getIsrExposure(self, dataRef, datasetType, dateObs=None, immediate=True)
Retrieve a calibration dataset for removing instrument signature.
Definition: isrTask.py:1578
lsst.pipe.base
Definition: __init__.py:1
lsst::meas::algorithms
Fit spatial kernel using approximate fluxes for candidates, and solving a linear system of equations.
Definition: CoaddBoundedField.h:34
lsst::ip::isr.isrTask.IsrTaskConnections.__init__
def __init__(self, *config=None)
Definition: isrTask.py:180
lsst::geom::Extent< int, 2 >
set
daf::base::PropertySet * set
Definition: fits.cc:912
lsst.pipe.base.connectionTypes
Definition: connectionTypes.py:1
lsst::afw::geom
Definition: frameSetUtils.h:40
lsst::ip::isr.isrTask.RunIsrConfig
Definition: isrTask.py:2457