lsst.pipe.drivers  18.1.0-3-g2e29e3d+8
background.py
Go to the documentation of this file.
1 from __future__ import absolute_import, division, print_function
2 
3 import numpy
4 import itertools
5 from scipy.ndimage import gaussian_filter
6 
7 import lsst.afw.math as afwMath
8 import lsst.afw.image as afwImage
9 import lsst.afw.geom as afwGeom
10 import lsst.afw.cameraGeom as afwCameraGeom
11 import lsst.geom as geom
12 import lsst.meas.algorithms as measAlg
13 import lsst.afw.table as afwTable
14 
15 from lsst.pex.config import Config, Field, ListField, ChoiceField, ConfigField, RangeField, ConfigurableField
16 from lsst.pipe.base import Task
17 
18 
19 def robustMean(array, rej=3.0):
20  """Measure a robust mean of an array
21 
22  Parameters
23  ----------
24  array : `numpy.ndarray`
25  Array for which to measure the mean.
26  rej : `float`
27  k-sigma rejection threshold.
28 
29  Returns
30  -------
31  mean : `array.dtype`
32  Robust mean of `array`.
33  """
34  q1, median, q3 = numpy.percentile(array, [25.0, 50.0, 100.0])
35  good = numpy.abs(array - median) < rej*0.74*(q3 - q1)
36  return array[good].mean()
37 
38 
39 class BackgroundConfig(Config):
40  """Configuration for background measurement"""
41  statistic = ChoiceField(dtype=str, default="MEANCLIP", doc="type of statistic to use for grid points",
42  allowed={"MEANCLIP": "clipped mean",
43  "MEAN": "unclipped mean",
44  "MEDIAN": "median"})
45  xBinSize = RangeField(dtype=int, default=32, min=1, doc="Superpixel size in x")
46  yBinSize = RangeField(dtype=int, default=32, min=1, doc="Superpixel size in y")
47  algorithm = ChoiceField(dtype=str, default="NATURAL_SPLINE", optional=True,
48  doc="How to interpolate the background values. "
49  "This maps to an enum; see afw::math::Background",
50  allowed={
51  "CONSTANT": "Use a single constant value",
52  "LINEAR": "Use linear interpolation",
53  "NATURAL_SPLINE": "cubic spline with zero second derivative at endpoints",
54  "AKIMA_SPLINE": "higher-level nonlinear spline that is more robust"
55  " to outliers",
56  "NONE": "No background estimation is to be attempted",
57  })
58  mask = ListField(dtype=str, default=["SAT", "BAD", "EDGE", "DETECTED", "DETECTED_NEGATIVE", "NO_DATA"],
59  doc="Names of mask planes to ignore while estimating the background")
60 
61 
62 class SkyStatsConfig(Config):
63  """Parameters controlling the measurement of sky statistics"""
64  statistic = ChoiceField(dtype=str, default="MEANCLIP", doc="type of statistic to use for grid points",
65  allowed={"MEANCLIP": "clipped mean",
66  "MEAN": "unclipped mean",
67  "MEDIAN": "median"})
68  clip = Field(doc="Clipping threshold for background", dtype=float, default=3.0)
69  nIter = Field(doc="Clipping iterations for background", dtype=int, default=3)
70  mask = ListField(doc="Mask planes to reject", dtype=str,
71  default=["SAT", "DETECTED", "DETECTED_NEGATIVE", "BAD", "NO_DATA"])
72 
73 
74 class SkyMeasurementConfig(Config):
75  """Configuration for SkyMeasurementTask"""
76  skyIter = Field(dtype=int, default=3, doc="k-sigma rejection iterations for sky scale")
77  skyRej = Field(dtype=float, default=3.0, doc="k-sigma rejection threshold for sky scale")
78  background = ConfigField(dtype=BackgroundConfig, doc="Background measurement")
79  xNumSamples = Field(dtype=int, default=4, doc="Number of samples in x for scaling sky frame")
80  yNumSamples = Field(dtype=int, default=4, doc="Number of samples in y for scaling sky frame")
81  stats = ConfigField(dtype=SkyStatsConfig, doc="Measurement of sky statistics in the samples")
82 
83 
84 class SkyMeasurementTask(Task):
85  """Task for creating, persisting and using sky frames
86 
87  A sky frame is like a fringe frame (the sum of many exposures of the night sky,
88  combined with rejection to remove astrophysical objects) except the structure
89  is on larger scales, and hence we bin the images and represent them as a
90  background model (a `lsst.afw.math.BackgroundMI`). The sky frame represents
91  the dominant response of the camera to the sky background.
92  """
93  ConfigClass = SkyMeasurementConfig
94 
95  def getSkyData(self, butler, calibId):
96  """Retrieve sky frame from the butler
97 
98  Parameters
99  ----------
100  butler : `lsst.daf.persistence.Butler`
101  Data butler
102  calibId : `dict`
103  Data identifier for calib
104 
105  Returns
106  -------
107  sky : `lsst.afw.math.BackgroundList`
108  Sky frame
109  """
110  exp = butler.get("sky", calibId)
111  return self.exposureToBackground(exp)
112 
113  @staticmethod
115  """Convert an exposure to background model
116 
117  Calibs need to be persisted as an Exposure, so we need to convert
118  the persisted Exposure to a background model.
119 
120  Parameters
121  ----------
122  bgExp : `lsst.afw.image.Exposure`
123  Background model in Exposure format.
124 
125  Returns
126  -------
127  bg : `lsst.afw.math.BackgroundList`
128  Background model
129  """
130  header = bgExp.getMetadata()
131  xMin = header.getScalar("BOX.MINX")
132  yMin = header.getScalar("BOX.MINY")
133  xMax = header.getScalar("BOX.MAXX")
134  yMax = header.getScalar("BOX.MAXY")
135  algorithm = header.getScalar("ALGORITHM")
136  bbox = geom.Box2I(geom.Point2I(xMin, yMin), geom.Point2I(xMax, yMax))
137  return afwMath.BackgroundList(
138  (afwMath.BackgroundMI(bbox, bgExp.getMaskedImage()),
139  afwMath.stringToInterpStyle(algorithm),
140  afwMath.stringToUndersampleStyle("REDUCE_INTERP_ORDER"),
141  afwMath.ApproximateControl.UNKNOWN,
142  0, 0, False))
143 
144  def backgroundToExposure(self, statsImage, bbox):
145  """Convert a background model to an exposure
146 
147  Calibs need to be persisted as an Exposure, so we need to convert
148  the background model to an Exposure.
149 
150  Parameters
151  ----------
152  statsImage : `lsst.afw.image.MaskedImageF`
153  Background model's statistics image.
154  bbox : `lsst.geom.Box2I`
155  Bounding box for image.
156 
157  Returns
158  -------
159  exp : `lsst.afw.image.Exposure`
160  Background model in Exposure format.
161  """
162  exp = afwImage.makeExposure(statsImage)
163  header = exp.getMetadata()
164  header.set("BOX.MINX", bbox.getMinX())
165  header.set("BOX.MINY", bbox.getMinY())
166  header.set("BOX.MAXX", bbox.getMaxX())
167  header.set("BOX.MAXY", bbox.getMaxY())
168  header.set("ALGORITHM", self.config.background.algorithm)
169  return exp
170 
171  def measureBackground(self, image):
172  """Measure a background model for image
173 
174  This doesn't use a full-featured background model (e.g., no Chebyshev
175  approximation) because we just want the binning behaviour. This will
176  allow us to average the bins later (`averageBackgrounds`).
177 
178  The `BackgroundMI` is wrapped in a `BackgroundList` so it can be
179  pickled and persisted.
180 
181  Parameters
182  ----------
183  image : `lsst.afw.image.MaskedImage`
184  Image for which to measure background.
185 
186  Returns
187  -------
188  bgModel : `lsst.afw.math.BackgroundList`
189  Background model.
190  """
191  stats = afwMath.StatisticsControl()
192  stats.setAndMask(image.getMask().getPlaneBitMask(self.config.background.mask))
193  stats.setNanSafe(True)
194  ctrl = afwMath.BackgroundControl(
195  self.config.background.algorithm,
196  max(int(image.getWidth()/self.config.background.xBinSize + 0.5), 1),
197  max(int(image.getHeight()/self.config.background.yBinSize + 0.5), 1),
198  "REDUCE_INTERP_ORDER",
199  stats,
200  self.config.background.statistic
201  )
202 
203  bg = afwMath.makeBackground(image, ctrl)
204 
205  return afwMath.BackgroundList((
206  bg,
207  afwMath.stringToInterpStyle(self.config.background.algorithm),
208  afwMath.stringToUndersampleStyle("REDUCE_INTERP_ORDER"),
209  afwMath.ApproximateControl.UNKNOWN,
210  0, 0, False
211  ))
212 
213  def averageBackgrounds(self, bgList):
214  """Average multiple background models
215 
216  The input background models should be a `BackgroundList` consisting
217  of a single `BackgroundMI`.
218 
219  Parameters
220  ----------
221  bgList : `list` of `lsst.afw.math.BackgroundList`
222  Background models to average.
223 
224  Returns
225  -------
226  bgExp : `lsst.afw.image.Exposure`
227  Background model in Exposure format.
228  """
229  assert all(len(bg) == 1 for bg in bgList), "Mixed bgList: %s" % ([len(bg) for bg in bgList],)
230  images = [bg[0][0].getStatsImage() for bg in bgList]
231  boxes = [bg[0][0].getImageBBox() for bg in bgList]
232  assert len(set((box.getMinX(), box.getMinY(), box.getMaxX(), box.getMaxY()) for box in boxes)) == 1, \
233  "Bounding boxes not all equal"
234  bbox = boxes.pop(0)
235 
236  # Ensure bad pixels are masked
237  maskVal = afwImage.Mask.getPlaneBitMask("BAD")
238  for img in images:
239  bad = numpy.isnan(img.getImage().getArray())
240  img.getMask().getArray()[bad] = maskVal
241 
242  stats = afwMath.StatisticsControl()
243  stats.setAndMask(maskVal)
244  stats.setNanSafe(True)
245  combined = afwMath.statisticsStack(images, afwMath.MEANCLIP, stats)
246 
247  # Set bad pixels to the median
248  # Specifically NOT going to attempt to interpolate the bad values because we're only working on a
249  # single CCD here and can't use the entire field-of-view to do the interpolation (which is what we
250  # would need to avoid introducing problems at the edge of CCDs).
251  array = combined.getImage().getArray()
252  bad = numpy.isnan(array)
253  median = numpy.median(array[~bad])
254  array[bad] = median
255 
256  # Put it into an exposure, which is required for calibs
257  return self.backgroundToExposure(combined, bbox)
258 
259  def measureScale(self, image, skyBackground):
260  """Measure scale of background model in image
261 
262  We treat the sky frame much as we would a fringe frame
263  (except the length scale of the variations is different):
264  we measure samples on the input image and the sky frame,
265  which we will use to determine the scaling factor in the
266  'solveScales` method.
267 
268  Parameters
269  ----------
270  image : `lsst.afw.image.Exposure` or `lsst.afw.image.MaskedImage`
271  Science image for which to measure scale.
272  skyBackground : `lsst.afw.math.BackgroundList`
273  Sky background model.
274 
275  Returns
276  -------
277  imageSamples : `numpy.ndarray`
278  Sample measurements on image.
279  skySamples : `numpy.ndarray`
280  Sample measurements on sky frame.
281  """
282  if isinstance(image, afwImage.Exposure):
283  image = image.getMaskedImage()
284  # Ensure more samples than pixels
285  xNumSamples = min(self.config.xNumSamples, image.getWidth())
286  yNumSamples = min(self.config.yNumSamples, image.getHeight())
287  xLimits = numpy.linspace(0, image.getWidth(), xNumSamples + 1, dtype=int)
288  yLimits = numpy.linspace(0, image.getHeight(), yNumSamples + 1, dtype=int)
289  sky = skyBackground.getImage()
290  maskVal = image.getMask().getPlaneBitMask(self.config.stats.mask)
291  ctrl = afwMath.StatisticsControl(self.config.stats.clip, self.config.stats.nIter, maskVal)
292  statistic = afwMath.stringToStatisticsProperty(self.config.stats.statistic)
293  imageSamples = []
294  skySamples = []
295  for xIndex, yIndex in itertools.product(range(xNumSamples), range(yNumSamples)):
296  # -1 on the stop because Box2I is inclusive of the end point and we don't want to overlap boxes
297  xStart, xStop = xLimits[xIndex], xLimits[xIndex + 1] - 1
298  yStart, yStop = yLimits[yIndex], yLimits[yIndex + 1] - 1
299  box = geom.Box2I(geom.Point2I(xStart, yStart), geom.Point2I(xStop, yStop))
300  subImage = image.Factory(image, box)
301  subSky = sky.Factory(sky, box)
302  imageSamples.append(afwMath.makeStatistics(subImage, statistic, ctrl).getValue())
303  skySamples.append(afwMath.makeStatistics(subSky, statistic, ctrl).getValue())
304  return imageSamples, skySamples
305 
306  def solveScales(self, scales):
307  """Solve multiple scales for a single scale factor
308 
309  Having measured samples from the image and sky frame, we
310  fit for the scaling factor.
311 
312  Parameters
313  ----------
314  scales : `list` of a `tuple` of two `numpy.ndarray` arrays
315  A `list` of the results from `measureScale` method.
316 
317  Returns
318  -------
319  scale : `float`
320  Scale factor.
321  """
322  imageSamples = []
323  skySamples = []
324  for ii, ss in scales:
325  imageSamples.extend(ii)
326  skySamples.extend(ss)
327  assert len(imageSamples) == len(skySamples)
328  imageSamples = numpy.array(imageSamples)
329  skySamples = numpy.array(skySamples)
330 
331  def solve(mask):
332  return afwMath.LeastSquares.fromDesignMatrix(skySamples[mask].reshape(mask.sum(), 1),
333  imageSamples[mask],
334  afwMath.LeastSquares.DIRECT_SVD).getSolution()
335 
336  mask = numpy.isfinite(imageSamples) & numpy.isfinite(skySamples)
337  for ii in range(self.config.skyIter):
338  solution = solve(mask)
339  residuals = imageSamples - solution*skySamples
340  lq, uq = numpy.percentile(residuals[mask], [25, 75])
341  stdev = 0.741*(uq - lq) # Robust stdev from IQR
342  with numpy.errstate(invalid="ignore"): # suppress NAN warnings
343  bad = numpy.abs(residuals) > self.config.skyRej*stdev
344  mask[bad] = False
345 
346  return solve(mask)
347 
348  def subtractSkyFrame(self, image, skyBackground, scale, bgList=None):
349  """Subtract sky frame from science image
350 
351  Parameters
352  ----------
353  image : `lsst.afw.image.Exposure` or `lsst.afw.image.MaskedImage`
354  Science image.
355  skyBackground : `lsst.afw.math.BackgroundList`
356  Sky background model.
357  scale : `float`
358  Scale to apply to background model.
359  bgList : `lsst.afw.math.BackgroundList`
360  List of backgrounds applied to image
361  """
362  if isinstance(image, afwImage.Exposure):
363  image = image.getMaskedImage()
364  if isinstance(image, afwImage.MaskedImage):
365  image = image.getImage()
366  image.scaledMinus(scale, skyBackground.getImage())
367  if bgList is not None:
368  # Append the sky frame to the list of applied background models
369  bgData = list(skyBackground[0])
370  bg = bgData[0]
371  statsImage = bg.getStatsImage().clone()
372  statsImage *= scale
373  newBg = afwMath.BackgroundMI(bg.getImageBBox(), statsImage)
374  newBgData = [newBg] + bgData[1:]
375  bgList.append(newBgData)
376 
377 
378 def interpolate1D(method, xSample, ySample, xInterp):
379  """Interpolate in one dimension
380 
381  Interpolates the curve provided by `xSample` and `ySample` at
382  the positions of `xInterp`. Automatically backs off the
383  interpolation method to achieve successful interpolation.
384 
385  Parameters
386  ----------
387  method : `lsst.afw.math.Interpolate.Style`
388  Interpolation method to use.
389  xSample : `numpy.ndarray`
390  Vector of ordinates.
391  ySample : `numpy.ndarray`
392  Vector of coordinates.
393  xInterp : `numpy.ndarray`
394  Vector of ordinates to which to interpolate.
395 
396  Returns
397  -------
398  yInterp : `numpy.ndarray`
399  Vector of interpolated coordinates.
400 
401  """
402  if len(xSample) == 0:
403  return numpy.ones_like(xInterp)*numpy.nan
404  try:
405  return afwMath.makeInterpolate(xSample.astype(float), ySample.astype(float),
406  method).interpolate(xInterp.astype(float))
407  except Exception:
408  if method == afwMath.Interpolate.CONSTANT:
409  # We've already tried the most basic interpolation and it failed
410  return numpy.ones_like(xInterp)*numpy.nan
411  newMethod = afwMath.lookupMaxInterpStyle(len(xSample))
412  if newMethod == method:
413  newMethod = afwMath.Interpolate.CONSTANT
414  return interpolate1D(newMethod, xSample, ySample, xInterp)
415 
416 
417 def interpolateBadPixels(array, isBad, interpolationStyle):
418  """Interpolate bad pixels in an image array
419 
420  The bad pixels are modified in the array.
421 
422  Parameters
423  ----------
424  array : `numpy.ndarray`
425  Image array with bad pixels.
426  isBad : `numpy.ndarray` of type `bool`
427  Boolean array indicating which pixels are bad.
428  interpolationStyle : `str`
429  Style for interpolation (see `lsst.afw.math.Background`);
430  supported values are CONSTANT, LINEAR, NATURAL_SPLINE,
431  AKIMA_SPLINE.
432  """
433  if numpy.all(isBad):
434  raise RuntimeError("No good pixels in image array")
435  height, width = array.shape
436  xIndices = numpy.arange(width, dtype=float)
437  yIndices = numpy.arange(height, dtype=float)
438  method = afwMath.stringToInterpStyle(interpolationStyle)
439  isGood = ~isBad
440  for y in range(height):
441  if numpy.any(isBad[y, :]) and numpy.any(isGood[y, :]):
442  array[y][isBad[y]] = interpolate1D(method, xIndices[isGood[y]], array[y][isGood[y]],
443  xIndices[isBad[y]])
444 
445  isBad = numpy.isnan(array)
446  isGood = ~isBad
447  for x in range(width):
448  if numpy.any(isBad[:, x]) and numpy.any(isGood[:, x]):
449  array[:, x][isBad[:, x]] = interpolate1D(method, yIndices[isGood[:, x]],
450  array[:, x][isGood[:, x]], yIndices[isBad[:, x]])
451 
452 
454  """Configuration for FocalPlaneBackground
455 
456  Note that `xSize` and `ySize` are floating-point values, as
457  the focal plane frame is usually defined in units of microns
458  or millimetres rather than pixels. As such, their values will
459  need to be revised according to each particular camera. For
460  this reason, no defaults are set for those.
461  """
462  xSize = Field(dtype=float, doc="Bin size in x")
463  ySize = Field(dtype=float, doc="Bin size in y")
464  minFrac = Field(dtype=float, default=0.1, doc="Minimum fraction of bin size for good measurement")
465  mask = ListField(dtype=str, doc="Mask planes to treat as bad",
466  default=["BAD", "SAT", "INTRP", "DETECTED", "DETECTED_NEGATIVE", "EDGE", "NO_DATA"])
467  interpolation = ChoiceField(
468  doc="how to interpolate the background values. This maps to an enum; see afw::math::Background",
469  dtype=str, default="AKIMA_SPLINE", optional=True,
470  allowed={
471  "CONSTANT": "Use a single constant value",
472  "LINEAR": "Use linear interpolation",
473  "NATURAL_SPLINE": "cubic spline with zero second derivative at endpoints",
474  "AKIMA_SPLINE": "higher-level nonlinear spline that is more robust to outliers",
475  "NONE": "No background estimation is to be attempted",
476  },
477  )
478  doSmooth = Field(dtype=bool, default=False, doc="Do smoothing?")
479  smoothScale = Field(dtype=float, default=2.0, doc="Smoothing scale, as a multiple of the bin size")
480  binning = Field(dtype=int, default=64, doc="Binning to use for CCD background model (pixels)")
481 
482 
483 class FocalPlaneBackground(object):
484  """Background model for a focal plane camera
485 
486  We model the background empirically with the "superpixel" method: we
487  measure the background in each superpixel and interpolate between
488  superpixels to yield the model.
489 
490  The principal difference between this and `lsst.afw.math.BackgroundMI`
491  is that here the superpixels are defined in the frame of the focal
492  plane of the camera which removes discontinuities across detectors.
493 
494  The constructor you probably want to use is the `fromCamera` classmethod.
495 
496  There are two use patterns for building a background model:
497 
498  * Serial: create a `FocalPlaneBackground`, then `addCcd` for each of the
499  CCDs in an exposure.
500 
501  * Parallel: create a `FocalPlaneBackground`, then `clone` it for each
502  of the CCDs in an exposure and use those to `addCcd` their respective
503  CCD image. Finally, `merge` all the clones into the original.
504 
505  Once you've built the background model, you can apply it to individual
506  CCDs with the `toCcdBackground` method.
507  """
508  @classmethod
509  def fromCamera(cls, config, camera):
510  """Construct from a camera object
511 
512  Parameters
513  ----------
514  config : `FocalPlaneBackgroundConfig`
515  Configuration for measuring backgrounds.
516  camera : `lsst.afw.cameraGeom.Camera`
517  Camera for which to measure backgrounds.
518  """
519  cameraBox = geom.Box2D()
520  for ccd in camera:
521  for point in ccd.getCorners(afwCameraGeom.FOCAL_PLANE):
522  cameraBox.include(point)
523 
524  width, height = cameraBox.getDimensions()
525  # Offset so that we run from zero
526  offset = geom.Extent2D(cameraBox.getMin())*-1
527  # Add an extra pixel buffer on either side
528  dims = geom.Extent2I(int(numpy.ceil(width/config.xSize)) + 2,
529  int(numpy.ceil(height/config.ySize)) + 2)
530  # Transform takes us from focal plane coordinates --> sample coordinates
531  transform = (geom.AffineTransform.makeTranslation(geom.Extent2D(1, 1))*
532  geom.AffineTransform.makeScaling(1.0/config.xSize, 1.0/config.ySize)*
533  geom.AffineTransform.makeTranslation(offset))
534 
535  return cls(config, dims, afwGeom.makeTransform(transform))
536 
537  def __init__(self, config, dims, transform, values=None, numbers=None):
538  """Constructor
539 
540  Developers should note that changes to the signature of this method
541  require coordinated changes to the `__reduce__` and `clone` methods.
542 
543  Parameters
544  ----------
545  config : `FocalPlaneBackgroundConfig`
546  Configuration for measuring backgrounds.
547  dims : `lsst.geom.Extent2I`
548  Dimensions for background samples.
549  transform : `lsst.afw.geom.TransformPoint2ToPoint2`
550  Transformation from focal plane coordinates to sample coordinates.
551  values : `lsst.afw.image.ImageF`
552  Measured background values.
553  numbers : `lsst.afw.image.ImageF`
554  Number of pixels in each background measurement.
555  """
556  self.config = config
557  self.dims = dims
558  self.transform = transform
559 
560  if values is None:
561  values = afwImage.ImageF(self.dims)
562  values.set(0.0)
563  else:
564  values = values.clone()
565  assert(values.getDimensions() == self.dims)
566  self._values = values
567  if numbers is None:
568  numbers = afwImage.ImageF(self.dims) # float for dynamic range and convenience
569  numbers.set(0.0)
570  else:
571  numbers = numbers.clone()
572  assert(numbers.getDimensions() == self.dims)
573  self._numbers = numbers
574 
575  def __reduce__(self):
576  return self.__class__, (self.config, self.dims, self.transform, self._values, self._numbers)
577 
578  def clone(self):
579  return self.__class__(self.config, self.dims, self.transform, self._values, self._numbers)
580 
581  def addCcd(self, exposure):
582  """Add CCD to model
583 
584  We measure the background on the CCD (clipped mean), and record
585  the results in the model. For simplicity, measurements are made
586  in a box on the CCD corresponding to the warped coordinates of the
587  superpixel rather than accounting for little rotations, etc.
588  We also record the number of pixels used in the measurement so we
589  can have a measure of confidence in each bin's value.
590 
591  Parameters
592  ----------
593  exposure : `lsst.afw.image.Exposure`
594  CCD exposure to measure
595  """
596  detector = exposure.getDetector()
597  transform = detector.getTransformMap().getTransform(detector.makeCameraSys(afwCameraGeom.PIXELS),
598  detector.makeCameraSys(afwCameraGeom.FOCAL_PLANE))
599  image = exposure.getMaskedImage()
600  maskVal = image.getMask().getPlaneBitMask(self.config.mask)
601 
602  # Warp the binned image to the focal plane
603  toSample = transform.then(self.transform) # CCD pixels --> focal plane --> sample
604 
605  warped = afwImage.ImageF(self._values.getBBox())
606  warpedCounts = afwImage.ImageF(self._numbers.getBBox())
607  width, height = warped.getDimensions()
608 
609  stats = afwMath.StatisticsControl()
610  stats.setAndMask(maskVal)
611  stats.setNanSafe(True)
612  # Iterating over individual pixels in python is usually bad because it's slow, but there aren't many.
613  pixels = itertools.product(range(width), range(height))
614  for xx, yy in pixels:
615  llc = toSample.applyInverse(geom.Point2D(xx - 0.5, yy - 0.5))
616  urc = toSample.applyInverse(geom.Point2D(xx + 0.5, yy + 0.5))
617  bbox = geom.Box2I(geom.Point2I(llc), geom.Point2I(urc))
618  bbox.clip(image.getBBox())
619  if bbox.isEmpty():
620  continue
621  subImage = image.Factory(image, bbox)
622  result = afwMath.makeStatistics(subImage, afwMath.MEANCLIP | afwMath.NPOINT, stats)
623  mean = result.getValue(afwMath.MEANCLIP)
624  num = result.getValue(afwMath.NPOINT)
625  if not numpy.isfinite(mean) or not numpy.isfinite(num):
626  continue
627  warped[xx, yy, afwImage.LOCAL] = mean*num
628  warpedCounts[xx, yy, afwImage.LOCAL] = num
629 
630  self._values += warped
631  self._numbers += warpedCounts
632 
633  def toCcdBackground(self, detector, bbox):
634  """Produce a background model for a CCD
635 
636  The superpixel background model is warped back to the
637  CCD frame, for application to the individual CCD.
638 
639  Parameters
640  ----------
641  detector : `lsst.afw.cameraGeom.Detector`
642  CCD for which to produce background model.
643  bbox : `lsst.geom.Box2I`
644  Bounding box of CCD exposure.
645 
646  Returns
647  -------
648  bg : `lsst.afw.math.BackgroundList`
649  Background model for CCD.
650  """
651  transform = detector.getTransformMap().getTransform(detector.makeCameraSys(afwCameraGeom.PIXELS),
652  detector.makeCameraSys(afwCameraGeom.FOCAL_PLANE))
653  binTransform = (geom.AffineTransform.makeScaling(self.config.binning)*
654  geom.AffineTransform.makeTranslation(geom.Extent2D(0.5, 0.5)))
655 
656  # Binned image on CCD --> unbinned image on CCD --> focal plane --> binned focal plane
657  toSample = afwGeom.makeTransform(binTransform).then(transform).then(self.transform)
658 
659  focalPlane = self.getStatsImage()
660  fpNorm = afwImage.ImageF(focalPlane.getBBox())
661  fpNorm.set(1.0)
662 
663  image = afwImage.ImageF(bbox.getDimensions()//self.config.binning)
664  norm = afwImage.ImageF(image.getBBox())
665  ctrl = afwMath.WarpingControl("bilinear")
666  afwMath.warpImage(image, focalPlane, toSample.inverted(), ctrl)
667  afwMath.warpImage(norm, fpNorm, toSample.inverted(), ctrl)
668  image /= norm
669 
670  mask = afwImage.Mask(image.getBBox())
671  isBad = numpy.isnan(image.getArray())
672  mask.getArray()[isBad] = mask.getPlaneBitMask("BAD")
673  image.getArray()[isBad] = image.getArray()[~isBad].mean()
674 
675  return afwMath.BackgroundList(
676  (afwMath.BackgroundMI(bbox, afwImage.makeMaskedImage(image, mask)),
677  afwMath.stringToInterpStyle(self.config.interpolation),
678  afwMath.stringToUndersampleStyle("REDUCE_INTERP_ORDER"),
679  afwMath.ApproximateControl.UNKNOWN,
680  0, 0, False)
681  )
682 
683  def merge(self, other):
684  """Merge with another FocalPlaneBackground
685 
686  This allows multiple background models to be constructed from
687  different CCDs, and then merged to form a single consistent
688  background model for the entire focal plane.
689 
690  Parameters
691  ----------
692  other : `FocalPlaneBackground`
693  Another background model to merge.
694 
695  Returns
696  -------
697  self : `FocalPlaneBackground`
698  The merged background model.
699  """
700  if (self.config.xSize, self.config.ySize) != (other.config.xSize, other.config.ySize):
701  raise RuntimeError("Size mismatch: %s vs %s" % ((self.config.xSize, self.config.ySize),
702  (other.config.xSize, other.config.ySize)))
703  if self.dims != other.dims:
704  raise RuntimeError("Dimensions mismatch: %s vs %s" % (self.dims, other.dims))
705  self._values += other._values
706  self._numbers += other._numbers
707  return self
708 
709  def __iadd__(self, other):
710  """Merge with another FocalPlaneBackground
711 
712  Parameters
713  ----------
714  other : `FocalPlaneBackground`
715  Another background model to merge.
716 
717  Returns
718  -------
719  self : `FocalPlaneBackground`
720  The merged background model.
721  """
722  return self.merge(other)
723 
724  def getStatsImage(self):
725  """Return the background model data
726 
727  This is the measurement of the background for each of the superpixels.
728  """
729  values = self._values.clone()
730  values /= self._numbers
731  thresh = self.config.minFrac*self.config.xSize*self.config.ySize
732  isBad = self._numbers.getArray() < thresh
733  if self.config.doSmooth:
734  array = values.getArray()
735  array[:] = smoothArray(array, isBad, self.config.smoothScale)
736  isBad = numpy.isnan(values.array)
737  if numpy.any(isBad):
738  interpolateBadPixels(values.getArray(), isBad, self.config.interpolation)
739  return values
740 
741 
742 class MaskObjectsConfig(Config):
743  """Configuration for MaskObjectsTask"""
744  nIter = Field(dtype=int, default=3, doc="Number of iterations")
745  subtractBackground = ConfigurableField(target=measAlg.SubtractBackgroundTask,
746  doc="Background subtraction")
747  detection = ConfigurableField(target=measAlg.SourceDetectionTask, doc="Source detection")
748  detectSigma = Field(dtype=float, default=5.0, doc="Detection threshold (standard deviations)")
749  doInterpolate = Field(dtype=bool, default=True, doc="Interpolate when removing objects?")
750  interpolate = ConfigurableField(target=measAlg.SubtractBackgroundTask, doc="Interpolation")
751 
752  def setDefaults(self):
753  self.detection.reEstimateBackground = False
754  self.detection.doTempLocalBackground = False
755  self.detection.doTempWideBackground = False
756  self.detection.thresholdValue = 2.5
757  self.subtractBackground.binSize = 1024
758  self.subtractBackground.useApprox = False
759  self.interpolate.binSize = 256
760  self.interpolate.useApprox = False
761 
762  def validate(self):
763  if (self.detection.reEstimateBackground or
764  self.detection.doTempLocalBackground or
765  self.detection.doTempWideBackground):
766  raise RuntimeError("Incorrect settings for object masking: reEstimateBackground, "
767  "doTempLocalBackground and doTempWideBackground must be False")
768 
769 
770 class MaskObjectsTask(Task):
771  """Iterative masking of objects on an Exposure
772 
773  This task makes more exhaustive object mask by iteratively doing detection
774  and background-subtraction. The purpose of this task is to get true
775  background removing faint tails of large objects. This is useful to get a
776  clean sky estimate from relatively small number of visits.
777 
778  We deliberately use the specified ``detectSigma`` instead of the PSF,
779  in order to better pick up the faint wings of objects.
780  """
781  ConfigClass = MaskObjectsConfig
782 
783  def __init__(self, *args, **kwargs):
784  super().__init__(*args, **kwargs)
785  # Disposable schema suppresses warning from SourceDetectionTask.__init__
786  self.makeSubtask("detection", schema=afwTable.Schema())
787  self.makeSubtask("interpolate")
788  self.makeSubtask("subtractBackground")
789 
790  def run(self, exposure, maskPlanes=None):
791  """Mask objects on Exposure
792 
793  Objects are found and removed.
794 
795  Parameters
796  ----------
797  exposure : `lsst.afw.image.Exposure`
798  Exposure on which to mask objects.
799  maskPlanes : iterable of `str`, optional
800  List of mask planes to remove.
801  """
802  self.findObjects(exposure)
803  self.removeObjects(exposure, maskPlanes)
804 
805  def findObjects(self, exposure):
806  """Iteratively find objects on an exposure
807 
808  Objects are masked with the ``DETECTED`` mask plane.
809 
810  Parameters
811  ----------
812  exposure : `lsst.afw.image.Exposure`
813  Exposure on which to mask objects.
814  """
815  for _ in range(self.config.nIter):
816  bg = self.subtractBackground.run(exposure).background
817  self.detection.detectFootprints(exposure, sigma=self.config.detectSigma, clearMask=True)
818  exposure.maskedImage += bg.getImage()
819 
820  def removeObjects(self, exposure, maskPlanes=None):
821  """Remove objects from exposure
822 
823  We interpolate over using a background model if ``doInterpolate`` is
824  set; otherwise we simply replace everything with the median.
825 
826  Parameters
827  ----------
828  exposure : `lsst.afw.image.Exposure`
829  Exposure on which to mask objects.
830  maskPlanes : iterable of `str`, optional
831  List of mask planes to remove. ``DETECTED`` will be added as well.
832  """
833  image = exposure.image
834  mask = exposure.mask
835  maskVal = mask.getPlaneBitMask("DETECTED")
836  if maskPlanes is not None:
837  maskVal |= mask.getPlaneBitMask(maskPlanes)
838  isBad = mask.array & maskVal > 0
839 
840  if self.config.doInterpolate:
841  smooth = self.interpolate.fitBackground(exposure.maskedImage)
842  replace = smooth.getImageF().array[isBad]
843  mask.array &= ~mask.getPlaneBitMask(["DETECTED"])
844  else:
845  replace = numpy.median(image.array[~isBad])
846  image.array[isBad] = replace
847 
848 
849 def smoothArray(array, bad, sigma):
850  """Gaussian-smooth an array while ignoring bad pixels
851 
852  It's not sufficient to set the bad pixels to zero, as then they're treated
853  as if they are zero, rather than being ignored altogether. We need to apply
854  a correction to that image that removes the effect of the bad pixels.
855 
856  Parameters
857  ----------
858  array : `numpy.ndarray` of floating-point
859  Array to smooth.
860  bad : `numpy.ndarray` of `bool`
861  Flag array indicating bad pixels.
862  sigma : `float`
863  Gaussian sigma.
864 
865  Returns
866  -------
867  convolved : `numpy.ndarray`
868  Smoothed image.
869  """
870  convolved = gaussian_filter(numpy.where(bad, 0.0, array), sigma, mode="constant", cval=0.0)
871  denominator = gaussian_filter(numpy.where(bad, 0.0, 1.0), sigma, mode="constant", cval=0.0)
872  return convolved/denominator
def removeObjects(self, exposure, maskPlanes=None)
Definition: background.py:820
def robustMean(array, rej=3.0)
Definition: background.py:19
def interpolateBadPixels(array, isBad, interpolationStyle)
Definition: background.py:417
def subtractSkyFrame(self, image, skyBackground, scale, bgList=None)
Definition: background.py:348
def run(self, exposure, maskPlanes=None)
Definition: background.py:790
def smoothArray(array, bad, sigma)
Definition: background.py:849
def interpolate1D(method, xSample, ySample, xInterp)
Definition: background.py:378
def __init__(self, config, dims, transform, values=None, numbers=None)
Definition: background.py:537
def measureScale(self, image, skyBackground)
Definition: background.py:259
def backgroundToExposure(self, statsImage, bbox)
Definition: background.py:144