lsst.pipe.tasks gcb6d5b21a0+3d108a17e1
Loading...
Searching...
No Matches
computeExposureSummaryStats.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# (https://www.lsst.org).
6# See the COPYRIGHT file at the top-level directory of this distribution
7# for details of code ownership.
8#
9# This program is free software: you can redistribute it and/or modify
10# it under the terms of the GNU General Public License as published by
11# the Free Software Foundation, either version 3 of the License, or
12# (at your option) any later version.
13#
14# This program is distributed in the hope that it will be useful,
15# but WITHOUT ANY WARRANTY; without even the implied warranty of
16# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17# GNU General Public License for more details.
18#
19# You should have received a copy of the GNU General Public License
20# along with this program. If not, see <https://www.gnu.org/licenses/>.
21
22__all__ = ["ComputeExposureSummaryStatsTask", "ComputeExposureSummaryStatsConfig"]
23
24import warnings
25import numpy as np
26from scipy.stats import median_abs_deviation as sigmaMad
27import astropy.units as units
28from astropy.time import Time
29from astropy.coordinates import AltAz, SkyCoord, EarthLocation
30from lsst.daf.base import DateTime
31
32import lsst.pipe.base as pipeBase
33import lsst.pex.config as pexConfig
34import lsst.afw.math as afwMath
35import lsst.afw.image as afwImage
36import lsst.geom as geom
37from lsst.utils.timer import timeMethod
38
39
40class ComputeExposureSummaryStatsConfig(pexConfig.Config):
41 """Config for ComputeExposureSummaryTask"""
42 sigmaClip = pexConfig.Field(
43 dtype=float,
44 doc="Sigma for outlier rejection for sky noise.",
45 default=3.0,
46 )
47 clipIter = pexConfig.Field(
48 dtype=int,
49 doc="Number of iterations of outlier rejection for sky noise.",
50 default=2,
51 )
52 badMaskPlanes = pexConfig.ListField(
53 dtype=str,
54 doc="Mask planes that, if set, the associated pixel should not be included sky noise calculation.",
55 default=("NO_DATA", "SUSPECT"),
56 )
57 starSelection = pexConfig.Field(
58 doc="Field to select sources to be used in the PSF statistics computation.",
59 dtype=str,
60 default="calib_psf_used"
61 )
62 starShape = pexConfig.Field(
63 doc="Base name of columns to use for the source shape in the PSF statistics computation.",
64 dtype=str,
65 default="base_SdssShape"
66 )
67 psfShape = pexConfig.Field(
68 doc="Base name of columns to use for the PSF shape in the PSF statistics computation.",
69 dtype=str,
70 default="base_SdssShape_psf"
71 )
72 psfSampling = pexConfig.Field(
73 dtype=int,
74 doc="Sampling rate in pixels in each dimension for the maxDistToNearestPsf metric "
75 "caclulation grid (the tradeoff is between adequate sampling versus speed).",
76 default=8,
77 )
78 psfGridSampling = pexConfig.Field(
79 dtype=int,
80 doc="Sampling rate in pixels in each dimension for PSF model robustness metric "
81 "caclulations grid (the tradeoff is between adequate sampling versus speed).",
82 default=96,
83 )
84 psfBadMaskPlanes = pexConfig.ListField(
85 dtype=str,
86 doc="Mask planes that, if set, the associated pixel should not be included in the PSF model "
87 "robutsness metric calculations (namely, maxDistToNearestPsf and psfTraceRadiusDelta).",
88 default=("BAD", "CR", "EDGE", "INTRP", "NO_DATA", "SAT", "SUSPECT"),
89 )
90
91
93 """Task to compute exposure summary statistics.
94
95 This task computes various quantities suitable for DPDD and other
96 downstream processing at the detector centers, including:
97 - psfSigma
98 - psfArea
99 - psfIxx
100 - psfIyy
101 - psfIxy
102 - ra
103 - decl
104 - zenithDistance
105 - zeroPoint
106 - skyBg
107 - skyNoise
108 - meanVar
109 - raCorners
110 - decCorners
111 - astromOffsetMean
112 - astromOffsetStd
113
114 These additional quantities are computed from the stars in the detector:
115 - psfStarDeltaE1Median
116 - psfStarDeltaE2Median
117 - psfStarDeltaE1Scatter
118 - psfStarDeltaE2Scatter
119 - psfStarDeltaSizeMedian
120 - psfStarDeltaSizeScatter
121 - psfStarScaledDeltaSizeScatter
122
123 These quantities are computed based on the PSF model and image mask
124 to assess the robustness of the PSF model across a given detector
125 (against, e.g., extrapolation instability):
126 - maxDistToNearestPsf
127 - psfTraceRadiusDelta
128 """
129 ConfigClass = ComputeExposureSummaryStatsConfig
130 _DefaultName = "computeExposureSummaryStats"
131
132 @timeMethod
133 def run(self, exposure, sources, background):
134 """Measure exposure statistics from the exposure, sources, and
135 background.
136
137 Parameters
138 ----------
139 exposure : `lsst.afw.image.ExposureF`
141 background : `lsst.afw.math.BackgroundList`
142
143 Returns
144 -------
145 summary : `lsst.afw.image.ExposureSummary`
146 """
147 self.log.info("Measuring exposure statistics")
148
149 summary = afwImage.ExposureSummaryStats()
150
151 bbox = exposure.getBBox()
152
153 psf = exposure.getPsf()
154 self.update_psf_stats(summary, psf, bbox, sources, image_mask=exposure.mask)
155
156 wcs = exposure.getWcs()
157 visitInfo = exposure.getInfo().getVisitInfo()
158 self.update_wcs_stats(summary, wcs, bbox, visitInfo)
159
160 photoCalib = exposure.getPhotoCalib()
161 self.update_photo_calib_stats(summary, photoCalib)
162
163 self.update_background_stats(summary, background)
164
165 self.update_masked_image_stats(summary, exposure.getMaskedImage())
166
167 md = exposure.getMetadata()
168 if 'SFM_ASTROM_OFFSET_MEAN' in md:
169 summary.astromOffsetMean = md['SFM_ASTROM_OFFSET_MEAN']
170 summary.astromOffsetStd = md['SFM_ASTROM_OFFSET_STD']
171
172 return summary
173
174 def update_psf_stats(self, summary, psf, bbox, sources=None, image_mask=None, sources_columns=None):
175 """Compute all summary-statistic fields that depend on the PSF model.
176
177 Parameters
178 ----------
179 summary : `lsst.afw.image.ExposureSummaryStats`
180 Summary object to update in-place.
181 psf : `lsst.afw.detection.Psf` or `None`
182 Point spread function model. If `None`, all fields that depend on
183 the PSF will be reset (generally to NaN).
184 bbox : `lsst.geom.Box2I`
185 Bounding box of the image for which summary stats are being
186 computed.
187 sources : `lsst.afw.table.SourceCatalog`, optional
188 Catalog for quantities that are computed from source table columns.
189 If `None`, these quantities will be reset (generally to NaN).
190 image_mask : `lsst.afw.image.Mask`, optional
191 Mask image that may be used to compute distance-to-nearest-star
192 metrics.
193 sources_columns : `collections.abc.Set` [ `str` ], optional
194 Set of all column names in ``sources``. If provided, ``sources``
195 may be any table type for which string indexes yield column arrays.
196 If not provided, ``sources`` is assumed to be an
198 """
199 nan = float("nan")
200 summary.psfSigma = nan
201 summary.psfIxx = nan
202 summary.psfIyy = nan
203 summary.psfIxy = nan
204 summary.psfArea = nan
205 summary.nPsfStar = 0
206 summary.psfStarDeltaE1Median = nan
207 summary.psfStarDeltaE2Median = nan
208 summary.psfStarDeltaE1Scatter = nan
209 summary.psfStarDeltaE2Scatter = nan
210 summary.psfStarDeltaSizeMedian = nan
211 summary.psfStarDeltaSizeScatter = nan
212 summary.psfStarScaledDeltaSizeScatter = nan
213 summary.maxDistToNearestPsf = nan
214 summary.psfTraceRadiusDelta = nan
215
216 if psf is None:
217 return
218 shape = psf.computeShape(bbox.getCenter())
219 summary.psfSigma = shape.getDeterminantRadius()
220 summary.psfIxx = shape.getIxx()
221 summary.psfIyy = shape.getIyy()
222 summary.psfIxy = shape.getIxy()
223 im = psf.computeKernelImage(bbox.getCenter())
224 # The calculation of effective psf area is taken from
225 # meas_base/src/PsfFlux.cc#L112. See
226 # https://github.com/lsst/meas_base/blob/
227 # 750bffe6620e565bda731add1509507f5c40c8bb/src/PsfFlux.cc#L112
228 summary.psfArea = float(np.sum(im.array)/np.sum(im.array**2.))
229
230 if image_mask is not None:
231 psfTraceRadiusDelta = psf_trace_radius_delta(
232 image_mask,
233 psf,
234 sampling=self.config.psfGridSampling,
235 bad_mask_bits=self.config.psfBadMaskPlanes
236 )
237 summary.psfTraceRadiusDelta = float(psfTraceRadiusDelta)
238
239 if sources is None:
240 # No sources are available (as in some tests)
241 return
242
243 if sources_columns is None:
244 sources_columns = sources.schema.getNames()
245 if (
246 self.config.starSelection not in sources_columns
247 or self.config.starShape + '_flag' not in sources_columns
248 ):
249 # The source catalog does not have the necessary fields (as in some
250 # tests).
251 return
252
253 psf_mask = sources[self.config.starSelection] & (~sources[self.config.starShape + '_flag'])
254 nPsfStar = psf_mask.sum()
255
256 if nPsfStar == 0:
257 # No stars to measure statistics, so we must return the defaults
258 # of 0 stars and NaN values.
259 return
260 psf_cat = sources[psf_mask].copy(deep=True)
261
262 starXX = psf_cat[self.config.starShape + '_xx']
263 starYY = psf_cat[self.config.starShape + '_yy']
264 starXY = psf_cat[self.config.starShape + '_xy']
265 psfXX = psf_cat[self.config.psfShape + '_xx']
266 psfYY = psf_cat[self.config.psfShape + '_yy']
267 psfXY = psf_cat[self.config.psfShape + '_xy']
268
269 starSize = (starXX*starYY - starXY**2.)**0.25
270 starE1 = (starXX - starYY)/(starXX + starYY)
271 starE2 = 2*starXY/(starXX + starYY)
272 starSizeMedian = np.median(starSize)
273
274 psfSize = (psfXX*psfYY - psfXY**2)**0.25
275 psfE1 = (psfXX - psfYY)/(psfXX + psfYY)
276 psfE2 = 2*psfXY/(psfXX + psfYY)
277
278 psfStarDeltaE1Median = np.median(starE1 - psfE1)
279 psfStarDeltaE1Scatter = sigmaMad(starE1 - psfE1, scale='normal')
280 psfStarDeltaE2Median = np.median(starE2 - psfE2)
281 psfStarDeltaE2Scatter = sigmaMad(starE2 - psfE2, scale='normal')
282
283 psfStarDeltaSizeMedian = np.median(starSize - psfSize)
284 psfStarDeltaSizeScatter = sigmaMad(starSize - psfSize, scale='normal')
285 psfStarScaledDeltaSizeScatter = psfStarDeltaSizeScatter/starSizeMedian**2.
286
287 summary.nPsfStar = int(nPsfStar)
288 summary.psfStarDeltaE1Median = float(psfStarDeltaE1Median)
289 summary.psfStarDeltaE2Median = float(psfStarDeltaE2Median)
290 summary.psfStarDeltaE1Scatter = float(psfStarDeltaE1Scatter)
291 summary.psfStarDeltaE2Scatter = float(psfStarDeltaE2Scatter)
292 summary.psfStarDeltaSizeMedian = float(psfStarDeltaSizeMedian)
293 summary.psfStarDeltaSizeScatter = float(psfStarDeltaSizeScatter)
294 summary.psfStarScaledDeltaSizeScatter = float(psfStarScaledDeltaSizeScatter)
295
296 if image_mask is not None:
297 maxDistToNearestPsf = maximum_nearest_psf_distance(
298 image_mask,
299 psf_cat,
300 sampling=self.config.psfSampling,
301 bad_mask_bits=self.config.psfBadMaskPlanes
302 )
303 summary.maxDistToNearestPsf = float(maxDistToNearestPsf)
304
305 def update_wcs_stats(self, summary, wcs, bbox, visitInfo):
306 """Compute all summary-statistic fields that depend on the WCS model.
307
308 Parameters
309 ----------
310 summary : `lsst.afw.image.ExposureSummaryStats`
311 Summary object to update in-place.
312 wcs : `lsst.afw.geom.SkyWcs` or `None`
313 Astrometric calibration model. If `None`, all fields that depend
314 on the WCS will be reset (generally to NaN).
315 bbox : `lsst.geom.Box2I`
316 Bounding box of the image for which summary stats are being
317 computed.
318 visitInfo : `lsst.afw.image.VisitInfo`
319 Observation information used in together with ``wcs`` to compute
320 the zenith distance.
321 """
322 nan = float("nan")
323 summary.raCorners = [nan]*4
324 summary.decCorners = [nan]*4
325 summary.ra = nan
326 summary.decl = nan
327 summary.zenithDistance = nan
328
329 if wcs is None:
330 return
331
332 sph_pts = wcs.pixelToSky(geom.Box2D(bbox).getCorners())
333 summary.raCorners = [float(sph.getRa().asDegrees()) for sph in sph_pts]
334 summary.decCorners = [float(sph.getDec().asDegrees()) for sph in sph_pts]
335
336 sph_pt = wcs.pixelToSky(bbox.getCenter())
337 summary.ra = sph_pt.getRa().asDegrees()
338 summary.decl = sph_pt.getDec().asDegrees()
339
340 date = visitInfo.getDate()
341
342 if date.isValid():
343 # We compute the zenithDistance at the center of the detector
344 # rather than use the boresight value available via the visitInfo,
345 # because the zenithDistance may vary significantly over a large
346 # field of view.
347 observatory = visitInfo.getObservatory()
348 loc = EarthLocation(lat=observatory.getLatitude().asDegrees()*units.deg,
349 lon=observatory.getLongitude().asDegrees()*units.deg,
350 height=observatory.getElevation()*units.m)
351 obstime = Time(visitInfo.getDate().get(system=DateTime.MJD),
352 location=loc, format='mjd')
353 coord = SkyCoord(
354 summary.ra*units.degree,
355 summary.decl*units.degree,
356 obstime=obstime,
357 location=loc,
358 )
359 with warnings.catch_warnings():
360 warnings.simplefilter('ignore')
361 altaz = coord.transform_to(AltAz)
362
363 summary.zenithDistance = float(90.0 - altaz.alt.degree)
364
365 def update_photo_calib_stats(self, summary, photo_calib):
366 """Compute all summary-statistic fields that depend on the photometric
367 calibration model.
368
369 Parameters
370 ----------
371 summary : `lsst.afw.image.ExposureSummaryStats`
372 Summary object to update in-place.
373 photo_calib : `lsst.afw.image.PhotoCalib` or `None`
374 Photometric calibration model. If `None`, all fields that depend
375 on the photometric calibration will be reset (generally to NaN).
376 """
377 if photo_calib is not None:
378 summary.zeroPoint = float(2.5*np.log10(photo_calib.getInstFluxAtZeroMagnitude()))
379 else:
380 summary.zeroPoint = float("nan")
381
382 def update_background_stats(self, summary, background):
383 """Compute summary-statistic fields that depend only on the
384 background model.
385
386 Parameters
387 ----------
388 summary : `lsst.afw.image.ExposureSummaryStats`
389 Summary object to update in-place.
390 background : `lsst.afw.math.BackgroundList` or `None`
391 Background model. If `None`, all fields that depend on the
392 background will be reset (generally to NaN).
393
394 Notes
395 -----
396 This does not include fields that depend on the background-subtracted
397 masked image; when the background changes, it should generally be
398 applied to the image and `update_masked_image_stats` should be called
399 as well.
400 """
401 if background is not None:
402 bgStats = (bg[0].getStatsImage().getImage().array
403 for bg in background)
404 summary.skyBg = float(sum(np.median(bg[np.isfinite(bg)]) for bg in bgStats))
405 else:
406 summary.skyBg = float("nan")
407
408 def update_masked_image_stats(self, summary, masked_image):
409 """Compute summary-statistic fields that depend on the masked image
410 itself.
411
412 Parameters
413 ----------
414 summary : `lsst.afw.image.ExposureSummaryStats`
415 Summary object to update in-place.
416 masked_image : `lsst.afw.image.MaskedImage` or `None`
417 Masked image. If `None`, all fields that depend
418 on the masked image will be reset (generally to NaN).
419 """
420 nan = float("nan")
421 if masked_image is None:
422 summary.skyNoise = nan
423 summary.meanVar = nan
424 return
425 statsCtrl = afwMath.StatisticsControl()
426 statsCtrl.setNumSigmaClip(self.config.sigmaClip)
427 statsCtrl.setNumIter(self.config.clipIter)
428 statsCtrl.setAndMask(afwImage.Mask.getPlaneBitMask(self.config.badMaskPlanes))
429 statsCtrl.setNanSafe(True)
430
431 statObj = afwMath.makeStatistics(masked_image, afwMath.STDEVCLIP, statsCtrl)
432 skyNoise, _ = statObj.getResult(afwMath.STDEVCLIP)
433 summary.skyNoise = skyNoise
434
435 statObj = afwMath.makeStatistics(masked_image.variance, masked_image.mask, afwMath.MEANCLIP,
436 statsCtrl)
437 meanVar, _ = statObj.getResult(afwMath.MEANCLIP)
438 summary.meanVar = meanVar
439
440
442 image_mask,
443 psf_cat,
444 sampling=8,
445 bad_mask_bits=["BAD", "CR", "INTRP", "SAT", "SUSPECT", "NO_DATA", "EDGE"],
446):
447 """Compute the maximum distance of an unmasked pixel to its nearest PSF.
448
449 Parameters
450 ----------
451 image_mask : `lsst.afw.image.Mask`
452 The mask plane assosiated with the exposure.
454 Catalog containing only the stars used in the PSF modeling.
455 sampling : `int`
456 Sampling rate in each dimension to create the grid of points on which
457 to evaluate the distance to the nearest PSF star. The tradeoff is
458 between adequate sampling versus speed.
459 bad_mask_bits : `list` [`str`]
460 Mask bits required to be absent for a pixel to be considered
461 "unmasked".
462
463 Returns
464 -------
465 max_dist_to_nearest_psf : `float`
466 The maximum distance (in pixels) of an unmasked pixel to its nearest
467 PSF model star.
468 """
469 mask_arr = image_mask.array[::sampling, ::sampling]
470 bitmask = image_mask.getPlaneBitMask(bad_mask_bits)
471 good = ((mask_arr & bitmask) == 0)
472
473 x = np.arange(good.shape[1]) * sampling
474 y = np.arange(good.shape[0]) * sampling
475 xx, yy = np.meshgrid(x, y)
476
477 dist_to_nearest_psf = np.full(good.shape, np.inf)
478 for psf in psf_cat:
479 x_psf = psf.getX()
480 y_psf = psf.getY()
481 dist_to_nearest_psf = np.minimum(dist_to_nearest_psf, np.hypot(xx - x_psf, yy - y_psf))
482 unmasked_dists = dist_to_nearest_psf * good
483 max_dist_to_nearest_psf = np.max(unmasked_dists)
484
485 return max_dist_to_nearest_psf
486
487
489 image_mask,
490 image_psf,
491 sampling=96,
492 bad_mask_bits=["BAD", "CR", "INTRP", "SAT", "SUSPECT", "NO_DATA", "EDGE"],
493):
494 """Compute the delta between the maximum and minimum model PSF trace radius
495 values evaluated on a grid of points lying in the unmasked region of the
496 image.
497
498 Parameters
499 ----------
500 image_mask : `lsst.afw.image.Mask`
501 The mask plane assosiated with the exposure.
502 image_psf : `lsst.afw.detection.Psf`
503 The PSF model assosiated with the exposure.
504 sampling : `int`
505 Sampling rate in each dimension to create the grid of points at which
506 to evaluate ``image_psf``s trace radius value. The tradeoff is between
507 adequate sampling versus speed.
508 bad_mask_bits : `list` [`str`]
509 Mask bits required to be absent for a pixel to be considered
510 "unmasked".
511
512 Returns
513 -------
514 psf_trace_radius_delta : `float`
515 The delta (in pixels) between the maximum and minimum model PSF trace
516 radius values evaluated on the x,y-grid subsampled on the unmasked
517 detector pixels by a factor of ``sampling``. If any model PSF trace
518 radius value on the grid evaluates to NaN, then NaN is returned
519 immediately.
520 """
521 psf_trace_radius_list = []
522 mask_arr = image_mask.array[::sampling, ::sampling]
523 bitmask = image_mask.getPlaneBitMask(bad_mask_bits)
524 good = ((mask_arr & bitmask) == 0)
525
526 x = np.arange(good.shape[1]) * sampling
527 y = np.arange(good.shape[0]) * sampling
528 xx, yy = np.meshgrid(x, y)
529
530 for x_mesh, y_mesh, good_mesh in zip(xx, yy, good):
531 for x_point, y_point, is_good in zip(x_mesh, y_mesh, good_mesh):
532 if is_good:
533 psf_trace_radius = image_psf.computeShape(geom.Point2D(x_point, y_point)).getTraceRadius()
534 if ~np.isfinite(psf_trace_radius):
535 return float("nan")
536 psf_trace_radius_list.append(psf_trace_radius)
537
538 psf_trace_radius_delta = np.max(psf_trace_radius_list) - np.min(psf_trace_radius_list)
539
540 return psf_trace_radius_delta
def update_psf_stats(self, summary, psf, bbox, sources=None, image_mask=None, sources_columns=None)
def psf_trace_radius_delta(image_mask, image_psf, sampling=96, bad_mask_bits=["BAD", "CR", "INTRP", "SAT", "SUSPECT", "NO_DATA", "EDGE"])
def maximum_nearest_psf_distance(image_mask, psf_cat, sampling=8, bad_mask_bits=["BAD", "CR", "INTRP", "SAT", "SUSPECT", "NO_DATA", "EDGE"])