lsst.ip.isr  20.0.0-18-g206d22f+d2859fcb1f
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 import numpy as np
23 import lsst.afw.math as afwMath
24 import lsst.afw.image as afwImage
25 import lsst.pipe.base as pipeBase
26 import lsst.pex.config as pexConfig
27 
28 __all__ = ["OverscanCorrectionTaskConfig", "OverscanCorrectionTask"]
29 
30 
31 class OverscanCorrectionTaskConfig(pexConfig.Config):
32  """Overscan correction options.
33  """
34  fitType = pexConfig.ChoiceField(
35  dtype=str,
36  doc="The method for fitting the overscan bias level.",
37  default='MEDIAN',
38  allowed={
39  "POLY": "Fit ordinary polynomial to the longest axis of the overscan region",
40  "CHEB": "Fit Chebyshev polynomial to the longest axis of the overscan region",
41  "LEG": "Fit Legendre polynomial to the longest axis of the overscan region",
42  "NATURAL_SPLINE": "Fit natural spline to the longest axis of the overscan region",
43  "CUBIC_SPLINE": "Fit cubic spline to the longest axis of the overscan region",
44  "AKIMA_SPLINE": "Fit Akima spline to the longest axis of the overscan region",
45  "MEAN": "Correct using the mean of the overscan region",
46  "MEANCLIP": "Correct using a clipped mean of the overscan region",
47  "MEDIAN": "Correct using the median of the overscan region",
48  "MEDIAN_PER_ROW": "Correct using the median per row of the overscan region",
49  },
50  )
51  order = pexConfig.Field(
52  dtype=int,
53  doc=("Order of polynomial to fit if overscan fit type is a polynomial, "
54  "or number of spline knots if overscan fit type is a spline."),
55  default=1,
56  )
57  numSigmaClip = pexConfig.Field(
58  dtype=float,
59  doc="Rejection threshold (sigma) for collapsing overscan before fit",
60  default=3.0,
61  )
62  maskPlanes = pexConfig.ListField(
63  dtype=str,
64  doc="Mask planes to reject when measuring overscan",
65  default=['SAT'],
66  )
67  overscanIsInt = pexConfig.Field(
68  dtype=bool,
69  doc="Treat overscan as an integer image for purposes of fitType=MEDIAN"
70  " and fitType=MEDIAN_PER_ROW.",
71  default=True,
72  )
73 
74 
75 class OverscanCorrectionTask(pipeBase.Task):
76  """Correction task for overscan.
77 
78  This class contains a number of utilities that are easier to
79  understand and use when they are not embedded in nested if/else
80  loops.
81 
82  Parameters
83  ----------
84  statControl : `lsst.afw.math.StatisticsControl`, optional
85  Statistics control object.
86  """
87  ConfigClass = OverscanCorrectionTaskConfig
88  _DefaultName = "overscan"
89 
90  def __init__(self, statControl=None, **kwargs):
91  super().__init__(**kwargs)
92  if statControl:
93  self.statControl = statControl
94  else:
95  self.statControl = afwMath.StatisticsControl()
96  self.statControl.setNumSigmaClip(self.config.numSigmaClip)
97  self.statControl.setAndMask(afwImage.Mask.getPlaneBitMask(self.config.maskPlanes))
98 
99  def run(self, ampImage, overscanImage):
100  """Measure and remove an overscan from an amplifier image.
101 
102  Parameters
103  ----------
104  ampImage : `lsst.afw.image.Image`
105  Image data that will have the overscan removed.
106  overscanImage : `lsst.afw.image.Image`
107  Overscan data that the overscan is measured from.
108 
109  Returns
110  -------
111  overscanResults : `lsst.pipe.base.Struct`
112  Result struct with components:
113 
114  ``imageFit``
115  Value or fit subtracted from the amplifier image data
116  (scalar or `lsst.afw.image.Image`).
117  ``overscanFit``
118  Value or fit subtracted from the overscan image data
119  (scalar or `lsst.afw.image.Image`).
120  ``overscanImage``
121  Image of the overscan region with the overscan
122  correction applied (`lsst.afw.image.Image`). This
123  quantity is used to estimate the amplifier read noise
124  empirically.
125 
126  Raises
127  ------
128  RuntimeError
129  Raised if an invalid overscan type is set.
130 
131  """
132  if self.config.fitType in ('MEAN', 'MEANCLIP', 'MEDIAN'):
133  overscanResult = self.measureConstantOverscan(overscanImage)
134  overscanValue = overscanResult.overscanValue
135  offImage = overscanValue
136  overscanModel = overscanValue
137  maskSuspect = None
138  elif self.config.fitType in ('MEDIAN_PER_ROW', 'POLY', 'CHEB', 'LEG',
139  'NATURAL_SPLINE', 'CUBIC_SPLINE', 'AKIMA_SPLINE'):
140  overscanResult = self.measureVectorOverscan(overscanImage)
141  overscanValue = overscanResult.overscanValue
142  maskArray = overscanResult.maskArray
143  isTransposed = overscanResult.isTransposed
144 
145  offImage = afwImage.ImageF(ampImage.getDimensions())
146  offArray = offImage.getArray()
147  overscanModel = afwImage.ImageF(overscanImage.getDimensions())
148  overscanArray = overscanModel.getArray()
149 
150  if hasattr(ampImage, 'getMask'):
151  maskSuspect = afwImage.Mask(ampImage.getDimensions())
152  else:
153  maskSuspect = None
154 
155  if isTransposed:
156  offArray[:, :] = overscanValue[np.newaxis, :]
157  overscanArray[:, :] = overscanValue[np.newaxis, :]
158  if maskSuspect:
159  maskSuspect.getArray()[:, maskArray] |= ampImage.getMask().getPlaneBitMask("SUSPECT")
160  else:
161  offArray[:, :] = overscanValue[:, np.newaxis]
162  overscanArray[:, :] = overscanValue[:, np.newaxis]
163  if maskSuspect:
164  maskSuspect.getArray()[maskArray, :] |= ampImage.getMask().getPlaneBitMask("SUSPECT")
165  else:
166  raise RuntimeError('%s : %s an invalid overscan type' %
167  ("overscanCorrection", self.config.fitType))
168 
169  self.debugView(overscanImage, overscanValue)
170 
171  ampImage -= offImage
172  if maskSuspect:
173  ampImage.getMask().getArray()[:, :] |= maskSuspect.getArray()[:, :]
174  overscanImage -= overscanModel
175  return pipeBase.Struct(imageFit=offImage,
176  overscanFit=overscanModel,
177  overscanImage=overscanImage,
178  edgeMask=maskSuspect)
179 
180  @staticmethod
181  def integerConvert(image):
182  """Return an integer version of the input image.
183 
184  Parameters
185  ----------
186  image : `numpy.ndarray`, `lsst.afw.image.Image` or `MaskedImage`
187  Image to convert to integers.
188 
189  Returns
190  -------
191  outI : `numpy.ndarray`, `lsst.afw.image.Image` or `MaskedImage`
192  The integer converted image.
193 
194  Raises
195  ------
196  RuntimeError
197  Raised if the input image could not be converted.
198  """
199  if hasattr(image, "image"):
200  # Is a maskedImage:
201  imageI = image.image.convertI()
202  outI = afwImage.MaskedImageI(imageI, image.mask, image.variance)
203  elif hasattr(image, "convertI"):
204  # Is an Image:
205  outI = image.convertI()
206  elif hasattr(image, "astype"):
207  # Is a numpy array:
208  outI = image.astype(int)
209  else:
210  raise RuntimeError("Could not convert this to integers: %s %s %s",
211  image, type(image), dir(image))
212  return outI
213 
214  # Constant methods
215  def measureConstantOverscan(self, image):
216  """Measure a constant overscan value.
217 
218  Parameters
219  ----------
220  image : `lsst.afw.image.Image` or `lsst.afw.image.MaskedImage`
221  Image data to measure the overscan from.
222 
223  Returns
224  -------
225  results : `lsst.pipe.base.Struct`
226  Overscan result with entries:
227  - ``overscanValue``: Overscan value to subtract (`float`)
228  - ``maskArray``: Placeholder for a mask array (`list`)
229  - ``isTransposed``: Orientation of the overscan (`bool`)
230  """
231  if self.config.fitType == 'MEDIAN':
232  calcImage = self.integerConvert(image)
233  else:
234  calcImage = image
235 
236  fitType = afwMath.stringToStatisticsProperty(self.config.fitType)
237  overscanValue = afwMath.makeStatistics(calcImage, fitType, self.statControl).getValue()
238 
239  return pipeBase.Struct(overscanValue=overscanValue,
240  maskArray=None,
241  isTransposed=False)
242 
243  # Vector correction utilities
244  def getImageArray(self, image):
245  """Extract the numpy array from the input image.
246 
247  Parameters
248  ----------
249  image : `lsst.afw.image.Image` or `lsst.afw.image.MaskedImage`
250  Image data to pull array from.
251 
252  calcImage : `numpy.ndarray`
253  Image data array for numpy operating.
254  """
255  if hasattr(image, "getImage"):
256  calcImage = image.getImage().getArray()
257  calcImage = np.ma.masked_where(image.getMask().getArray() & self.statControl.getAndMask(),
258  calcImage)
259  else:
260  calcImage = image.getArray()
261  return calcImage
262 
263  @staticmethod
264  def transpose(imageArray):
265  """Transpose input numpy array if necessary.
266 
267  Parameters
268  ----------
269  imageArray : `numpy.ndarray`
270  Image data to transpose.
271 
272  Returns
273  -------
274  imageArray : `numpy.ndarray`
275  Transposed image data.
276  isTransposed : `bool`
277  Indicates whether the input data was transposed.
278  """
279  if np.argmin(imageArray.shape) == 0:
280  return np.transpose(imageArray), True
281  else:
282  return imageArray, False
283 
284  def maskOutliers(self, imageArray):
285  """Mask outliers in a row of overscan data from a robust sigma
286  clipping procedure.
287 
288  Parameters
289  ----------
290  imageArray : `numpy.ndarray`
291  Image to filter along numpy axis=1.
292 
293  Returns
294  -------
295  maskedArray : `numpy.ma.masked_array`
296  Masked image marking outliers.
297  """
298  lq, median, uq = np.percentile(imageArray, [25.0, 50.0, 75.0], axis=1)
299  axisMedians = median
300  axisStdev = 0.74*(uq - lq) # robust stdev
301 
302  diff = np.abs(imageArray - axisMedians[:, np.newaxis])
303  return np.ma.masked_where(diff > self.statControl.getNumSigmaClip()
304  * axisStdev[:, np.newaxis], imageArray)
305 
306  @staticmethod
307  def collapseArray(maskedArray):
308  """Collapse overscan array (and mask) to a 1-D vector of values.
309 
310  Parameters
311  ----------
312  maskedArray : `numpy.ma.masked_array`
313  Masked array of input overscan data.
314 
315  Returns
316  -------
317  collapsed : `numpy.ma.masked_array`
318  Single dimensional overscan data, combined with the mean.
319  """
320  collapsed = np.mean(maskedArray, axis=1)
321  if collapsed.mask.sum() > 0:
322  collapsed.data[collapsed.mask] = np.mean(maskedArray.data[collapsed.mask], axis=1)
323  return collapsed
324 
325  def collapseArrayMedian(self, maskedArray):
326  """Collapse overscan array (and mask) to a 1-D vector of using the
327  correct integer median of row-values.
328 
329  Parameters
330  ----------
331  maskedArray : `numpy.ma.masked_array`
332  Masked array of input overscan data.
333 
334  Returns
335  -------
336  collapsed : `numpy.ma.masked_array`
337  Single dimensional overscan data, combined with the afwMath median.
338  """
339  integerMI = self.integerConvert(maskedArray)
340 
341  collapsed = []
342  fitType = afwMath.stringToStatisticsProperty('MEDIAN')
343  for row in integerMI:
344  newRow = row.compressed()
345  rowMedian = afwMath.makeStatistics(newRow, fitType, self.statControl).getValue()
346  collapsed.append(rowMedian)
347 
348  return np.array(collapsed)
349 
350  def splineFit(self, indices, collapsed, numBins):
351  """Wrapper function to match spline fit API to polynomial fit API.
352 
353  Parameters
354  ----------
355  indices : `numpy.ndarray`
356  Locations to evaluate the spline.
357  collapsed : `numpy.ndarray`
358  Collapsed overscan values corresponding to the spline
359  evaluation points.
360  numBins : `int`
361  Number of bins to use in constructing the spline.
362 
363  Returns
364  -------
365  interp : `lsst.afw.math.Interpolate`
366  Interpolation object for later evaluation.
367  """
368  if not np.ma.is_masked(collapsed):
369  collapsed.mask = np.array(len(collapsed)*[np.ma.nomask])
370 
371  numPerBin, binEdges = np.histogram(indices, bins=numBins,
372  weights=1 - collapsed.mask.astype(int))
373  with np.errstate(invalid="ignore"):
374  values = np.histogram(indices, bins=numBins,
375  weights=collapsed.data*~collapsed.mask)[0]/numPerBin
376  binCenters = np.histogram(indices, bins=numBins,
377  weights=indices*~collapsed.mask)[0]/numPerBin
378  interp = afwMath.makeInterpolate(binCenters.astype(float)[numPerBin > 0],
379  values.astype(float)[numPerBin > 0],
380  afwMath.stringToInterpStyle(self.config.fitType))
381  return interp
382 
383  @staticmethod
384  def splineEval(indices, interp):
385  """Wrapper function to match spline evaluation API to polynomial fit API.
386 
387  Parameters
388  ----------
389  indices : `numpy.ndarray`
390  Locations to evaluate the spline.
391  interp : `lsst.afw.math.interpolate`
392  Interpolation object to use.
393 
394  Returns
395  -------
396  values : `numpy.ndarray`
397  Evaluated spline values at each index.
398  """
399 
400  return interp.interpolate(indices.astype(float))
401 
402  @staticmethod
403  def maskExtrapolated(collapsed):
404  """Create mask if edges are extrapolated.
405 
406  Parameters
407  ----------
408  collapsed : `numpy.ma.masked_array`
409  Masked array to check the edges of.
410 
411  Returns
412  -------
413  maskArray : `numpy.ndarray`
414  Boolean numpy array of pixels to mask.
415  """
416  maskArray = np.full_like(collapsed, False, dtype=bool)
417  if np.ma.is_masked(collapsed):
418  num = len(collapsed)
419  for low in range(num):
420  if not collapsed.mask[low]:
421  break
422  if low > 0:
423  maskArray[:low] = True
424  for high in range(1, num):
425  if not collapsed.mask[-high]:
426  break
427  if high > 1:
428  maskArray[-high:] = True
429  return maskArray
430 
431  def measureVectorOverscan(self, image):
432  """Calculate the 1-d vector overscan from the input overscan image.
433 
434  Parameters
435  ----------
436  image : `lsst.afw.image.MaskedImage`
437  Image containing the overscan data.
438 
439  Returns
440  -------
441  results : `lsst.pipe.base.Struct`
442  Overscan result with entries:
443  - ``overscanValue``: Overscan value to subtract (`float`)
444  - ``maskArray`` : `list` [ `bool` ]
445  List of rows that should be masked as ``SUSPECT`` when the
446  overscan solution is applied.
447  - ``isTransposed`` : `bool`
448  Indicates if the overscan data was transposed during
449  calcuation, noting along which axis the overscan should be
450  subtracted.
451  """
452  calcImage = self.getImageArray(image)
453 
454  # operate on numpy-arrays from here
455  calcImage, isTransposed = self.transpose(calcImage)
456  masked = self.maskOutliers(calcImage)
457 
458  if self.config.fitType == 'MEDIAN_PER_ROW':
459  overscanVector = self.collapseArrayMedian(masked)
460  maskArray = self.maskExtrapolated(overscanVector)
461  else:
462  collapsed = self.collapseArray(masked)
463 
464  num = len(collapsed)
465  indices = 2.0*np.arange(num)/float(num) - 1.0
466 
467  poly = np.polynomial
468  fitter, evaler = {
469  'POLY': (poly.polynomial.polyfit, poly.polynomial.polyval),
470  'CHEB': (poly.chebyshev.chebfit, poly.chebyshev.chebval),
471  'LEG': (poly.legendre.legfit, poly.legendre.legval),
472  'NATURAL_SPLINE': (self.splineFit, self.splineEval),
473  'CUBIC_SPLINE': (self.splineFit, self.splineEval),
474  'AKIMA_SPLINE': (self.splineFit, self.splineEval)
475  }[self.config.fitType]
476 
477  coeffs = fitter(indices, collapsed, self.config.order)
478  overscanVector = evaler(indices, coeffs)
479  maskArray = self.maskExtrapolated(collapsed)
480  return pipeBase.Struct(overscanValue=np.array(overscanVector),
481  maskArray=maskArray,
482  isTransposed=isTransposed)
483 
484  def debugView(self, image, model):
485  """Debug display for the final overscan solution.
486 
487  Parameters
488  ----------
489  image : `lsst.afw.image.Image`
490  Input image the overscan solution was determined from.
491  model : `numpy.ndarray` or `float`
492  Overscan model determined for the image.
493  """
494  import lsstDebug
495  if not lsstDebug.Info(__name__).display:
496  return
497 
498  calcImage = self.getImageArray(image)
499  calcImage, isTransposed = self.transpose(calcImage)
500  masked = self.maskOutliers(calcImage)
501  collapsed = self.collapseArray(masked)
502 
503  num = len(collapsed)
504  indices = 2.0 * np.arange(num)/float(num) - 1.0
505 
506  if np.ma.is_masked(collapsed):
507  collapsedMask = collapsed.mask
508  else:
509  collapsedMask = np.array(num*[np.ma.nomask])
510 
511  import matplotlib.pyplot as plot
512  figure = plot.figure(1)
513  figure.clear()
514  axes = figure.add_axes((0.1, 0.1, 0.8, 0.8))
515  axes.plot(indices[~collapsedMask], collapsed[~collapsedMask], 'k+')
516  if collapsedMask.sum() > 0:
517  axes.plot(indices[collapsedMask], collapsed.data[collapsedMask], 'b+')
518  if isinstance(model, np.ndarray):
519  plotModel = model
520  else:
521  plotModel = np.zeros_like(indices)
522  plotModel += model
523  axes.plot(indices, plotModel, 'r-')
524  plot.xlabel("centered/scaled position along overscan region")
525  plot.ylabel("pixel value/fit value")
526  figure.show()
527  prompt = "Press Enter or c to continue [chp]..."
528  while True:
529  ans = input(prompt).lower()
530  if ans in ("", " ", "c",):
531  break
532  elif ans in ("p", ):
533  import pdb
534  pdb.set_trace()
535  elif ans in ("h", ):
536  print("[h]elp [c]ontinue [p]db")
537  plot.close()
lsst::ip::isr.overscan.OverscanCorrectionTask
Definition: overscan.py:75
lsst::afw::image
lsst::ip::isr.overscan.OverscanCorrectionTaskConfig
Definition: overscan.py:31
lsst::ip::isr.overscan.OverscanCorrectionTask.transpose
def transpose(imageArray)
Definition: overscan.py:264
lsst::ip::isr.overscan.OverscanCorrectionTask.collapseArray
def collapseArray(maskedArray)
Definition: overscan.py:307
lsst::ip::isr.overscan.OverscanCorrectionTask.integerConvert
def integerConvert(image)
Definition: overscan.py:181
lsst::ip::isr.overscan.OverscanCorrectionTask.maskExtrapolated
def maskExtrapolated(collapsed)
Definition: overscan.py:403
lsst::ip::isr.overscan.OverscanCorrectionTask.collapseArrayMedian
def collapseArrayMedian(self, maskedArray)
Definition: overscan.py:325
lsst::ip::isr.overscan.OverscanCorrectionTask.getImageArray
def getImageArray(self, image)
Definition: overscan.py:244
lsst::ip::isr.overscan.OverscanCorrectionTask.__init__
def __init__(self, statControl=None, **kwargs)
Definition: overscan.py:90
lsstDebug::Info
lsst::ip::isr.overscan.OverscanCorrectionTask.maskOutliers
def maskOutliers(self, imageArray)
Definition: overscan.py:284
lsst::ip::isr.overscan.OverscanCorrectionTask.statControl
statControl
Definition: overscan.py:93
lsst::ip::isr.overscan.OverscanCorrectionTask.measureVectorOverscan
def measureVectorOverscan(self, image)
Definition: overscan.py:431
lsst::pex::config
lsst::ip::isr.overscan.OverscanCorrectionTask.splineFit
def splineFit(self, indices, collapsed, numBins)
Definition: overscan.py:350
lsst::ip::isr.overscan.OverscanCorrectionTask.run
def run(self, ampImage, overscanImage)
Definition: overscan.py:99
lsst::afw::math
lsst::ip::isr.overscan.OverscanCorrectionTask.splineEval
def splineEval(indices, interp)
Definition: overscan.py:384
lsst::ip::isr.overscan.OverscanCorrectionTask.measureConstantOverscan
def measureConstantOverscan(self, image)
Definition: overscan.py:215
lsst::ip::isr.overscan.OverscanCorrectionTask.debugView
def debugView(self, image, model)
Definition: overscan.py:484
lsst::pipe::base