lsst.pipe.tasks  v23.0.x-g0162f352a6+f2fbba1123
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 import warnings
22 import numpy as np
23 from scipy.stats import median_abs_deviation as sigmaMad
24 import astropy.units as units
25 from astropy.time import Time
26 from astropy.coordinates import AltAz, SkyCoord, EarthLocation
27 from lsst.daf.base import DateTime
28 
29 import lsst.pipe.base as pipeBase
30 import lsst.pex.config as pexConfig
31 import lsst.afw.math as afwMath
32 import lsst.afw.image as afwImage
33 import lsst.geom
34 
35 
36 __all__ = ("ComputeExposureSummaryStatsTask", "ComputeExposureSummaryStatsConfig")
37 
38 
39 class ComputeExposureSummaryStatsConfig(pexConfig.Config):
40  """Config for ComputeExposureSummaryTask"""
41  sigmaClip = pexConfig.Field(
42  dtype=float,
43  doc="Sigma for outlier rejection for sky noise.",
44  default=3.0,
45  )
46  clipIter = pexConfig.Field(
47  dtype=int,
48  doc="Number of iterations of outlier rejection for sky noise.",
49  default=2,
50  )
51  badMaskPlanes = pexConfig.ListField(
52  dtype=str,
53  doc="Mask planes that, if set, the associated pixel should not be included sky noise calculation.",
54  default=("NO_DATA", "SUSPECT"),
55  )
56  starSelection = pexConfig.Field(
57  doc="Field to select sources to be used in the PSF statistics computation.",
58  dtype=str,
59  default="calib_psf_used"
60  )
61  starShape = pexConfig.Field(
62  doc="Base name of columns to use for the source shape in the PSF statistics computation.",
63  dtype=str,
64  default="base_SdssShape"
65  )
66  psfShape = pexConfig.Field(
67  doc="Base name of columns to use for the PSF shape in the PSF statistics computation.",
68  dtype=str,
69  default="base_SdssShape_psf"
70  )
71 
72 
73 class ComputeExposureSummaryStatsTask(pipeBase.Task):
74  """Task to compute exposure summary statistics.
75 
76  This task computes various quantities suitable for DPDD and other
77  downstream processing at the detector centers, including:
78  - psfSigma
79  - psfArea
80  - psfIxx
81  - psfIyy
82  - psfIxy
83  - ra
84  - decl
85  - zenithDistance
86  - zeroPoint
87  - skyBg
88  - skyNoise
89  - meanVar
90  - raCorners
91  - decCorners
92  - astromOffsetMean
93  - astromOffsetStd
94 
95  These additional quantities are computed from the stars in the detector:
96  - psfStarDeltaE1Median
97  - psfStarDeltaE2Median
98  - psfStarDeltaE1Scatter
99  - psfStarDeltaE2Scatter
100  - psfStarDeltaSizeMedian
101  - psfStarDeltaSizeScatter
102  - psfStarScaledDeltaSizeScatter
103  """
104  ConfigClass = ComputeExposureSummaryStatsConfig
105  _DefaultName = "computeExposureSummaryStats"
106 
107  @pipeBase.timeMethod
108  def run(self, exposure, sources, background):
109  """Measure exposure statistics from the exposure, sources, and background.
110 
111  Parameters
112  ----------
113  exposure : `lsst.afw.image.ExposureF`
114  sources : `lsst.afw.table.SourceCatalog`
115  background : `lsst.afw.math.BackgroundList`
116 
117  Returns
118  -------
119  summary : `lsst.afw.image.ExposureSummary`
120  """
121  self.log.info("Measuring exposure statistics")
122 
123  bbox = exposure.getBBox()
124 
125  psf = exposure.getPsf()
126  if psf is not None:
127  shape = psf.computeShape(bbox.getCenter())
128  psfSigma = shape.getDeterminantRadius()
129  psfIxx = shape.getIxx()
130  psfIyy = shape.getIyy()
131  psfIxy = shape.getIxy()
132  im = psf.computeKernelImage(bbox.getCenter())
133  # The calculation of effective psf area is taken from
134  # meas_base/src/PsfFlux.cc#L112. See
135  # https://github.com/lsst/meas_base/blob/
136  # 750bffe6620e565bda731add1509507f5c40c8bb/src/PsfFlux.cc#L112
137  psfArea = np.sum(im.array)/np.sum(im.array**2.)
138  else:
139  psfSigma = np.nan
140  psfIxx = np.nan
141  psfIyy = np.nan
142  psfIxy = np.nan
143  psfArea = np.nan
144 
145  wcs = exposure.getWcs()
146  if wcs is not None:
147  sph_pts = wcs.pixelToSky(lsst.geom.Box2D(bbox).getCorners())
148  raCorners = [float(sph.getRa().asDegrees()) for sph in sph_pts]
149  decCorners = [float(sph.getDec().asDegrees()) for sph in sph_pts]
150 
151  sph_pt = wcs.pixelToSky(bbox.getCenter())
152  ra = sph_pt.getRa().asDegrees()
153  decl = sph_pt.getDec().asDegrees()
154  else:
155  raCorners = [float(np.nan)]*4
156  decCorners = [float(np.nan)]*4
157  ra = np.nan
158  decl = np.nan
159 
160  photoCalib = exposure.getPhotoCalib()
161  if photoCalib is not None:
162  zeroPoint = 2.5*np.log10(photoCalib.getInstFluxAtZeroMagnitude())
163  else:
164  zeroPoint = np.nan
165 
166  visitInfo = exposure.getInfo().getVisitInfo()
167  date = visitInfo.getDate()
168 
169  if date.isValid():
170  # We compute the zenithDistance at the center of the detector rather
171  # than use the boresight value available via the visitInfo, because
172  # the zenithDistance may vary significantly over a large field of view.
173  observatory = visitInfo.getObservatory()
174  loc = EarthLocation(lat=observatory.getLatitude().asDegrees()*units.deg,
175  lon=observatory.getLongitude().asDegrees()*units.deg,
176  height=observatory.getElevation()*units.m)
177  obstime = Time(visitInfo.getDate().get(system=DateTime.MJD),
178  location=loc, format='mjd')
179  coord = SkyCoord(ra*units.degree, decl*units.degree, obstime=obstime, location=loc)
180  with warnings.catch_warnings():
181  warnings.simplefilter('ignore')
182  altaz = coord.transform_to(AltAz)
183 
184  zenithDistance = 90.0 - altaz.alt.degree
185  else:
186  zenithDistance = np.nan
187 
188  if background is not None:
189  bgStats = (bg[0].getStatsImage().getImage().array
190  for bg in background)
191  skyBg = sum(np.median(bg[np.isfinite(bg)]) for bg in bgStats)
192  else:
193  skyBg = np.nan
194 
195  statsCtrl = afwMath.StatisticsControl()
196  statsCtrl.setNumSigmaClip(self.config.sigmaClip)
197  statsCtrl.setNumIter(self.config.clipIter)
198  statsCtrl.setAndMask(afwImage.Mask.getPlaneBitMask(self.config.badMaskPlanes))
199  statsCtrl.setNanSafe(True)
200 
201  statObj = afwMath.makeStatistics(exposure.getMaskedImage(), afwMath.STDEVCLIP,
202  statsCtrl)
203  skyNoise, _ = statObj.getResult(afwMath.STDEVCLIP)
204 
205  statObj = afwMath.makeStatistics(exposure.getMaskedImage().getVariance(),
206  exposure.getMaskedImage().getMask(),
207  afwMath.MEANCLIP, statsCtrl)
208  meanVar, _ = statObj.getResult(afwMath.MEANCLIP)
209 
210  md = exposure.getMetadata()
211  if 'SFM_ASTROM_OFFSET_MEAN' in md:
212  astromOffsetMean = md['SFM_ASTROM_OFFSET_MEAN']
213  astromOffsetStd = md['SFM_ASTROM_OFFSET_STD']
214  else:
215  astromOffsetMean = np.nan
216  astromOffsetStd = np.nan
217 
218  psfStats = self._computePsfStats_computePsfStats(sources)
219 
220  # Note that all numpy values have to be explicitly converted to
221  # python floats for yaml serialization.
222  summary = afwImage.ExposureSummaryStats(
223  psfSigma=float(psfSigma),
224  psfArea=float(psfArea),
225  psfIxx=float(psfIxx),
226  psfIyy=float(psfIyy),
227  psfIxy=float(psfIxy),
228  ra=float(ra),
229  decl=float(decl),
230  zenithDistance=float(zenithDistance),
231  zeroPoint=float(zeroPoint),
232  skyBg=float(skyBg),
233  skyNoise=float(skyNoise),
234  meanVar=float(meanVar),
235  raCorners=raCorners,
236  decCorners=decCorners,
237  astromOffsetMean=astromOffsetMean,
238  astromOffsetStd=astromOffsetStd,
239  nPsfStar=psfStats.nPsfStar,
240  psfStarDeltaE1Median=psfStats.psfStarDeltaE1Median,
241  psfStarDeltaE2Median=psfStats.psfStarDeltaE2Median,
242  psfStarDeltaE1Scatter=psfStats.psfStarDeltaE1Scatter,
243  psfStarDeltaE2Scatter=psfStats.psfStarDeltaE2Scatter,
244  psfStarDeltaSizeMedian=psfStats.psfStarDeltaSizeMedian,
245  psfStarDeltaSizeScatter=psfStats.psfStarDeltaSizeScatter,
246  psfStarScaledDeltaSizeScatter=psfStats.psfStarScaledDeltaSizeScatter
247  )
248 
249  return summary
250 
251  def _computePsfStats(self, sources):
252  """Compute psf residual statistics.
253 
254  All residuals are computed using median statistics on the difference
255  between the sources and the models.
256 
257  Parameters
258  ----------
259  sources : `lsst.afw.table.SourceCatalog`
260  Source catalog on which to measure the PSF statistics.
261 
262  Returns
263  -------
264  psfStats : `lsst.pipe.base.Struct`
265  Struct with various psf stats.
266  """
267  psfStats = pipeBase.Struct(nPsfStar=0,
268  psfStarDeltaE1Median=float(np.nan),
269  psfStarDeltaE2Median=float(np.nan),
270  psfStarDeltaE1Scatter=float(np.nan),
271  psfStarDeltaE2Scatter=float(np.nan),
272  psfStarDeltaSizeMedian=float(np.nan),
273  psfStarDeltaSizeScatter=float(np.nan),
274  psfStarScaledDeltaSizeScatter=float(np.nan))
275 
276  if sources is None:
277  # No sources are available (as in some tests)
278  return psfStats
279 
280  names = sources.schema.getNames()
281  if self.config.starSelection not in names or self.config.starShape + '_flag' not in names:
282  # The source catalog does not have the necessary fields (as in some tests)
283  return psfStats
284 
285  mask = sources[self.config.starSelection] & (~sources[self.config.starShape + '_flag'])
286  nPsfStar = mask.sum()
287 
288  if nPsfStar == 0:
289  # No stars to measure statistics, so we must return the defaults
290  # of 0 stars and NaN values.
291  return psfStats
292 
293  starXX = sources[self.config.starShape + '_xx'][mask]
294  starYY = sources[self.config.starShape + '_yy'][mask]
295  starXY = sources[self.config.starShape + '_xy'][mask]
296  psfXX = sources[self.config.psfShape + '_xx'][mask]
297  psfYY = sources[self.config.psfShape + '_yy'][mask]
298  psfXY = sources[self.config.psfShape + '_xy'][mask]
299 
300  starSize = (starXX*starYY - starXY**2.)**0.25
301  starE1 = (starXX - starYY)/(starXX + starYY)
302  starE2 = 2*starXY/(starXX + starYY)
303  starSizeMedian = np.median(starSize)
304 
305  psfSize = (psfXX*psfYY - psfXY**2)**0.25
306  psfE1 = (psfXX - psfYY)/(psfXX + psfYY)
307  psfE2 = 2*psfXY/(psfXX + psfYY)
308 
309  psfStarDeltaE1Median = np.median(starE1 - psfE1)
310  psfStarDeltaE1Scatter = sigmaMad(starE1 - psfE1, scale='normal')
311  psfStarDeltaE2Median = np.median(starE2 - psfE2)
312  psfStarDeltaE2Scatter = sigmaMad(starE2 - psfE2, scale='normal')
313 
314  psfStarDeltaSizeMedian = np.median(starSize - psfSize)
315  psfStarDeltaSizeScatter = sigmaMad(starSize - psfSize, scale='normal')
316  psfStarScaledDeltaSizeScatter = psfStarDeltaSizeScatter/starSizeMedian**2.
317 
318  psfStats.nPsfStar = int(nPsfStar)
319  psfStats.psfStarDeltaE1Median = float(psfStarDeltaE1Median)
320  psfStats.psfStarDeltaE2Median = float(psfStarDeltaE2Median)
321  psfStats.psfStarDeltaE1Scatter = float(psfStarDeltaE1Scatter)
322  psfStats.psfStarDeltaE2Scatter = float(psfStarDeltaE2Scatter)
323  psfStats.psfStarDeltaSizeMedian = float(psfStarDeltaSizeMedian)
324  psfStats.psfStarDeltaSizeScatter = float(psfStarDeltaSizeScatter)
325  psfStats.psfStarScaledDeltaSizeScatter = float(psfStarScaledDeltaSizeScatter)
326 
327  return psfStats