LSST Applications 29.0.1,g0fba68d861+132dd21e0a,g107a963962+1bb9f809a9,g1fd858c14a+005be21cae,g21d47ad084+8a07b29876,g325378336f+5d73323c8f,g330003fc43+40b4eaffc6,g35bb328faa+fcb1d3bbc8,g36ff55ed5b+9c28a42a87,g4e0f332c67+5fbd1e3e73,g53246c7159+fcb1d3bbc8,g60b5630c4e+9c28a42a87,g67b6fd64d1+a38b34ea13,g78460c75b0+2f9a1b4bcd,g786e29fd12+cf7ec2a62a,g7b71ed6315+fcb1d3bbc8,g86c591e316+6b2b2d0295,g8852436030+bf14db0e33,g89139ef638+a38b34ea13,g8b8da53e10+e3777245af,g9125e01d80+fcb1d3bbc8,g989de1cb63+a38b34ea13,g9f1445be69+9c28a42a87,g9f33ca652e+52c8f07962,ga9baa6287d+9c28a42a87,ga9e4eb89a6+9f84bd6575,gabe3b4be73+1e0a283bba,gb037a4e798+f3cbcd26c0,gb1101e3267+e7be8da0f8,gb58c049af0+f03b321e39,gb89ab40317+a38b34ea13,gcf25f946ba+bf14db0e33,gd6cbbdb0b4+bce7f7457e,gd9a9a58781+fcb1d3bbc8,gde0f65d7ad+53d424b1ae,ge278dab8ac+222406d50a,ge410e46f29+a38b34ea13,ge80e9994a3+664d6357dc,gf67bdafdda+a38b34ea13
LSST Data Management Base Package
Loading...
Searching...
No Matches
lsst.ip.isr.overscan.SerialOverscanCorrectionTask Class Reference
Inheritance diagram for lsst.ip.isr.overscan.SerialOverscanCorrectionTask:
lsst.ip.isr.overscan.OverscanCorrectionTaskBase

Public Member Functions

 run (self, exposure, amp, isTransposed=False)
 
 correctOverscan (self, exposure, amp, imageBBox, overscanBBox, isTransposed=True, leadingToSkip=0, trailingToSkip=0, overscanFraction=1.0, imageThreshold=np.inf, maskedRowColumnGrowSize=0, medianSmoothingKernel=0, medianSmoothingOutlierThreshold=np.inf)
 
 broadcastFitToImage (self, overscanValue, imageArray, transpose=False)
 
 trimOverscan (self, exposure, amp, bbox, skipLeading, skipTrailing, transpose=False)
 
 fitOverscan (self, overscanImage, isTransposed=False)
 
 maskParallelOverscan (self, exposure, detector)
 
 measureConstantOverscan (self, image)
 
 getImageArray (self, image)
 
 maskOutliers (self, imageArray)
 
 fillMaskedPixels (self, overscanVector)
 
 collapseArray (self, maskedArray, fillMasked=True)
 
 splineFit (self, indices, collapsed, numBins)
 
 measureVectorOverscan (self, image, isTransposed=False)
 
 debugView (self, image, model, amp=None, isTransposed=True)
 

Static Public Member Functions

 splineEval (indices, interp)
 
 maskExtrapolated (collapsed)
 

Public Attributes

bool allowDebug = True
 
 statControl = statControl
 

Static Public Attributes

 ConfigClass = OverscanCorrectionTaskConfigBase
 

Static Protected Member Functions

 _maskRowsOrColumns (exposure, overscanBBox, overscanMaskedImage, overscanMask, maxDeviation, maskedRowColumnGrowSize, medianSmoothingKernel, medianSmoothingOutlierThreshold, doAbsoluteMaxDeviation, isTransposed)
 

Static Protected Attributes

str _DefaultName = "overscanBase"
 

Detailed Description

Correction task for serial overscan.

Parameters
----------
statControl : `lsst.afw.math.StatisticsControl`, optional
    Statistics control object.

Definition at line 1225 of file overscan.py.

Member Function Documentation

◆ _maskRowsOrColumns()

lsst.ip.isr.overscan.OverscanCorrectionTaskBase._maskRowsOrColumns ( exposure,
overscanBBox,
overscanMaskedImage,
overscanMask,
maxDeviation,
maskedRowColumnGrowSize,
medianSmoothingKernel,
medianSmoothingOutlierThreshold,
doAbsoluteMaxDeviation,
isTransposed )
staticprotectedinherited
Mask overscan rows (~serial) or columns (~parallel).

Parameters
----------
exposure : `lsst.afw.image.Exposure`
    Exposure containing the data.
overscanBBox: `lsst.geom.Box2I`
    Bounding box for the overscan data.
overscanMaskedImage : `lsst.afw.image.MaskedImage`
    Masked image containing the overscan data.
overscanMask : `np.ndarray`
    Numpy array of the mask bits, anded with appropriate
    mask planes.
maxDeviation : `float`
    Maximum deviation from median (overscan units) to mask in overscan
    correction. For parallel overscan this is a one-sided (positive
    only) cut.
maskedRowColumnGrowSize : `int`
    If a column (parallel overscan) or row (serial overscan) is
    completely masked, then grow the mask by this radius. If the
    value is <=0 then this will not be checked.
medianSmoothingKernel : `int`
    Kernel (pixels) to smooth rows/columns. If <=0, median smoothing
    is skipped. Otherwise must be odd.
medianSmoothingOutlierThreshold : `float`
    Outlier threshold after median smoothing (overscan units). This
    is applied only to positive outliers.
doAbsoluteMaxDeviation : `bool`
    If true, deviation comparisons will use the absolute value;
    otherwise it will cut positive outliers only.
isTransposed : `bool`
    If true, then the data will be transposed before fitting
    the overscan.

Returns
-------
badRowsOrColumns : `np.ndarray`
    Array of bad rows (serial) or columns (parallel) that were
    found, prior to dilation by maskedRowColumnGrowSize.

Definition at line 167 of file overscan.py.

178 ):
179 """Mask overscan rows (~serial) or columns (~parallel).
180
181 Parameters
182 ----------
183 exposure : `lsst.afw.image.Exposure`
184 Exposure containing the data.
185 overscanBBox: `lsst.geom.Box2I`
186 Bounding box for the overscan data.
187 overscanMaskedImage : `lsst.afw.image.MaskedImage`
188 Masked image containing the overscan data.
189 overscanMask : `np.ndarray`
190 Numpy array of the mask bits, anded with appropriate
191 mask planes.
192 maxDeviation : `float`
193 Maximum deviation from median (overscan units) to mask in overscan
194 correction. For parallel overscan this is a one-sided (positive
195 only) cut.
196 maskedRowColumnGrowSize : `int`
197 If a column (parallel overscan) or row (serial overscan) is
198 completely masked, then grow the mask by this radius. If the
199 value is <=0 then this will not be checked.
200 medianSmoothingKernel : `int`
201 Kernel (pixels) to smooth rows/columns. If <=0, median smoothing
202 is skipped. Otherwise must be odd.
203 medianSmoothingOutlierThreshold : `float`
204 Outlier threshold after median smoothing (overscan units). This
205 is applied only to positive outliers.
206 doAbsoluteMaxDeviation : `bool`
207 If true, deviation comparisons will use the absolute value;
208 otherwise it will cut positive outliers only.
209 isTransposed : `bool`
210 If true, then the data will be transposed before fitting
211 the overscan.
212
213 Returns
214 -------
215 badRowsOrColumns : `np.ndarray`
216 Array of bad rows (serial) or columns (parallel) that were
217 found, prior to dilation by maskedRowColumnGrowSize.
218 """
219 overscanArray = overscanMaskedImage.image.array
220
221 badRowsOrColumns = np.zeros(0, dtype=np.int64)
222
223 median = np.ma.median(np.ma.masked_where(overscanMask, overscanArray))
224 if doAbsoluteMaxDeviation:
225 delta = np.abs(overscanArray - median)
226 else:
227 delta = overscanArray - median
228
229 bad = np.where((delta > maxDeviation) & (~overscanMask))
230 # Mark the bad pixels as BAD
231 overscanMaskedImage.mask.array[bad] |= overscanMaskedImage.mask.getPlaneBitMask("BAD")
232
233 if isTransposed:
234 axis = 0
235 nComp = overscanArray.shape[0]
236 else:
237 axis = 1
238 nComp = overscanArray.shape[1]
239
240 # Check for completely masked row/column (from maxDeviation or
241 # previously applied SAT flag.)
242 if len(bad) > 0:
243 # We only need to look at the bad pixels set here for this
244 # mask growth.
245 overscanMaskTemp = np.zeros_like(overscanMask)
246 overscanMaskTemp[bad] = True
247
248 nMaskedArray = np.sum(overscanMaskTemp, axis=axis, dtype=np.int32)
249 badRowsOrColumns, = np.where(nMaskedArray == nComp)
250
251 # Perform median-smoothing outlier rejection if desired.
252 if medianSmoothingKernel > 0:
253 # We do a straight numpy median ignoring the mask.
254 # This will be fine because it avoids missing values,
255 # and very large deviations have already been flagged by
256 # maxDeviation or SAT.
257 rowsCols = np.median(overscanArray, axis=axis)
258 filtered = medfilt(rowsCols, kernel_size=medianSmoothingKernel)
259 delta = rowsCols - filtered
260
261 # We cannot reliably look for outliers within a kernel length
262 # of the edges.
263 high, = np.where(delta[medianSmoothingKernel: -medianSmoothingKernel]
264 >= medianSmoothingOutlierThreshold)
265 high += medianSmoothingKernel
266
267 if len(high) > 0:
268 badRowsOrColumns = np.unique(np.append(badRowsOrColumns, high))
269
270 # If we have any bad rows/columns, we need to dilate them
271 # and apply the mask to the parent overscan image.
272 if len(badRowsOrColumns) > 0:
273 dataView = afwImage.MaskedImageF(exposure.maskedImage,
274 overscanBBox,
275 afwImage.PARENT)
276 if isTransposed:
277 pixelsCopy = dataView.image.array[:, badRowsOrColumns].copy()
278 dataView.image.array[:, badRowsOrColumns] = 1e30
279 else:
280 pixelsCopy = dataView.image.array[badRowsOrColumns, :].copy()
281 dataView.image.array[badRowsOrColumns, :] = 1e30
282
283 makeThresholdMask(
284 maskedImage=dataView,
285 threshold=1e30,
286 growFootprints=maskedRowColumnGrowSize,
287 maskName="BAD",
288 )
289
290 if isTransposed:
291 dataView.image.array[:, badRowsOrColumns] = pixelsCopy
292 else:
293 dataView.image.array[badRowsOrColumns, :] = pixelsCopy
294
295 return badRowsOrColumns
296

◆ broadcastFitToImage()

lsst.ip.isr.overscan.OverscanCorrectionTaskBase.broadcastFitToImage ( self,
overscanValue,
imageArray,
transpose = False )
inherited
Broadcast 0 or 1 dimension fit to appropriate shape.

Parameters
----------
overscanValue : `numpy.ndarray`, (Nrows, ) or scalar
    Overscan fit to broadcast.
imageArray : `numpy.ndarray`, (Nrows, Ncols)
    Image array that we want to match.
transpose : `bool`, optional
    Switch order to broadcast along the other axis.

Returns
-------
overscanModel : `numpy.ndarray`, (Nrows, Ncols) or scalar
    Expanded overscan fit.

Raises
------
RuntimeError
    Raised if no axis has the appropriate dimension.

Definition at line 485 of file overscan.py.

485 def broadcastFitToImage(self, overscanValue, imageArray, transpose=False):
486 """Broadcast 0 or 1 dimension fit to appropriate shape.
487
488 Parameters
489 ----------
490 overscanValue : `numpy.ndarray`, (Nrows, ) or scalar
491 Overscan fit to broadcast.
492 imageArray : `numpy.ndarray`, (Nrows, Ncols)
493 Image array that we want to match.
494 transpose : `bool`, optional
495 Switch order to broadcast along the other axis.
496
497 Returns
498 -------
499 overscanModel : `numpy.ndarray`, (Nrows, Ncols) or scalar
500 Expanded overscan fit.
501
502 Raises
503 ------
504 RuntimeError
505 Raised if no axis has the appropriate dimension.
506 """
507 if isinstance(overscanValue, np.ndarray):
508 overscanModel = np.zeros_like(imageArray)
509
510 if transpose is False:
511 if imageArray.shape[0] == overscanValue.shape[0]:
512 overscanModel[:, :] = overscanValue[:, np.newaxis]
513 elif imageArray.shape[1] == overscanValue.shape[0]:
514 overscanModel[:, :] = overscanValue[np.newaxis, :]
515 elif imageArray.shape[0] == overscanValue.shape[1]:
516 overscanModel[:, :] = overscanValue[np.newaxis, :]
517 else:
518 raise RuntimeError(f"Could not broadcast {overscanValue.shape} to "
519 f"match {imageArray.shape}")
520 else:
521 if imageArray.shape[1] == overscanValue.shape[0]:
522 overscanModel[:, :] = overscanValue[np.newaxis, :]
523 elif imageArray.shape[0] == overscanValue.shape[0]:
524 overscanModel[:, :] = overscanValue[:, np.newaxis]
525 elif imageArray.shape[1] == overscanValue.shape[1]:
526 overscanModel[:, :] = overscanValue[:, np.newaxis]
527 else:
528 raise RuntimeError(f"Could not broadcast {overscanValue.shape} to "
529 f"match {imageArray.shape}")
530 else:
531 overscanModel = overscanValue
532
533 return overscanModel
534

◆ collapseArray()

lsst.ip.isr.overscan.OverscanCorrectionTaskBase.collapseArray ( self,
maskedArray,
fillMasked = True )
inherited
Collapse overscan array (and mask) to a 1-D vector of values.

Parameters
----------
maskedArray : `numpy.ma.masked_array`
    Masked array of input overscan data.
fillMasked : `bool`, optional
    If true, fill any pixels that are masked with a median of
    neighbors.

Returns
-------
collapsed : `numpy.ma.masked_array`
    Single dimensional overscan data, combined with the mean.

Definition at line 763 of file overscan.py.

763 def collapseArray(self, maskedArray, fillMasked=True):
764 """Collapse overscan array (and mask) to a 1-D vector of values.
765
766 Parameters
767 ----------
768 maskedArray : `numpy.ma.masked_array`
769 Masked array of input overscan data.
770 fillMasked : `bool`, optional
771 If true, fill any pixels that are masked with a median of
772 neighbors.
773
774 Returns
775 -------
776 collapsed : `numpy.ma.masked_array`
777 Single dimensional overscan data, combined with the mean.
778
779 """
780 collapsed = np.mean(maskedArray, axis=1)
781 if collapsed.mask.sum() > 0 and fillMasked:
782 collapsed = self.fillMaskedPixels(collapsed)
783
784 return collapsed
785

◆ correctOverscan()

lsst.ip.isr.overscan.OverscanCorrectionTaskBase.correctOverscan ( self,
exposure,
amp,
imageBBox,
overscanBBox,
isTransposed = True,
leadingToSkip = 0,
trailingToSkip = 0,
overscanFraction = 1.0,
imageThreshold = np.inf,
maskedRowColumnGrowSize = 0,
medianSmoothingKernel = 0,
medianSmoothingOutlierThreshold = np.inf )
inherited
Trim the exposure, fit the overscan, subtract the fit, and
calculate statistics.

Parameters
----------
exposure : `lsst.afw.image.Exposure`
    Exposure containing the data.
amp : `lsst.afw.cameraGeom.Amplifier`
    The amplifier that is to be corrected.
imageBBox: `lsst.geom.Box2I`
    Bounding box of the image data that will have the overscan
    subtracted.  If parallel overscan will be performed, that
    area is added to the image bounding box during serial
    overscan correction.
overscanBBox: `lsst.geom.Box2I`
    Bounding box for the overscan data.
isTransposed: `bool`
    If true, then the data will be transposed before fitting
    the overscan.
leadingToSkip : `int`, optional
    Leading rows/columns to skip.
trailingToSkip : `int`, optional
    Leading rows/columns to skip.
overscanFraction : `float`, optional
    If the overscan region median is greater than overscanFraction
    and the imaging region median is greater than imageThreshold
    then overscan correction will be skipped.
maxLevel : `float`, optional
    If the overscan region median is greater than overscanFraction
    and the imaging region median is greater than imageThreshold
    then overscan correction will be skipped.
maskedRowColumnGrowSize : `int`, optional
    If a column (parallel overscan) or row (serial overscan) is
    completely masked, then grow the mask by this radius. If the
    value is <=0 then this will not be checked.
medianSmoothingKernel : `int`, optional
    Kernel (pixels) to use to smooth rows/columns for row/column
    outlier rejection. Must be odd if positive; if <=0 median
    smoothing will not be used to find outliers.
medianSmoothingOutlierThreshold : `float`, optional
    Threshold to look for outliers after median smoothing (adu).

Returns
-------
results : `lsst.pipe.base.Struct`
    ``ampOverscanModel``
        Overscan model broadcast to the full image size.
        (`lsst.afw.image.Exposure`)
    ``overscanOverscanModel``
        Overscan model broadcast to the full overscan image
        size. (`lsst.afw.image.Exposure`)
    ``overscanImage``
        Overscan image with the overscan fit subtracted.
        (`lsst.afw.image.Exposure`)
    ``overscanValue``
        Overscan model. (`float` or `np.array`)
    ``overscanMean``
        Mean value of the overscan fit. (`float`)
    ``overscanMedian``
        Median value of the overscan fit. (`float`)
    ``overscanSigma``
        Standard deviation of the overscan fit. (`float`)
    ``overscanMeanResidual``
        Mean value of the overscan region after overscan
        subtraction. (`float`)
    ``overscanMedianResidual``
        Median value of the overscan region after overscan
        subtraction. (`float`)
    ``overscanSigmaResidual``
        Standard deviation of the overscan region after
        overscan subtraction. (`float`)

Definition at line 297 of file overscan.py.

311 ):
312 """Trim the exposure, fit the overscan, subtract the fit, and
313 calculate statistics.
314
315 Parameters
316 ----------
317 exposure : `lsst.afw.image.Exposure`
318 Exposure containing the data.
319 amp : `lsst.afw.cameraGeom.Amplifier`
320 The amplifier that is to be corrected.
321 imageBBox: `lsst.geom.Box2I`
322 Bounding box of the image data that will have the overscan
323 subtracted. If parallel overscan will be performed, that
324 area is added to the image bounding box during serial
325 overscan correction.
326 overscanBBox: `lsst.geom.Box2I`
327 Bounding box for the overscan data.
328 isTransposed: `bool`
329 If true, then the data will be transposed before fitting
330 the overscan.
331 leadingToSkip : `int`, optional
332 Leading rows/columns to skip.
333 trailingToSkip : `int`, optional
334 Leading rows/columns to skip.
335 overscanFraction : `float`, optional
336 If the overscan region median is greater than overscanFraction
337 and the imaging region median is greater than imageThreshold
338 then overscan correction will be skipped.
339 maxLevel : `float`, optional
340 If the overscan region median is greater than overscanFraction
341 and the imaging region median is greater than imageThreshold
342 then overscan correction will be skipped.
343 maskedRowColumnGrowSize : `int`, optional
344 If a column (parallel overscan) or row (serial overscan) is
345 completely masked, then grow the mask by this radius. If the
346 value is <=0 then this will not be checked.
347 medianSmoothingKernel : `int`, optional
348 Kernel (pixels) to use to smooth rows/columns for row/column
349 outlier rejection. Must be odd if positive; if <=0 median
350 smoothing will not be used to find outliers.
351 medianSmoothingOutlierThreshold : `float`, optional
352 Threshold to look for outliers after median smoothing (adu).
353
354 Returns
355 -------
356 results : `lsst.pipe.base.Struct`
357 ``ampOverscanModel``
358 Overscan model broadcast to the full image size.
359 (`lsst.afw.image.Exposure`)
360 ``overscanOverscanModel``
361 Overscan model broadcast to the full overscan image
362 size. (`lsst.afw.image.Exposure`)
363 ``overscanImage``
364 Overscan image with the overscan fit subtracted.
365 (`lsst.afw.image.Exposure`)
366 ``overscanValue``
367 Overscan model. (`float` or `np.array`)
368 ``overscanMean``
369 Mean value of the overscan fit. (`float`)
370 ``overscanMedian``
371 Median value of the overscan fit. (`float`)
372 ``overscanSigma``
373 Standard deviation of the overscan fit. (`float`)
374 ``overscanMeanResidual``
375 Mean value of the overscan region after overscan
376 subtraction. (`float`)
377 ``overscanMedianResidual``
378 Median value of the overscan region after overscan
379 subtraction. (`float`)
380 ``overscanSigmaResidual``
381 Standard deviation of the overscan region after
382 overscan subtraction. (`float`)
383 """
384 overscanBox = self.trimOverscan(exposure, amp, overscanBBox,
385 leadingToSkip,
386 trailingToSkip,
387 transpose=isTransposed)
388 overscanImage = exposure[overscanBox].getMaskedImage()
389
390 # Record the gain value if necessary to convert configs from
391 # electron to adu.
392 if exposure.metadata.get("LSST ISR UNITS", "adu") == "electron":
393 gain = exposure.metadata[f"LSST ISR GAIN {amp.getName()}"]
394 else:
395 gain = 1.0
396
397 # Mask pixels.
398 maskVal = overscanImage.mask.getPlaneBitMask(self.config.maskPlanes)
399 overscanMask = ~((overscanImage.mask.array & maskVal) == 0)
400
401 badResults = False
402 overscanMedian = np.nanmedian(overscanImage.image.array)
403 imageMedian = np.nanmedian(exposure[imageBBox].image.array)
404
405 if np.all(overscanMask):
406 self.log.warning(
407 "All overscan pixels masked when attempting overscan correction for %s",
408 amp.getName(),
409 )
410 badResults = True
411 elif overscanMedian/imageMedian > overscanFraction and imageMedian > imageThreshold:
412 self.log.warning(
413 "The level in the overscan region (%.2f) compared to the image region (%.2f) is "
414 "greater than the maximum fraction (%.2f) for %s",
415 overscanMedian,
416 imageMedian,
417 overscanFraction,
418 amp.getName(),
419 )
420 badResults = True
421
422 if badResults:
423 # Do not do overscan subtraction at all.
424 badRowsOrColumns = np.zeros(0, dtype=np.int64)
425 overscanResults = pipeBase.Struct(
426 overscanValue=0.0,
427 overscanMean=0.0,
428 overscanMedian=0.0,
429 overscanSigma=0.0,
430 )
431 else:
432 badRowsOrColumns = self._maskRowsOrColumns(
433 exposure,
434 overscanBBox,
435 overscanImage,
436 overscanMask,
437 gain * self.config.maxDeviation,
438 maskedRowColumnGrowSize,
439 medianSmoothingKernel,
440 gain * medianSmoothingOutlierThreshold,
441 self.config.doAbsoluteMaxDeviation,
442 isTransposed,
443 )
444 # Do overscan fit.
445 # CZW: Handle transposed correctly.
446 overscanResults = self.fitOverscan(overscanImage, isTransposed=isTransposed)
447
448 # Correct image region (and possibly parallel-overscan region).
449 ampImage = exposure[imageBBox]
450 ampOverscanModel = self.broadcastFitToImage(overscanResults.overscanValue,
451 ampImage.image.array,
452 transpose=isTransposed)
453 ampImage.image.array -= ampOverscanModel
454
455 # Correct overscan region (and possibly doubly-overscaned
456 # region).
457 overscanImage = exposure[overscanBBox]
458 # CZW: Transposed?
459 overscanOverscanModel = self.broadcastFitToImage(overscanResults.overscanValue,
460 overscanImage.image.array)
461 self.debugView(overscanImage, overscanResults.overscanValue, amp, isTransposed=isTransposed)
462 overscanImage.image.array -= overscanOverscanModel
463
464 # Find residual fit statistics.
465 stats = afwMath.makeStatistics(overscanImage.getMaskedImage(),
466 afwMath.MEAN | afwMath.MEDIAN | afwMath.STDEVCLIP, self.statControl)
467 residualMean = stats.getValue(afwMath.MEAN)
468 residualMedian = stats.getValue(afwMath.MEDIAN)
469 residualSigma = stats.getValue(afwMath.STDEVCLIP)
470
471 return pipeBase.Struct(
472 ampOverscanModel=ampOverscanModel,
473 overscanOverscanModel=overscanOverscanModel,
474 overscanImage=overscanImage,
475 overscanValue=overscanResults.overscanValue,
476 overscanMean=overscanResults.overscanMean,
477 overscanMedian=overscanResults.overscanMedian,
478 overscanSigma=overscanResults.overscanSigma,
479 overscanMeanResidual=residualMean,
480 overscanMedianResidual=residualMedian,
481 overscanSigmaResidual=residualSigma,
482 overscanBadRowsOrColumns=badRowsOrColumns,
483 )
484
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

◆ debugView()

lsst.ip.isr.overscan.OverscanCorrectionTaskBase.debugView ( self,
image,
model,
amp = None,
isTransposed = True )
inherited
Debug display for the final overscan solution.

Parameters
----------
image : `lsst.afw.image.Image`
    Input image the overscan solution was determined from.
model : `numpy.ndarray` or `float`
    Overscan model determined for the image.
amp : `lsst.afw.cameraGeom.Amplifier`, optional
    Amplifier to extract diagnostic information.
isTransposed : `bool`, optional
    Does the data need to be transposed before display?

Definition at line 967 of file overscan.py.

967 def debugView(self, image, model, amp=None, isTransposed=True):
968 """Debug display for the final overscan solution.
969
970 Parameters
971 ----------
972 image : `lsst.afw.image.Image`
973 Input image the overscan solution was determined from.
974 model : `numpy.ndarray` or `float`
975 Overscan model determined for the image.
976 amp : `lsst.afw.cameraGeom.Amplifier`, optional
977 Amplifier to extract diagnostic information.
978 isTransposed : `bool`, optional
979 Does the data need to be transposed before display?
980 """
981 import lsstDebug
982 if not lsstDebug.Info(__name__).display:
983 return
984 if not self.allowDebug:
985 return
986
987 calcImage = self.getImageArray(image)
988 # CZW: Check that this is ok
989 if isTransposed:
990 calcImage = np.transpose(calcImage)
991 masked = self.maskOutliers(calcImage)
992 collapsed = self.collapseArray(masked, fillMasked=False)
993
994 num = len(collapsed)
995 indices = 2.0 * np.arange(num)/float(num) - 1.0
996 indices = np.arange(num)
997
998 if np.ma.is_masked(collapsed):
999 collapsedMask = collapsed.mask
1000 else:
1001 collapsedMask = np.array(num*[np.ma.nomask])
1002
1003 import matplotlib.pyplot as plot
1004 figure = plot.figure(1)
1005 figure.clear()
1006 axes = figure.add_axes((0.1, 0.1, 0.8, 0.8))
1007 axes.plot(indices[~collapsedMask], collapsed[~collapsedMask], 'k+')
1008 if collapsedMask.sum() > 0:
1009 axes.plot(indices[collapsedMask], collapsed.data[collapsedMask], 'b+')
1010 if isinstance(model, np.ndarray):
1011 plotModel = model
1012 else:
1013 plotModel = np.zeros_like(indices)
1014 plotModel += model
1015
1016 axes.plot(indices, plotModel, 'r-')
1017 plot.xlabel("position along overscan region")
1018 plot.ylabel("pixel value/fit value")
1019 if amp:
1020 plot.title(f"{amp.getName()} DataX: "
1021 f"[{amp.getRawDataBBox().getBeginX()}:{amp.getRawBBox().getEndX()}]"
1022 f"OscanX: [{amp.getRawHorizontalOverscanBBox().getBeginX()}:"
1023 f"{amp.getRawHorizontalOverscanBBox().getEndX()}] {self.config.fitType}")
1024 else:
1025 plot.title("No amp supplied.")
1026 figure.show()
1027 prompt = "Press Enter or c to continue [chp]..."
1028 while True:
1029 ans = input(prompt).lower()
1030 if ans in ("", " ", "c",):
1031 break
1032 elif ans in ("p", ):
1033 import pdb
1034 pdb.set_trace()
1035 elif ans in ('x', ):
1036 self.allowDebug = False
1037 break
1038 elif ans in ("h", ):
1039 print("[h]elp [c]ontinue [p]db e[x]itDebug")
1040 plot.close()
1041
1042

◆ fillMaskedPixels()

lsst.ip.isr.overscan.OverscanCorrectionTaskBase.fillMaskedPixels ( self,
overscanVector )
inherited
Fill masked/NaN pixels in the overscan.

Parameters
----------
overscanVector : `np.array` or `np.ma.masked_array`
    Overscan vector to fill.

Returns
-------
overscanVector : `np.ma.masked_array`
    Filled vector.

Notes
-----
Each maskSlice is a section of overscan with contiguous masks.
Ideally this adds 5 pixels from the left and right of that
mask slice, and takes the median of those values to fill the
slice.  If this isn't possible, the median of all non-masked
values is used.  The mask is removed for the pixels filled.

Definition at line 722 of file overscan.py.

722 def fillMaskedPixels(self, overscanVector):
723 """Fill masked/NaN pixels in the overscan.
724
725 Parameters
726 ----------
727 overscanVector : `np.array` or `np.ma.masked_array`
728 Overscan vector to fill.
729
730 Returns
731 -------
732 overscanVector : `np.ma.masked_array`
733 Filled vector.
734
735 Notes
736 -----
737 Each maskSlice is a section of overscan with contiguous masks.
738 Ideally this adds 5 pixels from the left and right of that
739 mask slice, and takes the median of those values to fill the
740 slice. If this isn't possible, the median of all non-masked
741 values is used. The mask is removed for the pixels filled.
742 """
743 workingCopy = overscanVector
744 if not isinstance(overscanVector, np.ma.MaskedArray):
745 workingCopy = np.ma.masked_array(overscanVector,
746 mask=~np.isfinite(overscanVector))
747
748 defaultValue = np.median(workingCopy.data[~workingCopy.mask])
749 for maskSlice in np.ma.clump_masked(workingCopy):
750 neighborhood = []
751 if maskSlice.start > 5:
752 neighborhood.extend(workingCopy[maskSlice.start - 5:maskSlice.start].data)
753 if maskSlice.stop < workingCopy.size - 5:
754 neighborhood.extend(workingCopy[maskSlice.stop:maskSlice.stop+5].data)
755 if len(neighborhood) > 0:
756 workingCopy.data[maskSlice] = np.nanmedian(neighborhood)
757 workingCopy.mask[maskSlice] = False
758 else:
759 workingCopy.data[maskSlice] = defaultValue
760 workingCopy.mask[maskSlice] = False
761 return workingCopy
762

◆ fitOverscan()

lsst.ip.isr.overscan.OverscanCorrectionTaskBase.fitOverscan ( self,
overscanImage,
isTransposed = False )
inherited

Definition at line 582 of file overscan.py.

582 def fitOverscan(self, overscanImage, isTransposed=False):
583 if self.config.fitType in ('MEAN', 'MEANCLIP', 'MEDIAN'):
584 # Transposition has no effect here.
585 overscanResult = self.measureConstantOverscan(overscanImage)
586 overscanValue = overscanResult.overscanValue
587 overscanMean = overscanValue
588 overscanMedian = overscanValue
589 overscanSigma = 0.0
590 elif self.config.fitType in ('MEDIAN_PER_ROW', 'MEAN_PER_ROW', 'POLY', 'CHEB', 'LEG',
591 'NATURAL_SPLINE', 'CUBIC_SPLINE', 'AKIMA_SPLINE'):
592 # Force transposes as needed
593 overscanResult = self.measureVectorOverscan(overscanImage, isTransposed)
594 overscanValue = overscanResult.overscanValue
595
596 stats = afwMath.makeStatistics(overscanResult.overscanValue,
597 afwMath.MEAN | afwMath.MEDIAN | afwMath.STDEVCLIP,
598 self.statControl)
599 overscanMean = stats.getValue(afwMath.MEAN)
600 overscanMedian = stats.getValue(afwMath.MEDIAN)
601 overscanSigma = stats.getValue(afwMath.STDEVCLIP)
602 else:
603 raise ValueError('%s : %s an invalid overscan type' %
604 ("overscanCorrection", self.config.fitType))
605
606 return pipeBase.Struct(overscanValue=overscanValue,
607 overscanMean=overscanMean,
608 overscanMedian=overscanMedian,
609 overscanSigma=overscanSigma,
610 )
611

◆ getImageArray()

lsst.ip.isr.overscan.OverscanCorrectionTaskBase.getImageArray ( self,
image )
inherited
Extract the numpy array from the input image.

Parameters
----------
image : `lsst.afw.image.Image` or `lsst.afw.image.MaskedImage`
    Image data to pull array from.

calcImage : `numpy.ndarray`
    Image data array for numpy operating.

Definition at line 671 of file overscan.py.

671 def getImageArray(self, image):
672 """Extract the numpy array from the input image.
673
674 Parameters
675 ----------
676 image : `lsst.afw.image.Image` or `lsst.afw.image.MaskedImage`
677 Image data to pull array from.
678
679 calcImage : `numpy.ndarray`
680 Image data array for numpy operating.
681 """
682 if hasattr(image, "getImage"):
683 calcImage = image.getImage().getArray()
684 calcImage = np.ma.masked_where(image.getMask().getArray() & self.statControl.getAndMask(),
685 calcImage)
686 else:
687 calcImage = image.getArray()
688 return calcImage
689

◆ maskExtrapolated()

lsst.ip.isr.overscan.OverscanCorrectionTaskBase.maskExtrapolated ( collapsed)
staticinherited
Create mask if edges are extrapolated.

Parameters
----------
collapsed : `numpy.ma.masked_array`
    Masked array to check the edges of.

Returns
-------
maskArray : `numpy.ndarray`
    Boolean numpy array of pixels to mask.

Definition at line 852 of file overscan.py.

852 def maskExtrapolated(collapsed):
853 """Create mask if edges are extrapolated.
854
855 Parameters
856 ----------
857 collapsed : `numpy.ma.masked_array`
858 Masked array to check the edges of.
859
860 Returns
861 -------
862 maskArray : `numpy.ndarray`
863 Boolean numpy array of pixels to mask.
864 """
865 maskArray = np.full_like(collapsed, False, dtype=bool)
866 if np.ma.is_masked(collapsed):
867 num = len(collapsed)
868 for low in range(num):
869 if not collapsed.mask[low]:
870 break
871 if low > 0:
872 maskArray[:low] = True
873 for high in range(1, num):
874 if not collapsed.mask[-high]:
875 break
876 if high > 1:
877 maskArray[-high:] = True
878 return maskArray
879

◆ maskOutliers()

lsst.ip.isr.overscan.OverscanCorrectionTaskBase.maskOutliers ( self,
imageArray )
inherited
Mask  outliers in  a  row  of overscan  data  from  a robust  sigma
clipping procedure.

Parameters
----------
imageArray : `numpy.ndarray`
    Image to filter along numpy axis=1.

Returns
-------
maskedArray : `numpy.ma.masked_array`
    Masked image marking outliers.

Definition at line 690 of file overscan.py.

690 def maskOutliers(self, imageArray):
691 """Mask outliers in a row of overscan data from a robust sigma
692 clipping procedure.
693
694 Parameters
695 ----------
696 imageArray : `numpy.ndarray`
697 Image to filter along numpy axis=1.
698
699 Returns
700 -------
701 maskedArray : `numpy.ma.masked_array`
702 Masked image marking outliers.
703 """
704 lq, median, uq = np.percentile(np.ma.getdata(imageArray),
705 [25.0, 50.0, 75.0], axis=1)
706 axisMedians = median
707 axisStdev = 0.74*(uq - lq) # robust stdev
708
709 # Replace pixels that have excessively large stdev values
710 # with the median of stdev values. A large stdev likely
711 # indicates a bleed is spilling into the overscan.
712 axisStdev = np.where(axisStdev > 2.0 * np.median(axisStdev),
713 np.median(axisStdev), axisStdev)
714
715 # Mask pixels that are N-sigma away from their array medians.
716 diff = np.abs(imageArray - axisMedians[:, np.newaxis])
717 masked = np.ma.masked_where(diff > self.statControl.getNumSigmaClip()
718 * axisStdev[:, np.newaxis], imageArray)
719
720 return masked
721

◆ maskParallelOverscan()

lsst.ip.isr.overscan.OverscanCorrectionTaskBase.maskParallelOverscan ( self,
exposure,
detector )
inherited
Mask the union of high values on all amplifiers in the parallel
overscan.

This operates on the image in-place.

Parameters
----------
exposure : `lsst.afw.image.Exposure`
    An untrimmed raw exposure.
detector : `lsst.afw.cameraGeom.Detector`
    The detetor to use for amplifier geometry.

Definition at line 612 of file overscan.py.

612 def maskParallelOverscan(self, exposure, detector):
613 """Mask the union of high values on all amplifiers in the parallel
614 overscan.
615
616 This operates on the image in-place.
617
618 Parameters
619 ----------
620 exposure : `lsst.afw.image.Exposure`
621 An untrimmed raw exposure.
622 detector : `lsst.afw.cameraGeom.Detector`
623 The detetor to use for amplifier geometry.
624 """
625 parallelMask = None
626
627 for amp in detector:
628 dataView = afwImage.MaskedImageF(exposure.getMaskedImage(),
629 amp.getRawParallelOverscanBBox(),
630 afwImage.PARENT)
631 # This should mark all the saturated pixels as SAT.
632 makeThresholdMask(
633 maskedImage=dataView,
634 threshold=self.config.parallelOverscanMaskThreshold,
635 growFootprints=self.config.parallelOverscanMaskGrowSize,
636 maskName="SAT"
637 )
638 if parallelMask is None:
639 parallelMask = dataView.mask.array
640 else:
641 parallelMask |= dataView.mask.array
642 for amp in detector:
643 dataView = afwImage.MaskedImageF(exposure.getMaskedImage(),
644 amp.getRawParallelOverscanBBox(),
645 afwImage.PARENT)
646 dataView.mask.array |= parallelMask
647

◆ measureConstantOverscan()

lsst.ip.isr.overscan.OverscanCorrectionTaskBase.measureConstantOverscan ( self,
image )
inherited
Measure a constant overscan value.

Parameters
----------
image : `lsst.afw.image.Image` or `lsst.afw.image.MaskedImage`
    Image data to measure the overscan from.

Returns
-------
results : `lsst.pipe.base.Struct`
    Overscan result with entries:
    - ``overscanValue``: Overscan value to subtract (`float`)
    - ``isTransposed``: Orientation of the overscan (`bool`)

Definition at line 649 of file overscan.py.

649 def measureConstantOverscan(self, image):
650 """Measure a constant overscan value.
651
652 Parameters
653 ----------
654 image : `lsst.afw.image.Image` or `lsst.afw.image.MaskedImage`
655 Image data to measure the overscan from.
656
657 Returns
658 -------
659 results : `lsst.pipe.base.Struct`
660 Overscan result with entries:
661 - ``overscanValue``: Overscan value to subtract (`float`)
662 - ``isTransposed``: Orientation of the overscan (`bool`)
663 """
664 fitType = afwMath.stringToStatisticsProperty(self.config.fitType)
665 overscanValue = afwMath.makeStatistics(image, fitType, self.statControl).getValue()
666
667 return pipeBase.Struct(overscanValue=overscanValue,
668 isTransposed=False)
669
Property stringToStatisticsProperty(std::string const property)
Conversion function to switch a string to a Property (see Statistics.h)

◆ measureVectorOverscan()

lsst.ip.isr.overscan.OverscanCorrectionTaskBase.measureVectorOverscan ( self,
image,
isTransposed = False )
inherited
Calculate the 1-d vector overscan from the input overscan image.

Parameters
----------
image : `lsst.afw.image.MaskedImage`
    Image containing the overscan data.
isTransposed : `bool`
    If true, the image has been transposed.

Returns
-------
results : `lsst.pipe.base.Struct`
    Overscan result with entries:

    ``overscanValue``
        Overscan value to subtract (`float`)
    ``maskArray``
        List of rows that should be masked as ``SUSPECT`` when the
        overscan solution is applied. (`list` [ `bool` ])
    ``isTransposed``
       Indicates if the overscan data was transposed during
       calcuation, noting along which axis the overscan should be
       subtracted. (`bool`)

Definition at line 880 of file overscan.py.

880 def measureVectorOverscan(self, image, isTransposed=False):
881 """Calculate the 1-d vector overscan from the input overscan image.
882
883 Parameters
884 ----------
885 image : `lsst.afw.image.MaskedImage`
886 Image containing the overscan data.
887 isTransposed : `bool`
888 If true, the image has been transposed.
889
890 Returns
891 -------
892 results : `lsst.pipe.base.Struct`
893 Overscan result with entries:
894
895 ``overscanValue``
896 Overscan value to subtract (`float`)
897 ``maskArray``
898 List of rows that should be masked as ``SUSPECT`` when the
899 overscan solution is applied. (`list` [ `bool` ])
900 ``isTransposed``
901 Indicates if the overscan data was transposed during
902 calcuation, noting along which axis the overscan should be
903 subtracted. (`bool`)
904 """
905 calcImage = self.getImageArray(image)
906
907 # operate on numpy-arrays from here
908 if isTransposed:
909 calcImage = np.transpose(calcImage)
910 masked = self.maskOutliers(calcImage)
911
912 if self.config.fitType in ('MEDIAN_PER_ROW', "MEAN_PER_ROW"):
913 if self.config.overscanIsInt:
914 mi = afwImage.MaskedImageI(image.getBBox())
915 masked = masked.astype(int)
916 else:
917 mi = image.clone()
918
919 if isTransposed:
920 masked = masked.transpose()
921
922 mi.image.array[:, :] = masked.data[:, :]
923 if bool(masked.mask.shape):
924 mi.mask.array[:, :] = masked.mask[:, :]
925
926 if self.config.fitType == "MEDIAN_PER_ROW":
927 overscanVector = fitOverscanImage(mi, self.config.maskPlanes, isTransposed)
928 else:
929 overscanVector = fitOverscanImageMean(mi, self.config.maskPlanes, isTransposed)
930
931 overscanVector = self.fillMaskedPixels(overscanVector)
932 maskArray = self.maskExtrapolated(overscanVector)
933 else:
934 collapsed = self.collapseArray(masked)
935
936 num = len(collapsed)
937 indices = 2.0*np.arange(num)/float(num) - 1.0
938
939 poly = np.polynomial
940 fitter, evaler = {
941 'POLY': (poly.polynomial.polyfit, poly.polynomial.polyval),
942 'CHEB': (poly.chebyshev.chebfit, poly.chebyshev.chebval),
943 'LEG': (poly.legendre.legfit, poly.legendre.legval),
944 'NATURAL_SPLINE': (self.splineFit, self.splineEval),
945 'CUBIC_SPLINE': (self.splineFit, self.splineEval),
946 'AKIMA_SPLINE': (self.splineFit, self.splineEval)
947 }[self.config.fitType]
948
949 # These are the polynomial coefficients, or an
950 # interpolation object.
951 coeffs = fitter(indices, collapsed, self.config.order)
952
953 if isinstance(coeffs, float):
954 self.log.warning("Using fallback value %f due to fitter failure. Amplifier will be masked.",
955 coeffs)
956 overscanVector = np.full_like(indices, coeffs)
957 maskArray = np.full_like(collapsed, True, dtype=bool)
958 else:
959 # Otherwise we can just use things as normal.
960 overscanVector = evaler(indices, coeffs)
961 maskArray = self.maskExtrapolated(collapsed)
962
963 return pipeBase.Struct(overscanValue=np.array(overscanVector),
964 maskArray=maskArray,
965 isTransposed=isTransposed)
966

◆ run()

lsst.ip.isr.overscan.SerialOverscanCorrectionTask.run ( self,
exposure,
amp,
isTransposed = False )
Measure and remove serial overscan from an amplifier image.

Parameters
----------
exposure : `lsst.afw.image.Exposure`
    Image data that will have the overscan corrections applied.
amp : `lsst.afw.cameraGeom.Amplifier`
    Amplifier to use for debugging purposes.
isTransposed : `bool`, optional
    Is the image transposed, such that serial and parallel
    overscan regions are reversed?  Default is False.

Returns
-------
overscanResults : `lsst.pipe.base.Struct`
    Result struct with components:

    ``imageFit``
        Value or fit subtracted from the amplifier image data
        (scalar or `lsst.afw.image.Image`).
    ``overscanFit``
        Value or fit subtracted from the serial overscan image
        data (scalar or `lsst.afw.image.Image`).
    ``overscanImage``
        Image of the serial overscan region with the serial
        overscan correction applied
        (`lsst.afw.image.Image`). This quantity is used to
        estimate the amplifier read noise empirically.
    ``overscanMean``
        Mean of the fit serial overscan region.
    ``overscanMedian``
        Median of the fit serial overscan region.
    ``overscanSigma``
        Sigma of the fit serial overscan region.
    ``residualMean``
        Mean of the residual of the serial overscan region after
        correction.
    ``residualMedian``
        Median of the residual of the serial overscan region after
        correction.
    ``residualSigma``
        Mean of the residual of the serial overscan region after
        correction.

Raises
------
RuntimeError
    Raised if an invalid overscan type is set.

Reimplemented from lsst.ip.isr.overscan.OverscanCorrectionTaskBase.

Definition at line 1236 of file overscan.py.

1236 def run(self, exposure, amp, isTransposed=False):
1237 """Measure and remove serial overscan from an amplifier image.
1238
1239 Parameters
1240 ----------
1241 exposure : `lsst.afw.image.Exposure`
1242 Image data that will have the overscan corrections applied.
1243 amp : `lsst.afw.cameraGeom.Amplifier`
1244 Amplifier to use for debugging purposes.
1245 isTransposed : `bool`, optional
1246 Is the image transposed, such that serial and parallel
1247 overscan regions are reversed? Default is False.
1248
1249 Returns
1250 -------
1251 overscanResults : `lsst.pipe.base.Struct`
1252 Result struct with components:
1253
1254 ``imageFit``
1255 Value or fit subtracted from the amplifier image data
1256 (scalar or `lsst.afw.image.Image`).
1257 ``overscanFit``
1258 Value or fit subtracted from the serial overscan image
1259 data (scalar or `lsst.afw.image.Image`).
1260 ``overscanImage``
1261 Image of the serial overscan region with the serial
1262 overscan correction applied
1263 (`lsst.afw.image.Image`). This quantity is used to
1264 estimate the amplifier read noise empirically.
1265 ``overscanMean``
1266 Mean of the fit serial overscan region.
1267 ``overscanMedian``
1268 Median of the fit serial overscan region.
1269 ``overscanSigma``
1270 Sigma of the fit serial overscan region.
1271 ``residualMean``
1272 Mean of the residual of the serial overscan region after
1273 correction.
1274 ``residualMedian``
1275 Median of the residual of the serial overscan region after
1276 correction.
1277 ``residualSigma``
1278 Mean of the residual of the serial overscan region after
1279 correction.
1280
1281 Raises
1282 ------
1283 RuntimeError
1284 Raised if an invalid overscan type is set.
1285 """
1286 serialOverscanBBox = amp.getRawSerialOverscanBBox()
1287 imageBBox = amp.getRawDataBBox()
1288
1289 # We always want to extend the serial overscan bounding box to
1290 # the full size of the detector.
1291 parallelOverscanBBox = amp.getRawParallelOverscanBBox()
1292 imageBBox = imageBBox.expandedTo(parallelOverscanBBox)
1293
1294 if isTransposed:
1295 serialOverscanBBox = geom.Box2I(
1296 geom.Point2I(serialOverscanBBox.getMinX(), imageBBox.getEndY()),
1297 geom.Extent2I(imageBBox.getWidth(), serialOverscanBBox.getHeight()),
1298 )
1299 else:
1300 serialOverscanBBox = geom.Box2I(
1301 geom.Point2I(serialOverscanBBox.getMinX(),
1302 imageBBox.getMinY()),
1303 geom.Extent2I(serialOverscanBBox.getWidth(),
1304 imageBBox.getHeight()),
1305 )
1306
1307 results = self.correctOverscan(
1308 exposure,
1309 amp,
1310 imageBBox,
1311 serialOverscanBBox,
1312 isTransposed=isTransposed,
1313 leadingToSkip=self.config.leadingToSkip,
1314 trailingToSkip=self.config.trailingToSkip,
1315 )
1316 overscanMean = results.overscanMean
1317 overscanMedian = results.overscanMedian
1318 overscanSigma = results.overscanSigma
1319 residualMean = results.overscanMeanResidual
1320 residualMedian = results.overscanMedianResidual
1321 residualSigma = results.overscanSigmaResidual
1322 badRowsOrColumns = results.overscanBadRowsOrColumns
1323
1324 return pipeBase.Struct(
1325 imageFit=results.ampOverscanModel,
1326 overscanFit=results.overscanOverscanModel,
1327 overscanImage=results.overscanImage,
1328 overscanMean=overscanMean,
1329 overscanMedian=overscanMedian,
1330 overscanSigma=overscanSigma,
1331 residualMean=residualMean,
1332 residualMedian=residualMedian,
1333 residualSigma=residualSigma,
1334 overscanBadRowsOrColumns=badRowsOrColumns,
1335 )
1336
1337
An integer coordinate rectangle.
Definition Box.h:55

◆ splineEval()

lsst.ip.isr.overscan.OverscanCorrectionTaskBase.splineEval ( indices,
interp )
staticinherited
Wrapper function to match spline evaluation API to polynomial fit
API.

Parameters
----------
indices : `numpy.ndarray`
    Locations to evaluate the spline.
interp : `lsst.afw.math.interpolate`
    Interpolation object to use.

Returns
-------
values : `numpy.ndarray`
    Evaluated spline values at each index.

Definition at line 832 of file overscan.py.

832 def splineEval(indices, interp):
833 """Wrapper function to match spline evaluation API to polynomial fit
834 API.
835
836 Parameters
837 ----------
838 indices : `numpy.ndarray`
839 Locations to evaluate the spline.
840 interp : `lsst.afw.math.interpolate`
841 Interpolation object to use.
842
843 Returns
844 -------
845 values : `numpy.ndarray`
846 Evaluated spline values at each index.
847 """
848
849 return interp.interpolate(indices.astype(float))
850

◆ splineFit()

lsst.ip.isr.overscan.OverscanCorrectionTaskBase.splineFit ( self,
indices,
collapsed,
numBins )
inherited
Wrapper function to match spline fit API to polynomial fit API.

Parameters
----------
indices : `numpy.ndarray`
    Locations to evaluate the spline.
collapsed : `numpy.ndarray`
    Collapsed overscan values corresponding to the spline
    evaluation points.
numBins : `int`
    Number of bins to use in constructing the spline.

Returns
-------
interp : `lsst.afw.math.Interpolate`
    Interpolation object for later evaluation.

Definition at line 786 of file overscan.py.

786 def splineFit(self, indices, collapsed, numBins):
787 """Wrapper function to match spline fit API to polynomial fit API.
788
789 Parameters
790 ----------
791 indices : `numpy.ndarray`
792 Locations to evaluate the spline.
793 collapsed : `numpy.ndarray`
794 Collapsed overscan values corresponding to the spline
795 evaluation points.
796 numBins : `int`
797 Number of bins to use in constructing the spline.
798
799 Returns
800 -------
801 interp : `lsst.afw.math.Interpolate`
802 Interpolation object for later evaluation.
803 """
804 if not np.ma.is_masked(collapsed):
805 collapsed.mask = np.array(len(collapsed)*[np.ma.nomask])
806
807 numPerBin, binEdges = np.histogram(indices, bins=numBins,
808 weights=1 - collapsed.mask.astype(int))
809 with np.errstate(invalid="ignore"):
810 values = np.histogram(indices, bins=numBins,
811 weights=collapsed.data*~collapsed.mask)[0]/numPerBin
812 binCenters = np.histogram(indices, bins=numBins,
813 weights=indices*~collapsed.mask)[0]/numPerBin
814
815 if len(binCenters[numPerBin > 0]) < 5:
816 self.log.warning("Cannot do spline fitting for overscan: %s valid points.",
817 len(binCenters[numPerBin > 0]))
818 # Return a scalar value if we have one, otherwise
819 # return zero. This amplifier is hopefully already
820 # masked.
821 if len(values[numPerBin > 0]) != 0:
822 return float(values[numPerBin > 0][0])
823 else:
824 return 0.0
825
826 interp = afwMath.makeInterpolate(binCenters.astype(float)[numPerBin > 0],
827 values.astype(float)[numPerBin > 0],
828 afwMath.stringToInterpStyle(self.config.fitType))
829 return interp
830
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.

◆ trimOverscan()

lsst.ip.isr.overscan.OverscanCorrectionTaskBase.trimOverscan ( self,
exposure,
amp,
bbox,
skipLeading,
skipTrailing,
transpose = False )
inherited
Trim overscan region to remove edges.

Parameters
----------
exposure : `lsst.afw.image.Exposure`
    Exposure containing data.
amp : `lsst.afw.cameraGeom.Amplifier`
    Amplifier containing geometry information.
bbox : `lsst.geom.Box2I`
    Bounding box of the overscan region.
skipLeading : `int`
    Number of leading (towards data region) rows/columns to skip.
skipTrailing : `int`
    Number of trailing (away from data region) rows/columns to skip.
transpose : `bool`, optional
    Operate on the transposed array.

Returns
-------
overscanArray : `numpy.array`, (N, M)
    Data array to fit.
overscanMask : `numpy.array`, (N, M)
    Data mask.

Definition at line 535 of file overscan.py.

535 def trimOverscan(self, exposure, amp, bbox, skipLeading, skipTrailing, transpose=False):
536 """Trim overscan region to remove edges.
537
538 Parameters
539 ----------
540 exposure : `lsst.afw.image.Exposure`
541 Exposure containing data.
542 amp : `lsst.afw.cameraGeom.Amplifier`
543 Amplifier containing geometry information.
544 bbox : `lsst.geom.Box2I`
545 Bounding box of the overscan region.
546 skipLeading : `int`
547 Number of leading (towards data region) rows/columns to skip.
548 skipTrailing : `int`
549 Number of trailing (away from data region) rows/columns to skip.
550 transpose : `bool`, optional
551 Operate on the transposed array.
552
553 Returns
554 -------
555 overscanArray : `numpy.array`, (N, M)
556 Data array to fit.
557 overscanMask : `numpy.array`, (N, M)
558 Data mask.
559 """
560 dx0, dy0, dx1, dy1 = (0, 0, 0, 0)
561 dataBBox = amp.getRawDataBBox()
562 if transpose:
563 if dataBBox.getBeginY() < bbox.getBeginY():
564 dy0 += skipLeading
565 dy1 -= skipTrailing
566 else:
567 dy0 += skipTrailing
568 dy1 -= skipLeading
569 else:
570 if dataBBox.getBeginX() < bbox.getBeginX():
571 dx0 += skipLeading
572 dx1 -= skipTrailing
573 else:
574 dx0 += skipTrailing
575 dx1 -= skipLeading
576
577 overscanBBox = geom.Box2I(bbox.getBegin() + geom.Extent2I(dx0, dy0),
578 geom.Extent2I(bbox.getWidth() - dx0 + dx1,
579 bbox.getHeight() - dy0 + dy1))
580 return overscanBBox
581

Member Data Documentation

◆ _DefaultName

str lsst.ip.isr.overscan.OverscanCorrectionTaskBase._DefaultName = "overscanBase"
staticprotectedinherited

Definition at line 150 of file overscan.py.

◆ allowDebug

bool lsst.ip.isr.overscan.OverscanCorrectionTaskBase.allowDebug = True
inherited

Definition at line 154 of file overscan.py.

◆ ConfigClass

lsst.ip.isr.overscan.OverscanCorrectionTaskBase.ConfigClass = OverscanCorrectionTaskConfigBase
staticinherited

Definition at line 149 of file overscan.py.

◆ statControl

lsst.ip.isr.overscan.OverscanCorrectionTaskBase.statControl = statControl
inherited

Definition at line 157 of file overscan.py.


The documentation for this class was generated from the following file: