lsst.ip.isr g564bce9563+4d88945b8c
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
24import numpy as np
25import lsst.afw.math as afwMath
26import lsst.afw.image as afwImage
27import lsst.geom as geom
28import lsst.pipe.base as pipeBase
29import lsst.pex.config as pexConfig
30
31from .isr import fitOverscanImage
32
33
34class OverscanCorrectionTaskConfig(pexConfig.Config):
35 """Overscan correction options.
36 """
37 fitType = pexConfig.ChoiceField(
38 dtype=str,
39 doc="The method for fitting the overscan bias level.",
40 default='MEDIAN',
41 allowed={
42 "POLY": "Fit ordinary polynomial to the longest axis of the overscan region",
43 "CHEB": "Fit Chebyshev polynomial to the longest axis of the overscan region",
44 "LEG": "Fit Legendre polynomial to the longest axis of the overscan region",
45 "NATURAL_SPLINE": "Fit natural spline to the longest axis of the overscan region",
46 "CUBIC_SPLINE": "Fit cubic spline to the longest axis of the overscan region",
47 "AKIMA_SPLINE": "Fit Akima spline to the longest axis of the overscan region",
48 "MEAN": "Correct using the mean of the overscan region",
49 "MEANCLIP": "Correct using a clipped mean of the overscan region",
50 "MEDIAN": "Correct using the median of the overscan region",
51 "MEDIAN_PER_ROW": "Correct using the median per row of the overscan region",
52 },
53 )
54 order = pexConfig.Field(
55 dtype=int,
56 doc=("Order of polynomial to fit if overscan fit type is a polynomial, "
57 "or number of spline knots if overscan fit type is a spline."),
58 default=1,
59 )
60 numSigmaClip = pexConfig.Field(
61 dtype=float,
62 doc="Rejection threshold (sigma) for collapsing overscan before fit",
63 default=3.0,
64 )
65 maskPlanes = pexConfig.ListField(
66 dtype=str,
67 doc="Mask planes to reject when measuring overscan",
68 default=['BAD', 'SAT'],
69 )
70 overscanIsInt = pexConfig.Field(
71 dtype=bool,
72 doc="Treat overscan as an integer image for purposes of fitType=MEDIAN"
73 " and fitType=MEDIAN_PER_ROW.",
74 default=True,
75 )
76
77 doParallelOverscan = pexConfig.Field(
78 dtype=bool,
79 doc="Correct using parallel overscan after serial overscan correction?",
80 default=False,
81 )
82
83 leadingColumnsToSkip = pexConfig.Field(
84 dtype=int,
85 doc="Number of leading columns to skip in serial overscan correction.",
86 default=0,
87 )
88 trailingColumnsToSkip = pexConfig.Field(
89 dtype=int,
90 doc="Number of trailing columns to skip in serial overscan correction.",
91 default=0,
92 )
93 leadingRowsToSkip = pexConfig.Field(
94 dtype=int,
95 doc="Number of leading rows to skip in parallel overscan correction.",
96 default=0,
97 )
98 trailingRowsToSkip = pexConfig.Field(
99 dtype=int,
100 doc="Number of trailing rows to skip in parallel overscan correction.",
101 default=0,
102 )
103
104 maxDeviation = pexConfig.Field(
105 dtype=float,
106 doc="Maximum deviation from median (in ADU) to mask in overscan correction.",
107 default=1000.0, check=lambda x: x > 0,
108 )
109
110
111class OverscanCorrectionTask(pipeBase.Task):
112 """Correction task for overscan.
113
114 This class contains a number of utilities that are easier to
115 understand and use when they are not embedded in nested if/else
116 loops.
117
118 Parameters
119 ----------
120 statControl : `lsst.afw.math.StatisticsControl`, optional
121 Statistics control object.
122 """
123 ConfigClass = OverscanCorrectionTaskConfig
124 _DefaultName = "overscan"
125
126 def __init__(self, statControl=None, **kwargs):
127 super().__init__(**kwargs)
128 self.allowDebug = True
129
130 if statControl:
131 self.statControl = statControl
132 else:
133 self.statControl = afwMath.StatisticsControl()
134 self.statControl.setNumSigmaClip(self.config.numSigmaClip)
135 self.statControl.setAndMask(afwImage.Mask.getPlaneBitMask(self.config.maskPlanes))
136
137 def run(self, exposure, amp, isTransposed=False):
138 """Measure and remove an overscan from an amplifier image.
139
140 Parameters
141 ----------
142 exposure : `lsst.afw.image.Exposure`
143 Image data that will have the overscan corrections applied.
145 Amplifier to use for debugging purposes.
146 isTransposed : `bool`, optional
147 Is the image transposed, such that serial and parallel
148 overscan regions are reversed? Default is False.
149
150 Returns
151 -------
152 overscanResults : `lsst.pipe.base.Struct`
153 Result struct with components:
154
155 ``imageFit``
156 Value or fit subtracted from the amplifier image data
157 (scalar or `lsst.afw.image.Image`).
158 ``overscanFit``
159 Value or fit subtracted from the overscan image data
160 (scalar or `lsst.afw.image.Image`).
161 ``overscanImage``
162 Image of the overscan region with the overscan
163 correction applied (`lsst.afw.image.Image`). This
164 quantity is used to estimate the amplifier read noise
165 empirically.
166
167 Raises
168 ------
169 RuntimeError
170 Raised if an invalid overscan type is set.
171 """
172 # Do Serial overscan first.
173 serialOverscanBBox = amp.getRawSerialOverscanBBox()
174 imageBBox = amp.getRawDataBBox()
175
176 if self.config.doParallelOverscan:
177 # We need to extend the serial overscan BBox to the full
178 # size of the detector.
179 parallelOverscanBBox = amp.getRawParallelOverscanBBox()
180 imageBBox = imageBBox.expandedTo(parallelOverscanBBox)
181
182 serialOverscanBBox = geom.Box2I(geom.Point2I(serialOverscanBBox.getMinX(),
183 imageBBox.getMinY()),
184 geom.Extent2I(serialOverscanBBox.getWidth(),
185 imageBBox.getHeight()))
186 serialResults = self.correctOverscan(exposure, amp,
187 imageBBox, serialOverscanBBox, isTransposed=isTransposed)
188 overscanMean = serialResults.overscanMean
189 overscanSigma = serialResults.overscanSigma
190 residualMean = serialResults.overscanMeanResidual
191 residualSigma = serialResults.overscanSigmaResidual
192
193 # Do Parallel Overscan
194 if self.config.doParallelOverscan:
195 # This does not need any extensions, as we'll only
196 # subtract it from the data region.
197 parallelOverscanBBox = amp.getRawParallelOverscanBBox()
198 imageBBox = amp.getRawDataBBox()
199
200 parallelResults = self.correctOverscan(exposure, amp,
201 imageBBox, parallelOverscanBBox,
202 isTransposed=not isTransposed)
203
204 overscanMean = (overscanMean, parallelResults.overscanMean)
205 overscanSigma = (overscanSigma, parallelResults.overscanSigma)
206 residualMean = (residualMean, parallelResults.overscanMeanResidual)
207 residualSigma = (residualSigma, parallelResults.overscanSigmaResidual)
208
209 return pipeBase.Struct(imageFit=serialResults.ampOverscanModel,
210 overscanFit=serialResults.overscanOverscanModel,
211 overscanImage=serialResults.overscanImage,
212
213 overscanMean=overscanMean,
214 overscanSigma=overscanSigma,
215 residualMean=residualMean,
216 residualSigma=residualSigma)
217
218 def correctOverscan(self, exposure, amp, imageBBox, overscanBBox, isTransposed=True):
219 """
220 """
221 overscanBox = self.trimOverscan(exposure, amp, overscanBBox,
222 self.config.leadingColumnsToSkip,
223 self.config.trailingColumnsToSkip,
224 transpose=isTransposed)
225 overscanImage = exposure[overscanBox].getMaskedImage()
226 overscanArray = overscanImage.image.array
227
228 # Mask pixels.
229 maskVal = overscanImage.mask.getPlaneBitMask(self.config.maskPlanes)
230 overscanMask = ~((overscanImage.mask.array & maskVal) == 0)
231
232 median = np.ma.median(np.ma.masked_where(overscanMask, overscanArray))
233 bad = np.where(np.abs(overscanArray - median) > self.config.maxDeviation)
234 overscanMask[bad] = overscanImage.mask.getPlaneBitMask("SAT")
235
236 # Do overscan fit.
237 # CZW: Handle transposed correctly.
238 overscanResults = self.fitOverscan(overscanImage, isTransposed=isTransposed)
239
240 # Correct image region (and possibly parallel-overscan region).
241 ampImage = exposure[imageBBox]
242 ampOverscanModel = self.broadcastFitToImage(overscanResults.overscanValue,
243 ampImage.image.array,
244 transpose=isTransposed)
245 ampImage.image.array -= ampOverscanModel
246
247 # Correct overscan region (and possibly doubly-overscaned
248 # region).
249 overscanImage = exposure[overscanBBox]
250 # CZW: Transposed?
251 overscanOverscanModel = self.broadcastFitToImage(overscanResults.overscanValue,
252 overscanImage.image.array)
253 overscanImage.image.array -= overscanOverscanModel
254
255 self.debugView(overscanImage, overscanResults.overscanValue, amp)
256
257 # Find residual fit statistics.
258 stats = afwMath.makeStatistics(overscanImage.getMaskedImage(),
259 afwMath.MEDIAN | afwMath.STDEVCLIP, self.statControl)
260 residualMean = stats.getValue(afwMath.MEDIAN)
261 residualSigma = stats.getValue(afwMath.STDEVCLIP)
262
263 return pipeBase.Struct(ampOverscanModel=ampOverscanModel,
264 overscanOverscanModel=overscanOverscanModel,
265 overscanImage=overscanImage,
266 overscanValue=overscanResults.overscanValue,
267
268 overscanMean=overscanResults.overscanMean,
269 overscanSigma=overscanResults.overscanSigma,
270 overscanMeanResidual=residualMean,
271 overscanSigmaResidual=residualSigma
272 )
273
274 def broadcastFitToImage(self, overscanValue, imageArray, transpose=False):
275 """Broadcast 0 or 1 dimension fit to appropriate shape.
276
277 Parameters
278 ----------
279 overscanValue : `numpy.ndarray`, (Nrows, ) or scalar
280 Overscan fit to broadcast.
281 imageArray : `numpy.ndarray`, (Nrows, Ncols)
282 Image array that we want to match.
283 transpose : `bool`, optional
284 Switch order to broadcast along the other axis.
285
286 Returns
287 -------
288 overscanModel : `numpy.ndarray`, (Nrows, Ncols) or scalar
289 Expanded overscan fit.
290
291 Raises
292 ------
293 RuntimeError
294 Raised if no axis has the appropriate dimension.
295 """
296 if isinstance(overscanValue, np.ndarray):
297 overscanModel = np.zeros_like(imageArray)
298
299 if transpose is False:
300 if imageArray.shape[0] == overscanValue.shape[0]:
301 overscanModel[:, :] = overscanValue[:, np.newaxis]
302 elif imageArray.shape[1] == overscanValue.shape[0]:
303 overscanModel[:, :] = overscanValue[np.newaxis, :]
304 elif imageArray.shape[0] == overscanValue.shape[1]:
305 overscanModel[:, :] = overscanValue[np.newaxis, :]
306 else:
307 raise RuntimeError(f"Could not broadcast {overscanValue.shape} to "
308 f"match {imageArray.shape}")
309 else:
310 if imageArray.shape[1] == overscanValue.shape[0]:
311 overscanModel[:, :] = overscanValue[np.newaxis, :]
312 elif imageArray.shape[0] == overscanValue.shape[0]:
313 overscanModel[:, :] = overscanValue[:, np.newaxis]
314 elif imageArray.shape[1] == overscanValue.shape[1]:
315 overscanModel[:, :] = overscanValue[:, np.newaxis]
316 else:
317 raise RuntimeError(f"Could not broadcast {overscanValue.shape} to "
318 f"match {imageArray.shape}")
319 else:
320 overscanModel = overscanValue
321
322 return overscanModel
323
324 def trimOverscan(self, exposure, amp, bbox, skipLeading, skipTrailing, transpose=False):
325 """Trim overscan region to remove edges.
326
327 Parameters
328 ----------
329 exposure : `lsst.afw.image.Exposure`
330 Exposure containing data.
332 Amplifier containing geometry information.
333 bbox : `lsst.geom.Box2I`
334 Bounding box of the overscan region.
335 skipLeading : `int`
336 Number of leading (towards data region) rows/columns to skip.
337 skipTrailing : `int`
338 Number of trailing (away from data region) rows/columns to skip.
339 transpose : `bool`, optional
340 Operate on the transposed array.
341
342 Returns
343 -------
344 overscanArray : `numpy.array`, (N, M)
345 Data array to fit.
346 overscanMask : `numpy.array`, (N, M)
347 Data mask.
348 """
349 dx0, dy0, dx1, dy1 = (0, 0, 0, 0)
350 dataBBox = amp.getRawDataBBox()
351 if transpose:
352 if dataBBox.getBeginY() < bbox.getBeginY():
353 dy0 += skipLeading
354 dy1 -= skipTrailing
355 else:
356 dy0 += skipTrailing
357 dy1 -= skipLeading
358 else:
359 if dataBBox.getBeginX() < bbox.getBeginX():
360 dx0 += skipLeading
361 dx1 -= skipTrailing
362 else:
363 dx0 += skipTrailing
364 dx1 -= skipLeading
365
366 overscanBBox = geom.Box2I(bbox.getBegin() + geom.Extent2I(dx0, dy0),
367 geom.Extent2I(bbox.getWidth() - dx0 + dx1,
368 bbox.getHeight() - dy0 + dy1))
369 return overscanBBox
370
371 def fitOverscan(self, overscanImage, isTransposed=False):
372 if self.config.fitType in ('MEAN', 'MEANCLIP', 'MEDIAN'):
373 # Transposition has no effect here.
374 overscanResult = self.measureConstantOverscan(overscanImage)
375 overscanValue = overscanResult.overscanValue
376 overscanMean = overscanValue
377 overscanSigma = 0.0
378 elif self.config.fitType in ('MEDIAN_PER_ROW', 'POLY', 'CHEB', 'LEG',
379 'NATURAL_SPLINE', 'CUBIC_SPLINE', 'AKIMA_SPLINE'):
380 # Force transposes as needed
381 overscanResult = self.measureVectorOverscan(overscanImage, isTransposed)
382 overscanValue = overscanResult.overscanValue
383
384 stats = afwMath.makeStatistics(overscanResult.overscanValue,
385 afwMath.MEDIAN | afwMath.STDEVCLIP, self.statControl)
386 overscanMean = stats.getValue(afwMath.MEDIAN)
387 overscanSigma = stats.getValue(afwMath.STDEVCLIP)
388 else:
389 raise ValueError('%s : %s an invalid overscan type' %
390 ("overscanCorrection", self.config.fitType))
391
392 return pipeBase.Struct(overscanValue=overscanValue,
393 overscanMean=overscanMean,
394 overscanSigma=overscanSigma,
395 )
396
397 @staticmethod
398 def integerConvert(image):
399 """Return an integer version of the input image.
400
401 Parameters
402 ----------
403 image : `numpy.ndarray`, `lsst.afw.image.Image` or `MaskedImage`
404 Image to convert to integers.
405
406 Returns
407 -------
408 outI : `numpy.ndarray`, `lsst.afw.image.Image` or `MaskedImage`
409 The integer converted image.
410
411 Raises
412 ------
413 RuntimeError
414 Raised if the input image could not be converted.
415 """
416 if hasattr(image, "image"):
417 # Is a maskedImage:
418 imageI = image.image.convertI()
419 outI = afwImage.MaskedImageI(imageI, image.mask, image.variance)
420 elif hasattr(image, "convertI"):
421 # Is an Image:
422 outI = image.convertI()
423 elif hasattr(image, "astype"):
424 # Is a numpy array:
425 outI = image.astype(int)
426 else:
427 raise RuntimeError("Could not convert this to integers: %s %s %s",
428 image, type(image), dir(image))
429 return outI
430
431 # Constant methods
432 def measureConstantOverscan(self, image):
433 """Measure a constant overscan value.
434
435 Parameters
436 ----------
438 Image data to measure the overscan from.
439
440 Returns
441 -------
442 results : `lsst.pipe.base.Struct`
443 Overscan result with entries:
444 - ``overscanValue``: Overscan value to subtract (`float`)
445 - ``isTransposed``: Orientation of the overscan (`bool`)
446 """
447 if self.config.fitType == 'MEDIAN':
448 calcImage = self.integerConvert(image)
449 else:
450 calcImage = image
451 fitType = afwMath.stringToStatisticsProperty(self.config.fitType)
452 overscanValue = afwMath.makeStatistics(calcImage, fitType, self.statControl).getValue()
453
454 return pipeBase.Struct(overscanValue=overscanValue,
455 isTransposed=False)
456
457 # Vector correction utilities
458 def getImageArray(self, image):
459 """Extract the numpy array from the input image.
460
461 Parameters
462 ----------
464 Image data to pull array from.
465
466 calcImage : `numpy.ndarray`
467 Image data array for numpy operating.
468 """
469 if hasattr(image, "getImage"):
470 calcImage = image.getImage().getArray()
471 calcImage = np.ma.masked_where(image.getMask().getArray() & self.statControl.getAndMask(),
472 calcImage)
473 else:
474 calcImage = image.getArray()
475 return calcImage
476
477 def maskOutliers(self, imageArray):
478 """Mask outliers in a row of overscan data from a robust sigma
479 clipping procedure.
480
481 Parameters
482 ----------
483 imageArray : `numpy.ndarray`
484 Image to filter along numpy axis=1.
485
486 Returns
487 -------
488 maskedArray : `numpy.ma.masked_array`
489 Masked image marking outliers.
490 """
491 lq, median, uq = np.percentile(imageArray, [25.0, 50.0, 75.0], axis=1)
492 axisMedians = median
493 axisStdev = 0.74*(uq - lq) # robust stdev
494
495 diff = np.abs(imageArray - axisMedians[:, np.newaxis])
496 return np.ma.masked_where(diff > self.statControl.getNumSigmaClip()
497 * axisStdev[:, np.newaxis], imageArray)
498
499 @staticmethod
500 def collapseArray(maskedArray):
501 """Collapse overscan array (and mask) to a 1-D vector of values.
502
503 Parameters
504 ----------
505 maskedArray : `numpy.ma.masked_array`
506 Masked array of input overscan data.
507
508 Returns
509 -------
510 collapsed : `numpy.ma.masked_array`
511 Single dimensional overscan data, combined with the mean.
512 """
513 collapsed = np.mean(maskedArray, axis=1)
514 if collapsed.mask.sum() > 0:
515 collapsed.data[collapsed.mask] = np.mean(maskedArray.data[collapsed.mask], axis=1)
516 return collapsed
517
518 def collapseArrayMedian(self, maskedArray):
519 """Collapse overscan array (and mask) to a 1-D vector of using the
520 correct integer median of row-values.
521
522 Parameters
523 ----------
524 maskedArray : `numpy.ma.masked_array`
525 Masked array of input overscan data.
526
527 Returns
528 -------
529 collapsed : `numpy.ma.masked_array`
530 Single dimensional overscan data, combined with the afwMath median.
531 """
532 integerMI = self.integerConvert(maskedArray)
533
534 collapsed = []
535 fitType = afwMath.stringToStatisticsProperty('MEDIAN')
536 for row in integerMI:
537 newRow = row.compressed()
538 if len(newRow) > 0:
539 rowMedian = afwMath.makeStatistics(newRow, fitType, self.statControl).getValue()
540 else:
541 rowMedian = np.nan
542 collapsed.append(rowMedian)
543
544 return np.array(collapsed)
545
546 def splineFit(self, indices, collapsed, numBins):
547 """Wrapper function to match spline fit API to polynomial fit API.
548
549 Parameters
550 ----------
551 indices : `numpy.ndarray`
552 Locations to evaluate the spline.
553 collapsed : `numpy.ndarray`
554 Collapsed overscan values corresponding to the spline
555 evaluation points.
556 numBins : `int`
557 Number of bins to use in constructing the spline.
558
559 Returns
560 -------
562 Interpolation object for later evaluation.
563 """
564 if not np.ma.is_masked(collapsed):
565 collapsed.mask = np.array(len(collapsed)*[np.ma.nomask])
566
567 numPerBin, binEdges = np.histogram(indices, bins=numBins,
568 weights=1 - collapsed.mask.astype(int))
569 with np.errstate(invalid="ignore"):
570 values = np.histogram(indices, bins=numBins,
571 weights=collapsed.data*~collapsed.mask)[0]/numPerBin
572 binCenters = np.histogram(indices, bins=numBins,
573 weights=indices*~collapsed.mask)[0]/numPerBin
574
575 if len(binCenters[numPerBin > 0]) < 5:
576 self.log.warn("Cannot do spline fitting for overscan: %s valid points.",
577 len(binCenters[numPerBin > 0]))
578 # Return a scalar value if we have one, otherwise
579 # return zero. This amplifier is hopefully already
580 # masked.
581 if len(values[numPerBin > 0]) != 0:
582 return float(values[numPerBin > 0][0])
583 else:
584 return 0.0
585
586 interp = afwMath.makeInterpolate(binCenters.astype(float)[numPerBin > 0],
587 values.astype(float)[numPerBin > 0],
588 afwMath.stringToInterpStyle(self.config.fitType))
589 return interp
590
591 @staticmethod
592 def splineEval(indices, interp):
593 """Wrapper function to match spline evaluation API to polynomial fit
594 API.
595
596 Parameters
597 ----------
598 indices : `numpy.ndarray`
599 Locations to evaluate the spline.
600 interp : `lsst.afw.math.interpolate`
601 Interpolation object to use.
602
603 Returns
604 -------
605 values : `numpy.ndarray`
606 Evaluated spline values at each index.
607 """
608
609 return interp.interpolate(indices.astype(float))
610
611 @staticmethod
612 def maskExtrapolated(collapsed):
613 """Create mask if edges are extrapolated.
614
615 Parameters
616 ----------
617 collapsed : `numpy.ma.masked_array`
618 Masked array to check the edges of.
619
620 Returns
621 -------
622 maskArray : `numpy.ndarray`
623 Boolean numpy array of pixels to mask.
624 """
625 maskArray = np.full_like(collapsed, False, dtype=bool)
626 if np.ma.is_masked(collapsed):
627 num = len(collapsed)
628 for low in range(num):
629 if not collapsed.mask[low]:
630 break
631 if low > 0:
632 maskArray[:low] = True
633 for high in range(1, num):
634 if not collapsed.mask[-high]:
635 break
636 if high > 1:
637 maskArray[-high:] = True
638 return maskArray
639
640 def measureVectorOverscan(self, image, isTransposed=False):
641 """Calculate the 1-d vector overscan from the input overscan image.
642
643 Parameters
644 ----------
646 Image containing the overscan data.
647 isTransposed : `bool`
648 If true, the image has been transposed.
649
650 Returns
651 -------
652 results : `lsst.pipe.base.Struct`
653 Overscan result with entries:
654
655 ``overscanValue``
656 Overscan value to subtract (`float`)
657 ``maskArray``
658 List of rows that should be masked as ``SUSPECT`` when the
659 overscan solution is applied. (`list` [ `bool` ])
660 ``isTransposed``
661 Indicates if the overscan data was transposed during
662 calcuation, noting along which axis the overscan should be
663 subtracted. (`bool`)
664 """
665 calcImage = self.getImageArray(image)
666
667 # operate on numpy-arrays from here
668 if isTransposed:
669 calcImage = np.transpose(calcImage)
670 masked = self.maskOutliers(calcImage)
671
672 if self.config.fitType == 'MEDIAN_PER_ROW':
673 mi = afwImage.MaskedImageI(image.getBBox())
674 masked = masked.astype(int)
675 if isTransposed:
676 masked = masked.transpose()
677
678 mi.image.array[:, :] = masked.data[:, :]
679 if bool(masked.mask.shape):
680 mi.mask.array[:, :] = masked.mask[:, :]
681
682 overscanVector = fitOverscanImage(mi, self.config.maskPlanes, isTransposed)
683 maskArray = self.maskExtrapolated(overscanVector)
684 else:
685 collapsed = self.collapseArray(masked)
686
687 num = len(collapsed)
688 indices = 2.0*np.arange(num)/float(num) - 1.0
689
690 poly = np.polynomial
691 fitter, evaler = {
692 'POLY': (poly.polynomial.polyfit, poly.polynomial.polyval),
693 'CHEB': (poly.chebyshev.chebfit, poly.chebyshev.chebval),
694 'LEG': (poly.legendre.legfit, poly.legendre.legval),
695 'NATURAL_SPLINE': (self.splineFit, self.splineEval),
696 'CUBIC_SPLINE': (self.splineFit, self.splineEval),
697 'AKIMA_SPLINE': (self.splineFit, self.splineEval)
698 }[self.config.fitType]
699
700 # These are the polynomial coefficients, or an
701 # interpolation object.
702 coeffs = fitter(indices, collapsed, self.config.order)
703
704 if isinstance(coeffs, float):
705 self.log.warn("Using fallback value %f due to fitter failure. Amplifier will be masked.",
706 coeffs)
707 overscanVector = np.full_like(indices, coeffs)
708 maskArray = np.full_like(collapsed, True, dtype=bool)
709 else:
710 # Otherwise we can just use things as normal.
711 overscanVector = evaler(indices, coeffs)
712 maskArray = self.maskExtrapolated(collapsed)
713
714 return pipeBase.Struct(overscanValue=np.array(overscanVector),
715 maskArray=maskArray,
716 isTransposed=isTransposed)
717
718 def debugView(self, image, model, amp=None):
719 """Debug display for the final overscan solution.
720
721 Parameters
722 ----------
723 image : `lsst.afw.image.Image`
724 Input image the overscan solution was determined from.
725 model : `numpy.ndarray` or `float`
726 Overscan model determined for the image.
727 amp : `lsst.afw.cameraGeom.Amplifier`, optional
728 Amplifier to extract diagnostic information.
729 """
730 import lsstDebug
731 if not lsstDebug.Info(__name__).display:
732 return
733 if not self.allowDebug:
734 return
735
736 calcImage = self.getImageArray(image)
737 # CZW: Check that this is ok
738 calcImage = np.transpose(calcImage)
739 masked = self.maskOutliers(calcImage)
740 collapsed = self.collapseArray(masked)
741
742 num = len(collapsed)
743 indices = 2.0 * np.arange(num)/float(num) - 1.0
744
745 if np.ma.is_masked(collapsed):
746 collapsedMask = collapsed.mask
747 else:
748 collapsedMask = np.array(num*[np.ma.nomask])
749
750 import matplotlib.pyplot as plot
751 figure = plot.figure(1)
752 figure.clear()
753 axes = figure.add_axes((0.1, 0.1, 0.8, 0.8))
754 axes.plot(indices[~collapsedMask], collapsed[~collapsedMask], 'k+')
755 if collapsedMask.sum() > 0:
756 axes.plot(indices[collapsedMask], collapsed.data[collapsedMask], 'b+')
757 if isinstance(model, np.ndarray):
758 plotModel = model
759 else:
760 plotModel = np.zeros_like(indices)
761 plotModel += model
762 axes.plot(indices, plotModel, 'r-')
763 plot.xlabel("centered/scaled position along overscan region")
764 plot.ylabel("pixel value/fit value")
765 if amp:
766 plot.title(f"{amp.getName()} DataX: "
767 f"[{amp.getRawDataBBox().getBeginX()}:{amp.getRawBBox().getEndX()}]"
768 f"OscanX: [{amp.getRawHorizontalOverscanBBox().getBeginX()}:"
769 f"{amp.getRawHorizontalOverscanBBox().getEndX()}] {self.config.fitType}")
770 else:
771 plot.title("No amp supplied.")
772 figure.show()
773 prompt = "Press Enter or c to continue [chp]..."
774 while True:
775 ans = input(prompt).lower()
776 if ans in ("", " ", "c",):
777 break
778 elif ans in ("p", ):
779 import pdb
780 pdb.set_trace()
781 elif ans in ('x', ):
782 self.allowDebug = False
783 break
784 elif ans in ("h", ):
785 print("[h]elp [c]ontinue [p]db e[x]itDebug")
786 plot.close()
def collapseArrayMedian(self, maskedArray)
Definition: overscan.py:518
def fitOverscan(self, overscanImage, isTransposed=False)
Definition: overscan.py:371
def measureVectorOverscan(self, image, isTransposed=False)
Definition: overscan.py:640
def __init__(self, statControl=None, **kwargs)
Definition: overscan.py:126
def debugView(self, image, model, amp=None)
Definition: overscan.py:718
def splineFit(self, indices, collapsed, numBins)
Definition: overscan.py:546
def broadcastFitToImage(self, overscanValue, imageArray, transpose=False)
Definition: overscan.py:274
def correctOverscan(self, exposure, amp, imageBBox, overscanBBox, isTransposed=True)
Definition: overscan.py:218
def trimOverscan(self, exposure, amp, bbox, skipLeading, skipTrailing, transpose=False)
Definition: overscan.py:324
def run(self, exposure, amp, isTransposed=False)
Definition: overscan.py:137
std::vector< double > fitOverscanImage(lsst::afw::image::MaskedImage< ImagePixelT > const &overscan, std::vector< std::string > badPixelMask, bool isTransposed)
Definition: Isr.cc:53