29from .isr
import fitOverscanImage
31__all__ = [
"OverscanCorrectionTaskConfig",
"OverscanCorrectionTask"]
35 """Overscan correction options.
37 fitType = pexConfig.ChoiceField(
39 doc="The method for fitting the overscan bias level.",
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",
54 order = pexConfig.Field(
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."),
60 numSigmaClip = pexConfig.Field(
62 doc=
"Rejection threshold (sigma) for collapsing overscan before fit",
65 maskPlanes = pexConfig.ListField(
67 doc=
"Mask planes to reject when measuring overscan",
68 default=[
'BAD',
'SAT'],
70 overscanIsInt = pexConfig.Field(
72 doc=
"Treat overscan as an integer image for purposes of fitType=MEDIAN"
73 " and fitType=MEDIAN_PER_ROW.",
79 """Correction task for overscan.
81 This class contains a number of utilities that are easier
to
82 understand
and use when they are
not embedded
in nested
if/
else
88 Statistics control object.
90 ConfigClass = OverscanCorrectionTaskConfig
91 _DefaultName = "overscan"
93 def __init__(self, statControl=None, **kwargs):
101 self.
statControlstatControl.setNumSigmaClip(self.config.numSigmaClip)
102 self.
statControlstatControl.setAndMask(afwImage.Mask.getPlaneBitMask(self.config.maskPlanes))
104 def run(self, ampImage, overscanImage, amp=None):
105 """Measure and remove an overscan from an amplifier image.
110 Image data that will have the overscan removed.
112 Overscan data that the overscan is measured
from.
114 Amplifier to use
for debugging purposes.
118 overscanResults : `lsst.pipe.base.Struct`
119 Result struct
with components:
122 Value
or fit subtracted
from the amplifier image data
125 Value
or fit subtracted
from the overscan image data
128 Image of the overscan region
with the overscan
130 quantity
is used to estimate the amplifier read noise
136 Raised
if an invalid overscan type
is set.
139 if self.config.fitType
in (
'MEAN',
'MEANCLIP',
'MEDIAN'):
141 overscanValue = overscanResult.overscanValue
142 offImage = overscanValue
143 overscanModel = overscanValue
145 elif self.config.fitType
in (
'MEDIAN_PER_ROW',
'POLY',
'CHEB',
'LEG',
146 'NATURAL_SPLINE',
'CUBIC_SPLINE',
'AKIMA_SPLINE'):
148 overscanValue = overscanResult.overscanValue
149 maskArray = overscanResult.maskArray
150 isTransposed = overscanResult.isTransposed
152 offImage = afwImage.ImageF(ampImage.getDimensions())
153 offArray = offImage.getArray()
154 overscanModel = afwImage.ImageF(overscanImage.getDimensions())
155 overscanArray = overscanModel.getArray()
157 if hasattr(ampImage,
'getMask'):
163 offArray[:, :] = overscanValue[np.newaxis, :]
164 overscanArray[:, :] = overscanValue[np.newaxis, :]
166 maskSuspect.getArray()[:, maskArray] |= ampImage.getMask().getPlaneBitMask(
"SUSPECT")
168 offArray[:, :] = overscanValue[:, np.newaxis]
169 overscanArray[:, :] = overscanValue[:, np.newaxis]
171 maskSuspect.getArray()[maskArray, :] |= ampImage.getMask().getPlaneBitMask(
"SUSPECT")
173 raise RuntimeError(
'%s : %s an invalid overscan type' %
174 (
"overscanCorrection", self.config.fitType))
176 self.
debugViewdebugView(overscanImage, overscanValue, amp)
180 ampImage.getMask().getArray()[:, :] |= maskSuspect.getArray()[:, :]
181 overscanImage -= overscanModel
182 return pipeBase.Struct(imageFit=offImage,
183 overscanFit=overscanModel,
184 overscanImage=overscanImage,
185 edgeMask=maskSuspect)
189 """Return an integer version of the input image.
194 Image to convert to integers.
199 The integer converted image.
204 Raised
if the input image could
not be converted.
206 if hasattr(image,
"image"):
208 imageI = image.image.convertI()
209 outI = afwImage.MaskedImageI(imageI, image.mask, image.variance)
210 elif hasattr(image,
"convertI"):
212 outI = image.convertI()
213 elif hasattr(image,
"astype"):
215 outI = image.astype(int)
217 raise RuntimeError(
"Could not convert this to integers: %s %s %s",
218 image,
type(image), dir(image))
223 """Measure a constant overscan value.
228 Image data to measure the overscan
from.
232 results : `lsst.pipe.base.Struct`
233 Overscan result
with entries:
234 - ``overscanValue``: Overscan value to subtract (`float`)
235 - ``maskArray``: Placeholder
for a mask array (`list`)
236 - ``isTransposed``: Orientation of the overscan (`bool`)
238 if self.config.fitType ==
'MEDIAN':
246 return pipeBase.Struct(overscanValue=overscanValue,
252 """Extract the numpy array from the input image.
257 Image data to pull array
from.
259 calcImage : `numpy.ndarray`
260 Image data array
for numpy operating.
262 if hasattr(image,
"getImage"):
263 calcImage = image.getImage().getArray()
264 calcImage = np.ma.masked_where(image.getMask().getArray() & self.
statControlstatControl.getAndMask(),
267 calcImage = image.getArray()
272 """Transpose input numpy array if necessary.
276 imageArray : `numpy.ndarray`
277 Image data to transpose.
281 imageArray : `numpy.ndarray`
282 Transposed image data.
283 isTransposed : `bool`
284 Indicates whether the input data was transposed.
286 if np.argmin(imageArray.shape) == 0:
287 return np.transpose(imageArray),
True
289 return imageArray,
False
292 """Mask outliers in a row of overscan data
from a robust sigma
297 imageArray : `numpy.ndarray`
298 Image to filter along numpy axis=1.
302 maskedArray : `numpy.ma.masked_array`
303 Masked image marking outliers.
305 lq, median, uq = np.percentile(imageArray, [25.0, 50.0, 75.0], axis=1)
307 axisStdev = 0.74*(uq - lq)
309 diff = np.abs(imageArray - axisMedians[:, np.newaxis])
310 return np.ma.masked_where(diff > self.
statControlstatControl.getNumSigmaClip()
311 * axisStdev[:, np.newaxis], imageArray)
315 """Collapse overscan array (and mask) to a 1-D vector of values.
319 maskedArray : `numpy.ma.masked_array`
320 Masked array of input overscan data.
324 collapsed : `numpy.ma.masked_array`
325 Single dimensional overscan data, combined with the mean.
327 collapsed = np.mean(maskedArray, axis=1)
328 if collapsed.mask.sum() > 0:
329 collapsed.data[collapsed.mask] = np.mean(maskedArray.data[collapsed.mask], axis=1)
333 """Collapse overscan array (and mask) to a 1-D vector of using the
334 correct integer median of row-values.
338 maskedArray : `numpy.ma.masked_array`
339 Masked array of input overscan data.
343 collapsed : `numpy.ma.masked_array`
344 Single dimensional overscan data, combined with the afwMath median.
350 for row
in integerMI:
351 newRow = row.compressed()
356 collapsed.append(rowMedian)
358 return np.array(collapsed)
361 """Wrapper function to match spline fit API to polynomial fit API.
365 indices : `numpy.ndarray`
366 Locations to evaluate the spline.
367 collapsed : `numpy.ndarray`
368 Collapsed overscan values corresponding to the spline
371 Number of bins to use in constructing the spline.
376 Interpolation object
for later evaluation.
378 if not np.ma.is_masked(collapsed):
379 collapsed.mask = np.array(len(collapsed)*[np.ma.nomask])
381 numPerBin, binEdges = np.histogram(indices, bins=numBins,
382 weights=1 - collapsed.mask.astype(int))
383 with np.errstate(invalid=
"ignore"):
384 values = np.histogram(indices, bins=numBins,
385 weights=collapsed.data*~collapsed.mask)[0]/numPerBin
386 binCenters = np.histogram(indices, bins=numBins,
387 weights=indices*~collapsed.mask)[0]/numPerBin
389 if len(binCenters[numPerBin > 0]) < 5:
390 self.log.
warn(
"Cannot do spline fitting for overscan: %s valid points.",
391 len(binCenters[numPerBin > 0]))
395 if len(values[numPerBin > 0]) != 0:
396 return float(values[numPerBin > 0][0])
401 values.astype(float)[numPerBin > 0],
407 """Wrapper function to match spline evaluation API to polynomial fit
412 indices : `numpy.ndarray`
413 Locations to evaluate the spline.
414 interp : `lsst.afw.math.interpolate`
415 Interpolation object to use.
419 values : `numpy.ndarray`
420 Evaluated spline values at each index.
423 return interp.interpolate(indices.astype(float))
427 """Create mask if edges are extrapolated.
431 collapsed : `numpy.ma.masked_array`
432 Masked array to check the edges of.
436 maskArray : `numpy.ndarray`
437 Boolean numpy array of pixels to mask.
439 maskArray = np.full_like(collapsed, False, dtype=bool)
440 if np.ma.is_masked(collapsed):
442 for low
in range(num):
443 if not collapsed.mask[low]:
446 maskArray[:low] =
True
447 for high
in range(1, num):
448 if not collapsed.mask[-high]:
451 maskArray[-high:] =
True
455 """Calculate the 1-d vector overscan from the input overscan image.
460 Image containing the overscan data.
464 results : `lsst.pipe.base.Struct`
465 Overscan result with entries:
466 - ``overscanValue``: Overscan value to subtract (`float`)
467 - ``maskArray`` : `list` [ `bool` ]
468 List of rows that should be masked
as ``SUSPECT`` when the
469 overscan solution
is applied.
470 - ``isTransposed`` : `bool`
471 Indicates
if the overscan data was transposed during
472 calcuation, noting along which axis the overscan should be
478 calcImage, isTransposed = self.
transposetranspose(calcImage)
481 startTime = time.perf_counter()
483 if self.config.fitType ==
'MEDIAN_PER_ROW':
484 mi = afwImage.MaskedImageI(image.getBBox())
485 masked = masked.astype(int)
487 masked = masked.transpose()
489 mi.image.array[:, :] = masked.data[:, :]
490 if bool(masked.mask.shape):
491 mi.mask.array[:, :] = masked.mask[:, :]
499 indices = 2.0*np.arange(num)/
float(num) - 1.0
503 'POLY': (poly.polynomial.polyfit, poly.polynomial.polyval),
504 'CHEB': (poly.chebyshev.chebfit, poly.chebyshev.chebval),
505 'LEG': (poly.legendre.legfit, poly.legendre.legval),
509 }[self.config.fitType]
513 coeffs =
fitter(indices, collapsed, self.config.order)
515 if isinstance(coeffs, float):
516 self.log.
warn(
"Using fallback value %f due to fitter failure. Amplifier will be masked.",
518 overscanVector = np.full_like(indices, coeffs)
519 maskArray = np.full_like(collapsed,
True, dtype=bool)
522 overscanVector = evaler(indices, coeffs)
524 endTime = time.perf_counter()
525 self.log.
info(f
"Overscan measurement took {endTime - startTime}s for {self.config.fitType}")
526 return pipeBase.Struct(overscanValue=np.array(overscanVector),
528 isTransposed=isTransposed)
531 """Debug display for the final overscan solution.
536 Input image the overscan solution was determined from.
537 model : `numpy.ndarray`
or `float`
538 Overscan model determined
for the image.
540 Amplifier to extract diagnostic information.
549 calcImage, isTransposed = self.
transposetranspose(calcImage)
554 indices = 2.0 * np.arange(num)/
float(num) - 1.0
556 if np.ma.is_masked(collapsed):
557 collapsedMask = collapsed.mask
559 collapsedMask = np.array(num*[np.ma.nomask])
561 import matplotlib.pyplot
as plot
562 figure = plot.figure(1)
564 axes = figure.add_axes((0.1, 0.1, 0.8, 0.8))
565 axes.plot(indices[~collapsedMask], collapsed[~collapsedMask],
'k+')
566 if collapsedMask.sum() > 0:
567 axes.plot(indices[collapsedMask], collapsed.data[collapsedMask],
'b+')
568 if isinstance(model, np.ndarray):
571 plotModel = np.zeros_like(indices)
573 axes.plot(indices, plotModel,
'r-')
574 plot.xlabel(
"centered/scaled position along overscan region")
575 plot.ylabel(
"pixel value/fit value")
577 plot.title(f
"{amp.getName()} DataX: "
578 f
"[{amp.getRawDataBBox().getBeginX()}:{amp.getRawBBox().getEndX()}]"
579 f
"OscanX: [{amp.getRawHorizontalOverscanBBox().getBeginX()}:"
580 f
"{amp.getRawHorizontalOverscanBBox().getEndX()}] {self.config.fitType}")
582 plot.title(
"No amp supplied.")
584 prompt =
"Press Enter or c to continue [chp]..."
586 ans = input(prompt).lower()
587 if ans
in (
"",
" ",
"c",):
596 print(
"[h]elp [c]ontinue [p]db e[x]itDebug")
Geometry and electronic information about raw amplifier images.
A class to represent a 2-dimensional array of pixels.
Represent a 2-dimensional array of bitmask pixels.
A class to manipulate images, masks, and variance as a single object.
Pass parameters to a Statistics object.
def collapseArrayMedian(self, maskedArray)
def maskExtrapolated(collapsed)
def measureVectorOverscan(self, image)
def __init__(self, statControl=None, **kwargs)
def collapseArray(maskedArray)
def transpose(imageArray)
def maskOutliers(self, imageArray)
def debugView(self, image, model, amp=None)
def run(self, ampImage, overscanImage, amp=None)
def getImageArray(self, image)
def integerConvert(image)
def splineFit(self, indices, collapsed, numBins)
def splineEval(indices, interp)
def measureConstantOverscan(self, image)
Backwards-compatibility support for depersisting the old Calib (FluxMag0/FluxMag0Err) objects.
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)
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.
std::vector< double > fitOverscanImage(lsst::afw::image::MaskedImage< ImagePixelT > const &overscan, std::vector< std::string > badPixelMask, bool isTransposed)