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