LSST Applications g0f08755f38+82efc23009,g12f32b3c4e+e7bdf1200e,g1653933729+a8ce1bb630,g1a0ca8cf93+50eff2b06f,g28da252d5a+52db39f6a5,g2bbee38e9b+37c5a29d61,g2bc492864f+37c5a29d61,g2cdde0e794+c05ff076ad,g3156d2b45e+41e33cbcdc,g347aa1857d+37c5a29d61,g35bb328faa+a8ce1bb630,g3a166c0a6a+37c5a29d61,g3e281a1b8c+fb992f5633,g414038480c+7f03dfc1b0,g41af890bb2+11b950c980,g5fbc88fb19+17cd334064,g6b1c1869cb+12dd639c9a,g781aacb6e4+a8ce1bb630,g80478fca09+72e9651da0,g82479be7b0+04c31367b4,g858d7b2824+82efc23009,g9125e01d80+a8ce1bb630,g9726552aa6+8047e3811d,ga5288a1d22+e532dc0a0b,gae0086650b+a8ce1bb630,gb58c049af0+d64f4d3760,gc28159a63d+37c5a29d61,gcf0d15dbbd+2acd6d4d48,gd7358e8bfb+778a810b6e,gda3e153d99+82efc23009,gda6a2b7d83+2acd6d4d48,gdaeeff99f8+1711a396fd,ge2409df99d+6b12de1076,ge79ae78c31+37c5a29d61,gf0baf85859+d0a5978c5a,gf3967379c6+4954f8c433,gfb92a5be7c+82efc23009,gfec2e1e490+2aaed99252,w.2024.46
LSST Data Management Base Package
Loading...
Searching...
No Matches
overscan.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__all__ = ["OverscanCorrectionTaskConfig", "OverscanCorrectionTask",
23 "SerialOverscanCorrectionTaskConfig", "SerialOverscanCorrectionTask",
24 "ParallelOverscanCorrectionTaskConfig", "ParallelOverscanCorrectionTask",
25 ]
26
27import numpy as np
28import lsst.afw.math as afwMath
29import lsst.afw.image as afwImage
30import lsst.geom as geom
31import lsst.pipe.base as pipeBase
32import lsst.pex.config as pexConfig
33
34from .isr import fitOverscanImage, fitOverscanImageMean
35from .isrFunctions import makeThresholdMask
36
37
38class OverscanCorrectionTaskConfigBase(pexConfig.Config):
39 """Overscan correction options.
40 """
41 fitType = pexConfig.ChoiceField(
42 dtype=str,
43 doc="The method for fitting the overscan bias level.",
44 default='MEDIAN',
45 allowed={
46 "POLY": "Fit ordinary polynomial to the longest axis of the overscan region",
47 "CHEB": "Fit Chebyshev polynomial to the longest axis of the overscan region",
48 "LEG": "Fit Legendre polynomial to the longest axis of the overscan region",
49 "NATURAL_SPLINE": "Fit natural spline to the longest axis of the overscan region",
50 "CUBIC_SPLINE": "Fit cubic spline to the longest axis of the overscan region",
51 "AKIMA_SPLINE": "Fit Akima spline to the longest axis of the overscan region",
52 "MEAN": "Correct using the mean of the overscan region",
53 "MEANCLIP": "Correct using a clipped mean of the overscan region",
54 "MEDIAN": "Correct using the median of the overscan region",
55 "MEDIAN_PER_ROW": "Correct using the median per row of the overscan region",
56 "MEAN_PER_ROW": "Correct using the mean per row of the overscan region",
57 },
58 )
59 order = pexConfig.Field(
60 dtype=int,
61 doc=("Order of polynomial to fit if overscan fit type is a polynomial, "
62 "or number of spline knots if overscan fit type is a spline."),
63 default=1,
64 )
65 numSigmaClip = pexConfig.Field(
66 dtype=float,
67 doc="Rejection threshold (sigma) for collapsing overscan before fit",
68 default=3.0,
69 )
70 maskPlanes = pexConfig.ListField(
71 dtype=str,
72 doc="Mask planes to reject when measuring overscan",
73 default=['BAD', 'SAT'],
74 )
75 overscanIsInt = pexConfig.Field(
76 dtype=bool,
77 doc="Treat overscan as an integer image for purposes of fitType=MEDIAN"
78 " and fitType=MEDIAN_PER_ROW.",
79 default=True,
80 )
81 maxDeviation = pexConfig.Field(
82 dtype=float,
83 doc="Maximum deviation from median (in ADU) to mask in overscan correction.",
84 default=1000.0, check=lambda x: x > 0,
85 )
86
87
89 doParallelOverscan = pexConfig.Field(
90 dtype=bool,
91 doc="Correct using parallel overscan after serial overscan correction?",
92 default=False,
93 )
94 parallelOverscanMaskThreshold = pexConfig.Field(
95 dtype=int,
96 doc="Threshold above which pixels in the parallel overscan are masked as bleeds.",
97 default=100000,
98 )
99 parallelOverscanMaskGrowSize = pexConfig.Field(
100 dtype=int,
101 doc="Masks created from saturated bleeds should be grown by this many "
102 "pixels during construction of the parallel overscan mask. "
103 "This value determined from the ITL chip in the LATISS camera",
104 default=7,
105 )
106 leadingColumnsToSkip = pexConfig.Field(
107 dtype=int,
108 doc="Number of leading columns to skip in serial overscan correction.",
109 default=0,
110 )
111 trailingColumnsToSkip = pexConfig.Field(
112 dtype=int,
113 doc="Number of trailing columns to skip in serial overscan correction.",
114 default=0,
115 )
116 leadingRowsToSkip = pexConfig.Field(
117 dtype=int,
118 doc="Number of leading rows to skip in parallel overscan correction.",
119 default=0,
120 )
121 trailingRowsToSkip = pexConfig.Field(
122 dtype=int,
123 doc="Number of trailing rows to skip in parallel overscan correction.",
124 default=0,
125 )
126
127
128class OverscanCorrectionTaskBase(pipeBase.Task):
129 """Base Correction task for overscan.
130
131 This class contains a number of utilities that are easier to
132 understand and use when they are not embedded in nested if/else
133 loops.
134
135 Parameters
136 ----------
137 statControl : `lsst.afw.math.StatisticsControl`, optional
138 Statistics control object.
139 """
140 ConfigClass = OverscanCorrectionTaskConfigBase
141 _DefaultName = "overscanBase"
142
143 def __init__(self, statControl=None, **kwargs):
144 super().__init__(**kwargs)
145 self.allowDebug = True
146
147 if statControl:
148 self.statControl = statControl
149 else:
151 self.statControl.setNumSigmaClip(self.config.numSigmaClip)
152 self.statControl.setAndMask(afwImage.Mask.getPlaneBitMask(self.config.maskPlanes))
153
154 def run(self, exposure, amp, isTransposed=False):
155 raise NotImplementedError("run method is not defined for OverscanCorrectionTaskBase")
156
157 def correctOverscan(self, exposure, amp, imageBBox, overscanBBox,
158 isTransposed=True, leadingToSkip=0, trailingToSkip=0,
159 overscanFraction=1.0, imageThreshold=np.inf,
160 maskedRowColumnGrowSize=0):
161 """Trim the exposure, fit the overscan, subtract the fit, and
162 calculate statistics.
163
164 Parameters
165 ----------
166 exposure : `lsst.afw.image.Exposure`
167 Exposure containing the data.
168 amp : `lsst.afw.cameraGeom.Amplifier`
169 The amplifier that is to be corrected.
170 imageBBox: `lsst.geom.Box2I`
171 Bounding box of the image data that will have the overscan
172 subtracted. If parallel overscan will be performed, that
173 area is added to the image bounding box during serial
174 overscan correction.
175 overscanBBox: `lsst.geom.Box2I`
176 Bounding box for the overscan data.
177 isTransposed: `bool`
178 If true, then the data will be transposed before fitting
179 the overscan.
180 leadingToSkip : `int`, optional
181 Leading rows/columns to skip.
182 trailingToSkip : `int`, optional
183 Leading rows/columns to skip.
184 overscanFraction : `float`, optional
185 If the overscan region median is greater than overscanFraction
186 and the imaging region median is greater than imageThreshold
187 then overscan correction will be skipped.
188 maxLevel : `float`, optional
189 If the overscan region median is greater than overscanFraction
190 and the imaging region median is greater than imageThreshold
191 then overscan correction will be skipped.
192 maskedRowColumnGrowSize : `int`, optional
193 If a column (parallel overscan) or row (serial overscan) is
194 completely masked, then grow the mask by this radius. If the
195 value is <=0 then this will not be checked.
196
197 Returns
198 -------
199 results : `lsst.pipe.base.Struct`
200 ``ampOverscanModel``
201 Overscan model broadcast to the full image size.
202 (`lsst.afw.image.Exposure`)
203 ``overscanOverscanModel``
204 Overscan model broadcast to the full overscan image
205 size. (`lsst.afw.image.Exposure`)
206 ``overscanImage``
207 Overscan image with the overscan fit subtracted.
208 (`lsst.afw.image.Exposure`)
209 ``overscanValue``
210 Overscan model. (`float` or `np.array`)
211 ``overscanMean``
212 Mean value of the overscan fit. (`float`)
213 ``overscanMedian``
214 Median value of the overscan fit. (`float`)
215 ``overscanSigma``
216 Standard deviation of the overscan fit. (`float`)
217 ``overscanMeanResidual``
218 Mean value of the overscan region after overscan
219 subtraction. (`float`)
220 ``overscanMedianResidual``
221 Median value of the overscan region after overscan
222 subtraction. (`float`)
223 ``overscanSigmaResidual``
224 Standard deviation of the overscan region after
225 overscan subtraction. (`float`)
226 """
227 overscanBox = self.trimOverscan(exposure, amp, overscanBBox,
228 leadingToSkip,
229 trailingToSkip,
230 transpose=isTransposed)
231 overscanImage = exposure[overscanBox].getMaskedImage()
232 overscanArray = overscanImage.image.array
233
234 # Mask pixels.
235 maskVal = overscanImage.mask.getPlaneBitMask(self.config.maskPlanes)
236 overscanMask = ~((overscanImage.mask.array & maskVal) == 0)
237
238 badResults = False
239 overscanMedian = np.nanmedian(overscanImage.image.array)
240 imageMedian = np.nanmedian(exposure[imageBBox].image.array)
241
242 if np.all(overscanMask):
243 self.log.warning(
244 "All overscan pixels masked when attempting overscan correction for %s",
245 amp.getName(),
246 )
247 badResults = True
248 elif overscanMedian/imageMedian > overscanFraction and imageMedian > imageThreshold:
249 self.log.warning(
250 "The level in the overscan region (%.2f) compared to the image region (%.2f) is "
251 "greater than the maximum fraction (%.2f) for %s",
252 overscanMedian,
253 imageMedian,
254 overscanFraction,
255 amp.getName(),
256 )
257 badResults = True
258
259 if badResults:
260 # Do not do overscan subtraction at all.
261 overscanResults = pipeBase.Struct(
262 overscanValue=0.0,
263 overscanMean=0.0,
264 overscanMedian=0.0,
265 overscanSigma=0.0,
266 )
267 else:
268 median = np.ma.median(np.ma.masked_where(overscanMask, overscanArray))
269 bad = np.where((np.abs(overscanArray - median) > self.config.maxDeviation) & (~overscanMask))
270 # Mark the bad pixels as BAD.
271 overscanImage.mask.array[bad] = overscanImage.mask.getPlaneBitMask("BAD")
272
273 # Check for completely masked row/column.
274 if maskedRowColumnGrowSize > 0 and len(bad) > 0:
275 if isTransposed:
276 axis = 0
277 nComp = overscanArray.shape[0]
278 else:
279 axis = 1
280 nComp = overscanArray.shape[1]
281
282 # We only need to look at the bad pixels set here for this
283 # mask growth.
284 overscanMaskTemp = np.zeros_like(overscanMask)
285 overscanMaskTemp[bad] = True
286
287 nMaskedArray = np.sum(overscanMaskTemp, axis=axis, dtype=np.int32)
288 badRowsColumns, = np.where(nMaskedArray == nComp)
289 if len(badRowsColumns) > 0:
290 dataView = afwImage.MaskedImageF(exposure.maskedImage,
291 overscanBBox,
292 afwImage.PARENT)
293 if isTransposed:
294 pixelsCopy = dataView.image.array[:, badRowsColumns].copy()
295 dataView.image.array[:, badRowsColumns] = 1e30
296 else:
297 pixelsCopy = dataView.image.array[badRowsColumns, :].copy()
298 dataView.image.array[badRowsColumns, :] = 1e30
299
300 makeThresholdMask(
301 maskedImage=dataView,
302 threshold=1e30,
303 growFootprints=maskedRowColumnGrowSize,
304 maskName="BAD",
305 )
306
307 if isTransposed:
308 dataView.image.array[:, badRowsColumns] = pixelsCopy
309 else:
310 dataView.image.array[badRowsColumns, :] = 1e30
311
312 # Do overscan fit.
313 # CZW: Handle transposed correctly.
314 overscanResults = self.fitOverscan(overscanImage, isTransposed=isTransposed)
315
316 # Correct image region (and possibly parallel-overscan region).
317 ampImage = exposure[imageBBox]
318 ampOverscanModel = self.broadcastFitToImage(overscanResults.overscanValue,
319 ampImage.image.array,
320 transpose=isTransposed)
321 ampImage.image.array -= ampOverscanModel
322
323 # Correct overscan region (and possibly doubly-overscaned
324 # region).
325 overscanImage = exposure[overscanBBox]
326 # CZW: Transposed?
327 overscanOverscanModel = self.broadcastFitToImage(overscanResults.overscanValue,
328 overscanImage.image.array)
329 self.debugView(overscanImage, overscanResults.overscanValue, amp, isTransposed=isTransposed)
330 overscanImage.image.array -= overscanOverscanModel
331
332 # Find residual fit statistics.
333 stats = afwMath.makeStatistics(overscanImage.getMaskedImage(),
334 afwMath.MEAN | afwMath.MEDIAN | afwMath.STDEVCLIP, self.statControl)
335 residualMean = stats.getValue(afwMath.MEAN)
336 residualMedian = stats.getValue(afwMath.MEDIAN)
337 residualSigma = stats.getValue(afwMath.STDEVCLIP)
338
339 return pipeBase.Struct(ampOverscanModel=ampOverscanModel,
340 overscanOverscanModel=overscanOverscanModel,
341 overscanImage=overscanImage,
342 overscanValue=overscanResults.overscanValue,
343
344 overscanMean=overscanResults.overscanMean,
345 overscanMedian=overscanResults.overscanMedian,
346 overscanSigma=overscanResults.overscanSigma,
347 overscanMeanResidual=residualMean,
348 overscanMedianResidual=residualMedian,
349 overscanSigmaResidual=residualSigma
350 )
351
352 def broadcastFitToImage(self, overscanValue, imageArray, transpose=False):
353 """Broadcast 0 or 1 dimension fit to appropriate shape.
354
355 Parameters
356 ----------
357 overscanValue : `numpy.ndarray`, (Nrows, ) or scalar
358 Overscan fit to broadcast.
359 imageArray : `numpy.ndarray`, (Nrows, Ncols)
360 Image array that we want to match.
361 transpose : `bool`, optional
362 Switch order to broadcast along the other axis.
363
364 Returns
365 -------
366 overscanModel : `numpy.ndarray`, (Nrows, Ncols) or scalar
367 Expanded overscan fit.
368
369 Raises
370 ------
371 RuntimeError
372 Raised if no axis has the appropriate dimension.
373 """
374 if isinstance(overscanValue, np.ndarray):
375 overscanModel = np.zeros_like(imageArray)
376
377 if transpose is False:
378 if imageArray.shape[0] == overscanValue.shape[0]:
379 overscanModel[:, :] = overscanValue[:, np.newaxis]
380 elif imageArray.shape[1] == overscanValue.shape[0]:
381 overscanModel[:, :] = overscanValue[np.newaxis, :]
382 elif imageArray.shape[0] == overscanValue.shape[1]:
383 overscanModel[:, :] = overscanValue[np.newaxis, :]
384 else:
385 raise RuntimeError(f"Could not broadcast {overscanValue.shape} to "
386 f"match {imageArray.shape}")
387 else:
388 if imageArray.shape[1] == overscanValue.shape[0]:
389 overscanModel[:, :] = overscanValue[np.newaxis, :]
390 elif imageArray.shape[0] == overscanValue.shape[0]:
391 overscanModel[:, :] = overscanValue[:, np.newaxis]
392 elif imageArray.shape[1] == overscanValue.shape[1]:
393 overscanModel[:, :] = overscanValue[:, np.newaxis]
394 else:
395 raise RuntimeError(f"Could not broadcast {overscanValue.shape} to "
396 f"match {imageArray.shape}")
397 else:
398 overscanModel = overscanValue
399
400 return overscanModel
401
402 def trimOverscan(self, exposure, amp, bbox, skipLeading, skipTrailing, transpose=False):
403 """Trim overscan region to remove edges.
404
405 Parameters
406 ----------
407 exposure : `lsst.afw.image.Exposure`
408 Exposure containing data.
409 amp : `lsst.afw.cameraGeom.Amplifier`
410 Amplifier containing geometry information.
411 bbox : `lsst.geom.Box2I`
412 Bounding box of the overscan region.
413 skipLeading : `int`
414 Number of leading (towards data region) rows/columns to skip.
415 skipTrailing : `int`
416 Number of trailing (away from data region) rows/columns to skip.
417 transpose : `bool`, optional
418 Operate on the transposed array.
419
420 Returns
421 -------
422 overscanArray : `numpy.array`, (N, M)
423 Data array to fit.
424 overscanMask : `numpy.array`, (N, M)
425 Data mask.
426 """
427 dx0, dy0, dx1, dy1 = (0, 0, 0, 0)
428 dataBBox = amp.getRawDataBBox()
429 if transpose:
430 if dataBBox.getBeginY() < bbox.getBeginY():
431 dy0 += skipLeading
432 dy1 -= skipTrailing
433 else:
434 dy0 += skipTrailing
435 dy1 -= skipLeading
436 else:
437 if dataBBox.getBeginX() < bbox.getBeginX():
438 dx0 += skipLeading
439 dx1 -= skipTrailing
440 else:
441 dx0 += skipTrailing
442 dx1 -= skipLeading
443
444 overscanBBox = geom.Box2I(bbox.getBegin() + geom.Extent2I(dx0, dy0),
445 geom.Extent2I(bbox.getWidth() - dx0 + dx1,
446 bbox.getHeight() - dy0 + dy1))
447 return overscanBBox
448
449 def fitOverscan(self, overscanImage, isTransposed=False):
450 if self.config.fitType in ('MEAN', 'MEANCLIP', 'MEDIAN'):
451 # Transposition has no effect here.
452 overscanResult = self.measureConstantOverscan(overscanImage)
453 overscanValue = overscanResult.overscanValue
454 overscanMean = overscanValue
455 overscanMedian = overscanValue
456 overscanSigma = 0.0
457 elif self.config.fitType in ('MEDIAN_PER_ROW', 'MEAN_PER_ROW', 'POLY', 'CHEB', 'LEG',
458 'NATURAL_SPLINE', 'CUBIC_SPLINE', 'AKIMA_SPLINE'):
459 # Force transposes as needed
460 overscanResult = self.measureVectorOverscan(overscanImage, isTransposed)
461 overscanValue = overscanResult.overscanValue
462
463 stats = afwMath.makeStatistics(overscanResult.overscanValue,
464 afwMath.MEAN | afwMath.MEDIAN | afwMath.STDEVCLIP,
465 self.statControl)
466 overscanMean = stats.getValue(afwMath.MEAN)
467 overscanMedian = stats.getValue(afwMath.MEDIAN)
468 overscanSigma = stats.getValue(afwMath.STDEVCLIP)
469 else:
470 raise ValueError('%s : %s an invalid overscan type' %
471 ("overscanCorrection", self.config.fitType))
472
473 return pipeBase.Struct(overscanValue=overscanValue,
474 overscanMean=overscanMean,
475 overscanMedian=overscanMedian,
476 overscanSigma=overscanSigma,
477 )
478
479 def maskParallelOverscan(self, exposure, detector):
480 """Mask the union of high values on all amplifiers in the parallel
481 overscan.
482
483 This operates on the image in-place.
484
485 Parameters
486 ----------
487 exposure : `lsst.afw.image.Exposure`
488 An untrimmed raw exposure.
489 detector : `lsst.afw.cameraGeom.Detector`
490 The detetor to use for amplifier geometry.
491 """
492 parallelMask = None
493
494 for amp in detector:
495 dataView = afwImage.MaskedImageF(exposure.getMaskedImage(),
496 amp.getRawParallelOverscanBBox(),
497 afwImage.PARENT)
498 # This should mark all the saturated pixels as SAT.
499 makeThresholdMask(
500 maskedImage=dataView,
501 threshold=self.config.parallelOverscanMaskThreshold,
502 growFootprints=self.config.parallelOverscanMaskGrowSize,
503 maskName="SAT"
504 )
505 if parallelMask is None:
506 parallelMask = dataView.mask.array
507 else:
508 parallelMask |= dataView.mask.array
509 for amp in detector:
510 dataView = afwImage.MaskedImageF(exposure.getMaskedImage(),
511 amp.getRawParallelOverscanBBox(),
512 afwImage.PARENT)
513 dataView.mask.array |= parallelMask
514
515 # Constant methods
516 def measureConstantOverscan(self, image):
517 """Measure a constant overscan value.
518
519 Parameters
520 ----------
521 image : `lsst.afw.image.Image` or `lsst.afw.image.MaskedImage`
522 Image data to measure the overscan from.
523
524 Returns
525 -------
526 results : `lsst.pipe.base.Struct`
527 Overscan result with entries:
528 - ``overscanValue``: Overscan value to subtract (`float`)
529 - ``isTransposed``: Orientation of the overscan (`bool`)
530 """
531 fitType = afwMath.stringToStatisticsProperty(self.config.fitType)
532 overscanValue = afwMath.makeStatistics(image, fitType, self.statControl).getValue()
533
534 return pipeBase.Struct(overscanValue=overscanValue,
535 isTransposed=False)
536
537 # Vector correction utilities
538 def getImageArray(self, image):
539 """Extract the numpy array from the input image.
540
541 Parameters
542 ----------
543 image : `lsst.afw.image.Image` or `lsst.afw.image.MaskedImage`
544 Image data to pull array from.
545
546 calcImage : `numpy.ndarray`
547 Image data array for numpy operating.
548 """
549 if hasattr(image, "getImage"):
550 calcImage = image.getImage().getArray()
551 calcImage = np.ma.masked_where(image.getMask().getArray() & self.statControl.getAndMask(),
552 calcImage)
553 else:
554 calcImage = image.getArray()
555 return calcImage
556
557 def maskOutliers(self, imageArray):
558 """Mask outliers in a row of overscan data from a robust sigma
559 clipping procedure.
560
561 Parameters
562 ----------
563 imageArray : `numpy.ndarray`
564 Image to filter along numpy axis=1.
565
566 Returns
567 -------
568 maskedArray : `numpy.ma.masked_array`
569 Masked image marking outliers.
570 """
571 lq, median, uq = np.percentile(np.ma.getdata(imageArray),
572 [25.0, 50.0, 75.0], axis=1)
573 axisMedians = median
574 axisStdev = 0.74*(uq - lq) # robust stdev
575
576 # Replace pixels that have excessively large stdev values
577 # with the median of stdev values. A large stdev likely
578 # indicates a bleed is spilling into the overscan.
579 axisStdev = np.where(axisStdev > 2.0 * np.median(axisStdev),
580 np.median(axisStdev), axisStdev)
581
582 # Mask pixels that are N-sigma away from their array medians.
583 diff = np.abs(imageArray - axisMedians[:, np.newaxis])
584 masked = np.ma.masked_where(diff > self.statControl.getNumSigmaClip()
585 * axisStdev[:, np.newaxis], imageArray)
586
587 return masked
588
589 def fillMaskedPixels(self, overscanVector):
590 """Fill masked/NaN pixels in the overscan.
591
592 Parameters
593 ----------
594 overscanVector : `np.array` or `np.ma.masked_array`
595 Overscan vector to fill.
596
597 Returns
598 -------
599 overscanVector : `np.ma.masked_array`
600 Filled vector.
601
602 Notes
603 -----
604 Each maskSlice is a section of overscan with contiguous masks.
605 Ideally this adds 5 pixels from the left and right of that
606 mask slice, and takes the median of those values to fill the
607 slice. If this isn't possible, the median of all non-masked
608 values is used. The mask is removed for the pixels filled.
609 """
610 workingCopy = overscanVector
611 if not isinstance(overscanVector, np.ma.MaskedArray):
612 workingCopy = np.ma.masked_array(overscanVector,
613 mask=~np.isfinite(overscanVector))
614
615 defaultValue = np.median(workingCopy.data[~workingCopy.mask])
616 for maskSlice in np.ma.clump_masked(workingCopy):
617 neighborhood = []
618 if maskSlice.start > 5:
619 neighborhood.extend(workingCopy[maskSlice.start - 5:maskSlice.start].data)
620 if maskSlice.stop < workingCopy.size - 5:
621 neighborhood.extend(workingCopy[maskSlice.stop:maskSlice.stop+5].data)
622 if len(neighborhood) > 0:
623 workingCopy.data[maskSlice] = np.nanmedian(neighborhood)
624 workingCopy.mask[maskSlice] = False
625 else:
626 workingCopy.data[maskSlice] = defaultValue
627 workingCopy.mask[maskSlice] = False
628 return workingCopy
629
630 def collapseArray(self, maskedArray, fillMasked=True):
631 """Collapse overscan array (and mask) to a 1-D vector of values.
632
633 Parameters
634 ----------
635 maskedArray : `numpy.ma.masked_array`
636 Masked array of input overscan data.
637 fillMasked : `bool`, optional
638 If true, fill any pixels that are masked with a median of
639 neighbors.
640
641 Returns
642 -------
643 collapsed : `numpy.ma.masked_array`
644 Single dimensional overscan data, combined with the mean.
645
646 """
647 collapsed = np.mean(maskedArray, axis=1)
648 if collapsed.mask.sum() > 0 and fillMasked:
649 collapsed = self.fillMaskedPixels(collapsed)
650
651 return collapsed
652
653 def splineFit(self, indices, collapsed, numBins):
654 """Wrapper function to match spline fit API to polynomial fit API.
655
656 Parameters
657 ----------
658 indices : `numpy.ndarray`
659 Locations to evaluate the spline.
660 collapsed : `numpy.ndarray`
661 Collapsed overscan values corresponding to the spline
662 evaluation points.
663 numBins : `int`
664 Number of bins to use in constructing the spline.
665
666 Returns
667 -------
668 interp : `lsst.afw.math.Interpolate`
669 Interpolation object for later evaluation.
670 """
671 if not np.ma.is_masked(collapsed):
672 collapsed.mask = np.array(len(collapsed)*[np.ma.nomask])
673
674 numPerBin, binEdges = np.histogram(indices, bins=numBins,
675 weights=1 - collapsed.mask.astype(int))
676 with np.errstate(invalid="ignore"):
677 values = np.histogram(indices, bins=numBins,
678 weights=collapsed.data*~collapsed.mask)[0]/numPerBin
679 binCenters = np.histogram(indices, bins=numBins,
680 weights=indices*~collapsed.mask)[0]/numPerBin
681
682 if len(binCenters[numPerBin > 0]) < 5:
683 self.log.warning("Cannot do spline fitting for overscan: %s valid points.",
684 len(binCenters[numPerBin > 0]))
685 # Return a scalar value if we have one, otherwise
686 # return zero. This amplifier is hopefully already
687 # masked.
688 if len(values[numPerBin > 0]) != 0:
689 return float(values[numPerBin > 0][0])
690 else:
691 return 0.0
692
693 interp = afwMath.makeInterpolate(binCenters.astype(float)[numPerBin > 0],
694 values.astype(float)[numPerBin > 0],
695 afwMath.stringToInterpStyle(self.config.fitType))
696 return interp
697
698 @staticmethod
699 def splineEval(indices, interp):
700 """Wrapper function to match spline evaluation API to polynomial fit
701 API.
702
703 Parameters
704 ----------
705 indices : `numpy.ndarray`
706 Locations to evaluate the spline.
707 interp : `lsst.afw.math.interpolate`
708 Interpolation object to use.
709
710 Returns
711 -------
712 values : `numpy.ndarray`
713 Evaluated spline values at each index.
714 """
715
716 return interp.interpolate(indices.astype(float))
717
718 @staticmethod
719 def maskExtrapolated(collapsed):
720 """Create mask if edges are extrapolated.
721
722 Parameters
723 ----------
724 collapsed : `numpy.ma.masked_array`
725 Masked array to check the edges of.
726
727 Returns
728 -------
729 maskArray : `numpy.ndarray`
730 Boolean numpy array of pixels to mask.
731 """
732 maskArray = np.full_like(collapsed, False, dtype=bool)
733 if np.ma.is_masked(collapsed):
734 num = len(collapsed)
735 for low in range(num):
736 if not collapsed.mask[low]:
737 break
738 if low > 0:
739 maskArray[:low] = True
740 for high in range(1, num):
741 if not collapsed.mask[-high]:
742 break
743 if high > 1:
744 maskArray[-high:] = True
745 return maskArray
746
747 def measureVectorOverscan(self, image, isTransposed=False):
748 """Calculate the 1-d vector overscan from the input overscan image.
749
750 Parameters
751 ----------
752 image : `lsst.afw.image.MaskedImage`
753 Image containing the overscan data.
754 isTransposed : `bool`
755 If true, the image has been transposed.
756
757 Returns
758 -------
759 results : `lsst.pipe.base.Struct`
760 Overscan result with entries:
761
762 ``overscanValue``
763 Overscan value to subtract (`float`)
764 ``maskArray``
765 List of rows that should be masked as ``SUSPECT`` when the
766 overscan solution is applied. (`list` [ `bool` ])
767 ``isTransposed``
768 Indicates if the overscan data was transposed during
769 calcuation, noting along which axis the overscan should be
770 subtracted. (`bool`)
771 """
772 calcImage = self.getImageArray(image)
773
774 # operate on numpy-arrays from here
775 if isTransposed:
776 calcImage = np.transpose(calcImage)
777 masked = self.maskOutliers(calcImage)
778
779 if self.config.fitType in ('MEDIAN_PER_ROW', "MEAN_PER_ROW"):
780 if self.config.overscanIsInt:
781 mi = afwImage.MaskedImageI(image.getBBox())
782 masked = masked.astype(int)
783 else:
784 mi = image.clone()
785
786 if isTransposed:
787 masked = masked.transpose()
788
789 mi.image.array[:, :] = masked.data[:, :]
790 if bool(masked.mask.shape):
791 mi.mask.array[:, :] = masked.mask[:, :]
792
793 if self.config.fitType == "MEDIAN_PER_ROW":
794 overscanVector = fitOverscanImage(mi, self.config.maskPlanes, isTransposed)
795 else:
796 overscanVector = fitOverscanImageMean(mi, self.config.maskPlanes, isTransposed)
797
798 overscanVector = self.fillMaskedPixels(overscanVector)
799 maskArray = self.maskExtrapolated(overscanVector)
800 else:
801 collapsed = self.collapseArray(masked)
802
803 num = len(collapsed)
804 indices = 2.0*np.arange(num)/float(num) - 1.0
805
806 poly = np.polynomial
807 fitter, evaler = {
808 'POLY': (poly.polynomial.polyfit, poly.polynomial.polyval),
809 'CHEB': (poly.chebyshev.chebfit, poly.chebyshev.chebval),
810 'LEG': (poly.legendre.legfit, poly.legendre.legval),
811 'NATURAL_SPLINE': (self.splineFitsplineFit, self.splineEvalsplineEval),
812 'CUBIC_SPLINE': (self.splineFitsplineFit, self.splineEvalsplineEval),
813 'AKIMA_SPLINE': (self.splineFitsplineFit, self.splineEvalsplineEval)
814 }[self.config.fitType]
815
816 # These are the polynomial coefficients, or an
817 # interpolation object.
818 coeffs = fitter(indices, collapsed, self.config.order)
819
820 if isinstance(coeffs, float):
821 self.log.warning("Using fallback value %f due to fitter failure. Amplifier will be masked.",
822 coeffs)
823 overscanVector = np.full_like(indices, coeffs)
824 maskArray = np.full_like(collapsed, True, dtype=bool)
825 else:
826 # Otherwise we can just use things as normal.
827 overscanVector = evaler(indices, coeffs)
828 maskArray = self.maskExtrapolated(collapsed)
829
830 return pipeBase.Struct(overscanValue=np.array(overscanVector),
831 maskArray=maskArray,
832 isTransposed=isTransposed)
833
834 def debugView(self, image, model, amp=None, isTransposed=True):
835 """Debug display for the final overscan solution.
836
837 Parameters
838 ----------
839 image : `lsst.afw.image.Image`
840 Input image the overscan solution was determined from.
841 model : `numpy.ndarray` or `float`
842 Overscan model determined for the image.
843 amp : `lsst.afw.cameraGeom.Amplifier`, optional
844 Amplifier to extract diagnostic information.
845 isTransposed : `bool`, optional
846 Does the data need to be transposed before display?
847 """
848 import lsstDebug
849 if not lsstDebug.Info(__name__).display:
850 return
851 if not self.allowDebug:
852 return
853
854 calcImage = self.getImageArray(image)
855 # CZW: Check that this is ok
856 if isTransposed:
857 calcImage = np.transpose(calcImage)
858 masked = self.maskOutliers(calcImage)
859 collapsed = self.collapseArray(masked, fillMasked=False)
860
861 num = len(collapsed)
862 indices = 2.0 * np.arange(num)/float(num) - 1.0
863 indices = np.arange(num)
864
865 if np.ma.is_masked(collapsed):
866 collapsedMask = collapsed.mask
867 else:
868 collapsedMask = np.array(num*[np.ma.nomask])
869
870 import matplotlib.pyplot as plot
871 figure = plot.figure(1)
872 figure.clear()
873 axes = figure.add_axes((0.1, 0.1, 0.8, 0.8))
874 axes.plot(indices[~collapsedMask], collapsed[~collapsedMask], 'k+')
875 if collapsedMask.sum() > 0:
876 axes.plot(indices[collapsedMask], collapsed.data[collapsedMask], 'b+')
877 if isinstance(model, np.ndarray):
878 plotModel = model
879 else:
880 plotModel = np.zeros_like(indices)
881 plotModel += model
882
883 axes.plot(indices, plotModel, 'r-')
884 plot.xlabel("position along overscan region")
885 plot.ylabel("pixel value/fit value")
886 if amp:
887 plot.title(f"{amp.getName()} DataX: "
888 f"[{amp.getRawDataBBox().getBeginX()}:{amp.getRawBBox().getEndX()}]"
889 f"OscanX: [{amp.getRawHorizontalOverscanBBox().getBeginX()}:"
890 f"{amp.getRawHorizontalOverscanBBox().getEndX()}] {self.config.fitType}")
891 else:
892 plot.title("No amp supplied.")
893 figure.show()
894 prompt = "Press Enter or c to continue [chp]..."
895 while True:
896 ans = input(prompt).lower()
897 if ans in ("", " ", "c",):
898 break
899 elif ans in ("p", ):
900 import pdb
901 pdb.set_trace()
902 elif ans in ('x', ):
903 self.allowDebug = False
904 break
905 elif ans in ("h", ):
906 print("[h]elp [c]ontinue [p]db e[x]itDebug")
907 plot.close()
908
909
911 """Correction task for serial/parallel overscan.
912
913 (Will be deprecated)
914
915 This class contains a number of utilities that are easier to
916 understand and use when they are not embedded in nested if/else
917 loops.
918
919 Parameters
920 ----------
921 statControl : `lsst.afw.math.StatisticsControl`, optional
922 Statistics control object.
923 """
924 ConfigClass = OverscanCorrectionTaskConfig
925 _DefaultName = "overscan"
926
927 def run(self, exposure, amp, isTransposed=False):
928 """Measure and remove serial/parallel overscan from an amplifier image.
929
930 This will be deprecated.
931
932 Parameters
933 ----------
934 exposure : `lsst.afw.image.Exposure`
935 Image data that will have the overscan corrections applied.
936 amp : `lsst.afw.cameraGeom.Amplifier`
937 Amplifier to use for debugging purposes.
938 isTransposed : `bool`, optional
939 Is the image transposed, such that serial and parallel
940 overscan regions are reversed? Default is False.
941
942 Returns
943 -------
944 overscanResults : `lsst.pipe.base.Struct`
945 Result struct with components:
946
947 ``imageFit``
948 Value or fit subtracted from the amplifier image data
949 (scalar or `lsst.afw.image.Image`).
950 ``overscanFit``
951 Value or fit subtracted from the serial overscan image
952 data (scalar or `lsst.afw.image.Image`).
953 ``overscanImage``
954 Image of the serial overscan region with the serial
955 overscan correction applied
956 (`lsst.afw.image.Image`). This quantity is used to
957 estimate the amplifier read noise empirically.
958 ``parallelOverscanFit``
959 Value or fit subtracted from the parallel overscan
960 image data (scalar, `lsst.afw.image.Image`, or None).
961 ``parallelOverscanImage``
962 Image of the parallel overscan region with the
963 parallel overscan correction applied
964 (`lsst.afw.image.Image` or None).
965 ``overscanMean``
966 Mean of the fit serial overscan region.
967 This and the following values will be tuples of
968 (serial, parallel) if doParallelOverscan=True.
969 ``overscanMedian``
970 Median of the fit serial overscan region.
971 ``overscanSigma``
972 Sigma of the fit serial overscan region.
973 ``residualMean``
974 Mean of the residual of the serial overscan region after
975 correction.
976 ``residualMedian``
977 Median of the residual of the serial overscan region after
978 correction.
979 ``residualSigma``
980 Mean of the residual of the serial overscan region after
981 correction.
982
983
984 Raises
985 ------
986 RuntimeError
987 Raised if an invalid overscan type is set.
988 """
989 # Do Serial overscan first.
990 serialOverscanBBox = amp.getRawSerialOverscanBBox()
991 imageBBox = amp.getRawDataBBox()
992
993 if self.config.doParallelOverscan:
994 # We need to extend the serial overscan BBox to the full
995 # size of the detector.
996 parallelOverscanBBox = amp.getRawParallelOverscanBBox()
997 imageBBox = imageBBox.expandedTo(parallelOverscanBBox)
998
999 if isTransposed:
1000 serialOverscanBBox = geom.Box2I(
1001 geom.Point2I(serialOverscanBBox.getMinX(), imageBBox.getEndY()),
1002 geom.Extent2I(imageBBox.getWidth(), serialOverscanBBox.getHeight()),
1003 )
1004 else:
1005 serialOverscanBBox = geom.Box2I(
1006 geom.Point2I(serialOverscanBBox.getMinX(),
1007 imageBBox.getMinY()),
1008 geom.Extent2I(serialOverscanBBox.getWidth(),
1009 imageBBox.getHeight()),
1010 )
1011
1012 serialResults = self.correctOverscan(
1013 exposure,
1014 amp,
1015 imageBBox,
1016 serialOverscanBBox,
1017 isTransposed=isTransposed,
1018 leadingToSkip=self.config.leadingColumnsToSkip,
1019 trailingToSkip=self.config.trailingColumnsToSkip,
1020 )
1021 overscanMean = serialResults.overscanMean
1022 overscanMedian = serialResults.overscanMedian
1023 overscanSigma = serialResults.overscanSigma
1024 residualMean = serialResults.overscanMeanResidual
1025 residualMedian = serialResults.overscanMedianResidual
1026 residualSigma = serialResults.overscanSigmaResidual
1027
1028 # Do Parallel Overscan
1029 parallelResults = None
1030 if self.config.doParallelOverscan:
1031 # This does not need any extensions, as we'll only
1032 # subtract it from the data region.
1033 parallelOverscanBBox = amp.getRawParallelOverscanBBox()
1034 imageBBox = amp.getRawDataBBox()
1035
1036 # The serial overscan correction has removed some signal
1037 # from the parallel overscan region, but that is largely a
1038 # constant offset. The collapseArray method now attempts
1039 # to fill fully masked columns with the median of
1040 # neighboring values, with a fallback to the median of the
1041 # correction in all other columns. Filling with neighbor
1042 # values ensures that large variations in the parallel
1043 # overscan do not create new outlier points. The
1044 # MEDIAN_PER_ROW method does this filling as a separate
1045 # operation, using the same method.
1046 parallelResults = self.correctOverscan(
1047 exposure,
1048 amp,
1049 imageBBox,
1050 parallelOverscanBBox,
1051 isTransposed=not isTransposed,
1052 leadingToSkip=self.config.leadingRowsToSkip,
1053 trailingToSkip=self.config.trailingRowsToSkip,
1054 )
1055 overscanMean = (overscanMean, parallelResults.overscanMean)
1056 overscanMedian = (overscanMedian, parallelResults.overscanMedian)
1057 overscanSigma = (overscanSigma, parallelResults.overscanSigma)
1058 residualMean = (residualMean, parallelResults.overscanMeanResidual)
1059 residualMedian = (residualMedian, parallelResults.overscanMedianResidual)
1060 residualSigma = (residualSigma, parallelResults.overscanSigmaResidual)
1061
1062 parallelOverscanFit = parallelResults.overscanOverscanModel if parallelResults else None
1063 parallelOverscanImage = parallelResults.overscanImage if parallelResults else None
1064
1065 return pipeBase.Struct(imageFit=serialResults.ampOverscanModel,
1066 overscanFit=serialResults.overscanOverscanModel,
1067 overscanImage=serialResults.overscanImage,
1068
1069 parallelOverscanFit=parallelOverscanFit,
1070 parallelOverscanImage=parallelOverscanImage,
1071 overscanMean=overscanMean,
1072 overscanMedian=overscanMedian,
1073 overscanSigma=overscanSigma,
1074 residualMean=residualMean,
1075 residualMedian=residualMedian,
1076 residualSigma=residualSigma)
1077
1078
1080 leadingToSkip = pexConfig.Field(
1081 dtype=int,
1082 doc="Number of leading values to skip in serial overscan correction.",
1083 default=0,
1084 )
1085 trailingToSkip = pexConfig.Field(
1086 dtype=int,
1087 doc="Number of trailing values to skip in serial overscan correction.",
1088 default=0,
1089 )
1090
1091
1093 """Correction task for serial overscan.
1094
1095 Parameters
1096 ----------
1097 statControl : `lsst.afw.math.StatisticsControl`, optional
1098 Statistics control object.
1099 """
1100 ConfigClass = SerialOverscanCorrectionTaskConfig
1101 _DefaultName = "serialOverscan"
1102
1103 def run(self, exposure, amp, isTransposed=False):
1104 """Measure and remove serial overscan from an amplifier image.
1105
1106 Parameters
1107 ----------
1108 exposure : `lsst.afw.image.Exposure`
1109 Image data that will have the overscan corrections applied.
1110 amp : `lsst.afw.cameraGeom.Amplifier`
1111 Amplifier to use for debugging purposes.
1112 isTransposed : `bool`, optional
1113 Is the image transposed, such that serial and parallel
1114 overscan regions are reversed? Default is False.
1115
1116 Returns
1117 -------
1118 overscanResults : `lsst.pipe.base.Struct`
1119 Result struct with components:
1120
1121 ``imageFit``
1122 Value or fit subtracted from the amplifier image data
1123 (scalar or `lsst.afw.image.Image`).
1124 ``overscanFit``
1125 Value or fit subtracted from the serial overscan image
1126 data (scalar or `lsst.afw.image.Image`).
1127 ``overscanImage``
1128 Image of the serial overscan region with the serial
1129 overscan correction applied
1130 (`lsst.afw.image.Image`). This quantity is used to
1131 estimate the amplifier read noise empirically.
1132 ``overscanMean``
1133 Mean of the fit serial overscan region.
1134 ``overscanMedian``
1135 Median of the fit serial overscan region.
1136 ``overscanSigma``
1137 Sigma of the fit serial overscan region.
1138 ``residualMean``
1139 Mean of the residual of the serial overscan region after
1140 correction.
1141 ``residualMedian``
1142 Median of the residual of the serial overscan region after
1143 correction.
1144 ``residualSigma``
1145 Mean of the residual of the serial overscan region after
1146 correction.
1147
1148 Raises
1149 ------
1150 RuntimeError
1151 Raised if an invalid overscan type is set.
1152 """
1153 serialOverscanBBox = amp.getRawSerialOverscanBBox()
1154 imageBBox = amp.getRawDataBBox()
1155
1156 # We always want to extend the serial overscan bounding box to
1157 # the full size of the detector.
1158 parallelOverscanBBox = amp.getRawParallelOverscanBBox()
1159 imageBBox = imageBBox.expandedTo(parallelOverscanBBox)
1160
1161 if isTransposed:
1162 serialOverscanBBox = geom.Box2I(
1163 geom.Point2I(serialOverscanBBox.getMinX(), imageBBox.getEndY()),
1164 geom.Extent2I(imageBBox.getWidth(), serialOverscanBBox.getHeight()),
1165 )
1166 else:
1167 serialOverscanBBox = geom.Box2I(
1168 geom.Point2I(serialOverscanBBox.getMinX(),
1169 imageBBox.getMinY()),
1170 geom.Extent2I(serialOverscanBBox.getWidth(),
1171 imageBBox.getHeight()),
1172 )
1173
1174 results = self.correctOverscan(
1175 exposure,
1176 amp,
1177 imageBBox,
1178 serialOverscanBBox,
1179 isTransposed=isTransposed,
1180 leadingToSkip=self.config.leadingToSkip,
1181 trailingToSkip=self.config.trailingToSkip,
1182 )
1183 overscanMean = results.overscanMean
1184 overscanMedian = results.overscanMedian
1185 overscanSigma = results.overscanSigma
1186 residualMean = results.overscanMeanResidual
1187 residualMedian = results.overscanMedianResidual
1188 residualSigma = results.overscanSigmaResidual
1189
1190 return pipeBase.Struct(
1191 imageFit=results.ampOverscanModel,
1192 overscanFit=results.overscanOverscanModel,
1193 overscanImage=results.overscanImage,
1194 overscanMean=overscanMean,
1195 overscanMedian=overscanMedian,
1196 overscanSigma=overscanSigma,
1197 residualMean=residualMean,
1198 residualMedian=residualMedian,
1199 residualSigma=residualSigma,
1200 )
1201
1202
1204 doParallelOverscanSaturation = pexConfig.Field(
1205 dtype=bool,
1206 doc="Mask saturated pixels in parallel overscan region?",
1207 default=True,
1208 )
1209 parallelOverscanSaturationLevel = pexConfig.Field(
1210 dtype=float,
1211 doc="The saturation level (adu) to use if not specified in call to "
1212 "maskParallelOverscanAmp. This should be low enough to capture "
1213 "all possible amplifiers for defect detection.",
1214 default=20000.,
1215 )
1216 parallelOverscanSaturationLevelAdjustmentFactor = pexConfig.Field(
1217 dtype=float,
1218 doc="The parallel overscan saturation level may be below that of "
1219 "the data region. This factor is applied to the amplifier "
1220 "saturation value when evaluating saturation in the parallel "
1221 "overscan region.",
1222 default=0.75,
1223 )
1224 parallelOverscanMaskGrowSize = pexConfig.Field(
1225 dtype=int,
1226 doc="Masks created from saturated bleeds should be grown by this many "
1227 "pixels during construction of the parallel overscan mask. "
1228 "This value determined from the ITL chip in the LATISS camera",
1229 default=7,
1230 )
1231 parallelOverscanMaskedColumnGrowSize = pexConfig.Field(
1232 dtype=int,
1233 doc="When a full column is masked in the parallel overscan (at less "
1234 "than saturation) the mask should be grown by this many pixels. "
1235 "This value is determined from ITL chips in LATISS and LSSTCam.",
1236 default=2,
1237 )
1238 leadingToSkip = pexConfig.Field(
1239 dtype=int,
1240 doc="Number of leading values to skip in parallel overscan correction.",
1241 default=0,
1242 )
1243 trailingToSkip = pexConfig.Field(
1244 dtype=int,
1245 doc="Number of trailing values to skip in parallel overscan correction.",
1246 default=0,
1247 )
1248 parallelOverscanFraction = pexConfig.Field(
1249 dtype=float,
1250 doc="When the parallel overscan region median is greater than parallelOverscanFraction "
1251 "and the imaging region median is greater than parallelOverscanImageThreshold "
1252 "then parallel overscan subtraction will be turned off, as this is usually "
1253 "due to the region being flooded with spillover from a super-saturated flat.",
1254 default=0.5,
1255 )
1256 parallelOverscanImageThreshold = pexConfig.Field(
1257 dtype=float,
1258 doc="When the parallel overscan region median is greater than parallelOverscanFraction "
1259 "and the imaging region median is greater than parallelOverscanImageThreshold "
1260 "then parallel overscan subtraction will be turned off, as this is usually "
1261 "due to the region being flooded with spillover from a super-saturated flat.",
1262 default=10000.0,
1263 )
1264
1265
1267 """Correction task for parallel overscan.
1268
1269 Parameters
1270 ----------
1271 statControl : `lsst.afw.math.StatisticsControl`, optional
1272 Statistics control object.
1273 """
1274 ConfigClass = ParallelOverscanCorrectionTaskConfig
1275 _DefaultName = "parallelOverscan"
1276
1277 def run(self, exposure, amp, isTransposed=False):
1278 """Measure and remove parallel overscan from an amplifier image.
1279
1280 This method assumes that serial overscan has already been
1281 removed from the amplifier.
1282
1283 Parameters
1284 ----------
1285 exposure : `lsst.afw.image.Exposure`
1286 Image data that will have the overscan corrections applied.
1287 amp : `lsst.afw.cameraGeom.Amplifier`
1288 Amplifier to use for debugging purposes.
1289 isTransposed : `bool`, optional
1290 Is the image transposed, such that serial and parallel
1291 overscan regions are reversed? Default is False.
1292
1293 Returns
1294 -------
1295 overscanResults : `lsst.pipe.base.Struct`
1296 Result struct with components:
1297
1298 ``imageFit``
1299 Value or fit subtracted from the amplifier image data
1300 (scalar or `lsst.afw.image.Image`).
1301 ``overscanFit``
1302 Value or fit subtracted from the parallel overscan image
1303 data (scalar or `lsst.afw.image.Image`).
1304 ``overscanImage``
1305 Image of the parallel overscan region with the parallel
1306 overscan correction applied
1307 (`lsst.afw.image.Image`). This quantity is used to
1308 estimate the amplifier read noise empirically.
1309 ``overscanMean``
1310 Mean of the fit parallel overscan region.
1311 ``overscanMedian``
1312 Median of the fit parallel overscan region.
1313 ``overscanSigma``
1314 Sigma of the fit parallel overscan region.
1315 ``residualMean``
1316 Mean of the residual of the parallel overscan region after
1317 correction.
1318 ``residualMedian``
1319 Median of the residual of the parallel overscan region after
1320 correction.
1321 ``residualSigma``
1322 Mean of the residual of the parallel overscan region after
1323 correction.
1324
1325 Raises
1326 ------
1327 RuntimeError
1328 Raised if an invalid overscan type is set.
1329 """
1330 # This does not need any extending, as we only subtract
1331 # from the data region.
1332 parallelOverscanBBox = amp.getRawParallelOverscanBBox()
1333 imageBBox = amp.getRawDataBBox()
1334
1335 # The serial overscan correction has removed some signal
1336 # from the parallel overscan region, but that is largely a
1337 # constant offset. The collapseArray method now attempts
1338 # to fill fully masked columns with the median of
1339 # neighboring values, with a fallback to the median of the
1340 # correction in all other columns. Filling with neighbor
1341 # values ensures that large variations in the parallel
1342 # overscan do not create new outlier points. The
1343 # MEDIAN_PER_ROW method does this filling as a separate
1344 # operation, using the same method.
1345 results = self.correctOverscan(
1346 exposure,
1347 amp,
1348 imageBBox,
1349 parallelOverscanBBox,
1350 isTransposed=not isTransposed,
1351 leadingToSkip=self.config.leadingToSkip,
1352 trailingToSkip=self.config.trailingToSkip,
1353 overscanFraction=self.config.parallelOverscanFraction,
1354 imageThreshold=self.config.parallelOverscanImageThreshold,
1355 maskedRowColumnGrowSize=self.config.parallelOverscanMaskedColumnGrowSize,
1356 )
1357 overscanMean = results.overscanMean
1358 overscanMedian = results.overscanMedian
1359 overscanSigma = results.overscanSigma
1360 residualMean = results.overscanMeanResidual
1361 residualMedian = results.overscanMedianResidual
1362 residualSigma = results.overscanSigmaResidual
1363
1364 return pipeBase.Struct(
1365 imageFit=results.ampOverscanModel,
1366 overscanFit=results.overscanOverscanModel,
1367 overscanImage=results.overscanImage,
1368 overscanMean=overscanMean,
1369 overscanMedian=overscanMedian,
1370 overscanSigma=overscanSigma,
1371 residualMean=residualMean,
1372 residualMedian=residualMedian,
1373 residualSigma=residualSigma,
1374 )
1375
1376 def maskParallelOverscanAmp(self, exposure, amp, saturationLevel=None):
1377 """Mask parallel overscan, growing saturated pixels.
1378
1379 This operates on the image in-place.
1380
1381 Parameters
1382 ----------
1383 exposure : `lsst.afw.image.Exposure`
1384 An untrimmed raw exposure.
1385 amp : `lsst.afw.cameraGeom.Amplifier`
1386 The amplifier to use for masking.
1387 saturationLevel : `float`, optional
1388 Saturation level to use for masking.
1389 """
1390 if not self.config.doParallelOverscanSaturation:
1391 # This is a no-op.
1392 return
1393
1394 if saturationLevel is None:
1395 saturationLevel = self.config.parallelOverscanSaturationLevel
1396
1397 dataView = afwImage.MaskedImageF(exposure.getMaskedImage(),
1398 amp.getRawParallelOverscanBBox(),
1399 afwImage.PARENT)
1400 # This should mark all of these saturated pixels as SAT.
1401 makeThresholdMask(
1402 maskedImage=dataView,
1403 threshold=saturationLevel,
1404 growFootprints=self.config.parallelOverscanMaskGrowSize,
1405 maskName="SAT",
1406 )
Pass parameters to a Statistics object.
Definition Statistics.h:83
An integer coordinate rectangle.
Definition Box.h:55
debugView(self, image, model, amp=None, isTransposed=True)
Definition overscan.py:834
fitOverscan(self, overscanImage, isTransposed=False)
Definition overscan.py:449
correctOverscan(self, exposure, amp, imageBBox, overscanBBox, isTransposed=True, leadingToSkip=0, trailingToSkip=0, overscanFraction=1.0, imageThreshold=np.inf, maskedRowColumnGrowSize=0)
Definition overscan.py:160
collapseArray(self, maskedArray, fillMasked=True)
Definition overscan.py:630
trimOverscan(self, exposure, amp, bbox, skipLeading, skipTrailing, transpose=False)
Definition overscan.py:402
maskParallelOverscan(self, exposure, detector)
Definition overscan.py:479
broadcastFitToImage(self, overscanValue, imageArray, transpose=False)
Definition overscan.py:352
splineFit(self, indices, collapsed, numBins)
Definition overscan.py:653
run(self, exposure, amp, isTransposed=False)
Definition overscan.py:154
__init__(self, statControl=None, **kwargs)
Definition overscan.py:143
measureVectorOverscan(self, image, isTransposed=False)
Definition overscan.py:747
run(self, exposure, amp, isTransposed=False)
Definition overscan.py:927
maskParallelOverscanAmp(self, exposure, amp, saturationLevel=None)
Definition overscan.py:1376
run(self, exposure, amp, isTransposed=False)
Definition overscan.py:1277
run(self, exposure, amp, isTransposed=False)
Definition overscan.py:1103
Statistics makeStatistics(lsst::afw::image::Image< Pixel > const &img, lsst::afw::image::Mask< image::MaskPixel > const &msk, int const flags, StatisticsControl const &sctrl=StatisticsControl())
Handle a watered-down front-end to the constructor (no variance)
Definition Statistics.h:361
Property stringToStatisticsProperty(std::string const property)
Conversion function to switch a string to a Property (see Statistics.h)
Interpolate::Style stringToInterpStyle(std::string const &style)
Conversion function to switch a string to an Interpolate::Style.
std::shared_ptr< Interpolate > makeInterpolate(std::vector< double > const &x, std::vector< double > const &y, Interpolate::Style const style=Interpolate::AKIMA_SPLINE)
A factory function to make Interpolate objects.