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