lsst.pipe.tasks  21.0.0-167-gb57ef11d+e02347d8c1
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 import astropy.units as units
24 from astropy.time import Time
25 from astropy.coordinates import AltAz, SkyCoord, EarthLocation
26 from lsst.daf.base import DateTime
27 
28 import lsst.pipe.base as pipeBase
29 import lsst.pex.config as pexConfig
30 import lsst.afw.math as afwMath
31 import lsst.afw.image as afwImage
32 import lsst.geom
33 from lsst.utils.timer import timeMethod
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 
57 
58 class ComputeExposureSummaryStatsTask(pipeBase.Task):
59  """Task to compute exposure summary statistics.
60 
61  This task computes various quantities suitable for DPDD and other
62  downstream processing at the detector centers, including:
63  - psfSigma
64  - psfArea
65  - psfIxx
66  - psfIyy
67  - psfIxy
68  - ra
69  - decl
70  - zenithDistance
71  - zeroPoint
72  - skyBg
73  - skyNoise
74  - meanVar
75  - raCorners
76  - decCorners
77  - astromOffsetMean
78  - astromOffsetStd
79  """
80  ConfigClass = ComputeExposureSummaryStatsConfig
81  _DefaultName = "computeExposureSummaryStats"
82 
83  @timeMethod
84  def run(self, exposure, sources, background):
85  """Measure exposure statistics from the exposure, sources, and background.
86 
87  Parameters
88  ----------
89  exposure : `lsst.afw.image.ExposureF`
90  sources : `lsst.afw.table.SourceCatalog`
91  background : `lsst.afw.math.BackgroundList`
92 
93  Returns
94  -------
95  summary : `lsst.afw.image.ExposureSummary`
96  """
97  self.log.info("Measuring exposure statistics")
98 
99  bbox = exposure.getBBox()
100 
101  psf = exposure.getPsf()
102  if psf is not None:
103  shape = psf.computeShape(bbox.getCenter())
104  psfSigma = shape.getDeterminantRadius()
105  psfIxx = shape.getIxx()
106  psfIyy = shape.getIyy()
107  psfIxy = shape.getIxy()
108  im = psf.computeKernelImage(bbox.getCenter())
109  # The calculation of effective psf area is taken from
110  # meas_base/src/PsfFlux.cc#L112. See
111  # https://github.com/lsst/meas_base/blob/
112  # 750bffe6620e565bda731add1509507f5c40c8bb/src/PsfFlux.cc#L112
113  psfArea = np.sum(im.array)/np.sum(im.array**2.)
114  else:
115  psfSigma = np.nan
116  psfIxx = np.nan
117  psfIyy = np.nan
118  psfIxy = np.nan
119  psfArea = np.nan
120 
121  wcs = exposure.getWcs()
122  if wcs is not None:
123  sph_pts = wcs.pixelToSky(lsst.geom.Box2D(bbox).getCorners())
124  raCorners = [float(sph.getRa().asDegrees()) for sph in sph_pts]
125  decCorners = [float(sph.getDec().asDegrees()) for sph in sph_pts]
126 
127  sph_pt = wcs.pixelToSky(bbox.getCenter())
128  ra = sph_pt.getRa().asDegrees()
129  decl = sph_pt.getDec().asDegrees()
130  else:
131  raCorners = [float(np.nan)]*4
132  decCorners = [float(np.nan)]*4
133  ra = np.nan
134  decl = np.nan
135 
136  photoCalib = exposure.getPhotoCalib()
137  if photoCalib is not None:
138  zeroPoint = 2.5*np.log10(photoCalib.getInstFluxAtZeroMagnitude())
139  else:
140  zeroPoint = np.nan
141 
142  visitInfo = exposure.getInfo().getVisitInfo()
143  date = visitInfo.getDate()
144 
145  if date.isValid():
146  # We compute the zenithDistance at the center of the detector rather
147  # than use the boresight value available via the visitInfo, because
148  # the zenithDistance may vary significantly over a large field of view.
149  observatory = visitInfo.getObservatory()
150  loc = EarthLocation(lat=observatory.getLatitude().asDegrees()*units.deg,
151  lon=observatory.getLongitude().asDegrees()*units.deg,
152  height=observatory.getElevation()*units.m)
153  obstime = Time(visitInfo.getDate().get(system=DateTime.MJD),
154  location=loc, format='mjd')
155  coord = SkyCoord(ra*units.degree, decl*units.degree, obstime=obstime, location=loc)
156  with warnings.catch_warnings():
157  warnings.simplefilter('ignore')
158  altaz = coord.transform_to(AltAz)
159 
160  zenithDistance = 90.0 - altaz.alt.degree
161  else:
162  zenithDistance = np.nan
163 
164  if background is not None:
165  bgStats = (bg[0].getStatsImage().getImage().array
166  for bg in background)
167  skyBg = sum(np.median(bg[np.isfinite(bg)]) for bg in bgStats)
168  else:
169  skyBg = np.nan
170 
171  statsCtrl = afwMath.StatisticsControl()
172  statsCtrl.setNumSigmaClip(self.config.sigmaClip)
173  statsCtrl.setNumIter(self.config.clipIter)
174  statsCtrl.setAndMask(afwImage.Mask.getPlaneBitMask(self.config.badMaskPlanes))
175  statsCtrl.setNanSafe(True)
176 
177  statObj = afwMath.makeStatistics(exposure.getMaskedImage(), afwMath.STDEVCLIP,
178  statsCtrl)
179  skyNoise, _ = statObj.getResult(afwMath.STDEVCLIP)
180 
181  statObj = afwMath.makeStatistics(exposure.getMaskedImage().getVariance(),
182  exposure.getMaskedImage().getMask(),
183  afwMath.MEANCLIP, statsCtrl)
184  meanVar, _ = statObj.getResult(afwMath.MEANCLIP)
185 
186  md = exposure.getMetadata()
187  if 'SFM_ASTROM_OFFSET_MEAN' in md:
188  astromOffsetMean = md['SFM_ASTROM_OFFSET_MEAN']
189  astromOffsetStd = md['SFM_ASTROM_OFFSET_STD']
190  else:
191  astromOffsetMean = np.nan
192  astromOffsetStd = np.nan
193 
194  summary = afwImage.ExposureSummaryStats(psfSigma=float(psfSigma),
195  psfArea=float(psfArea),
196  psfIxx=float(psfIxx),
197  psfIyy=float(psfIyy),
198  psfIxy=float(psfIxy),
199  ra=float(ra),
200  decl=float(decl),
201  zenithDistance=float(zenithDistance),
202  zeroPoint=float(zeroPoint),
203  skyBg=float(skyBg),
204  skyNoise=float(skyNoise),
205  meanVar=float(meanVar),
206  raCorners=raCorners,
207  decCorners=decCorners,
208  astromOffsetMean=astromOffsetMean,
209  astromOffsetStd=astromOffsetStd)
210 
211  return summary