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