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