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
36 __all__ = (
"ComputeExposureSummaryStatsTask",
"ComputeExposureSummaryStatsConfig")
40 """Config for ComputeExposureSummaryTask"""
41 sigmaClip = pexConfig.Field(
43 doc=
"Sigma for outlier rejection for sky noise.",
46 clipIter = pexConfig.Field(
48 doc=
"Number of iterations of outlier rejection for sky noise.",
51 badMaskPlanes = pexConfig.ListField(
53 doc=
"Mask planes that, if set, the associated pixel should not be included sky noise calculation.",
54 default=(
"NO_DATA",
"SUSPECT"),
56 starSelection = pexConfig.Field(
57 doc=
"Field to select sources to be used in the PSF statistics computation.",
59 default=
"calib_psf_used"
61 starShape = pexConfig.Field(
62 doc=
"Base name of columns to use for the source shape in the PSF statistics computation.",
64 default=
"base_SdssShape"
66 psfShape = pexConfig.Field(
67 doc=
"Base name of columns to use for the PSF shape in the PSF statistics computation.",
69 default=
"base_SdssShape_psf"
74 """Task to compute exposure summary statistics.
76 This task computes various quantities suitable for DPDD and other
77 downstream processing at the detector centers, including:
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
104 ConfigClass = ComputeExposureSummaryStatsConfig
105 _DefaultName =
"computeExposureSummaryStats"
108 def run(self, exposure, sources, background):
109 """Measure exposure statistics from the exposure, sources, and background.
113 exposure : `lsst.afw.image.ExposureF`
114 sources : `lsst.afw.table.SourceCatalog`
115 background : `lsst.afw.math.BackgroundList`
119 summary : `lsst.afw.image.ExposureSummary`
121 self.log.info(
"Measuring exposure statistics")
123 bbox = exposure.getBBox()
125 psf = exposure.getPsf()
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())
137 psfArea = np.sum(im.array)/np.sum(im.array**2.)
145 wcs = exposure.getWcs()
148 raCorners = [float(sph.getRa().asDegrees())
for sph
in sph_pts]
149 decCorners = [float(sph.getDec().asDegrees())
for sph
in sph_pts]
151 sph_pt = wcs.pixelToSky(bbox.getCenter())
152 ra = sph_pt.getRa().asDegrees()
153 decl = sph_pt.getDec().asDegrees()
155 raCorners = [float(np.nan)]*4
156 decCorners = [float(np.nan)]*4
160 photoCalib = exposure.getPhotoCalib()
161 if photoCalib
is not None:
162 zeroPoint = 2.5*np.log10(photoCalib.getInstFluxAtZeroMagnitude())
166 visitInfo = exposure.getInfo().getVisitInfo()
167 date = visitInfo.getDate()
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)
184 zenithDistance = 90.0 - altaz.alt.degree
186 zenithDistance = np.nan
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)
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)
201 statObj = afwMath.makeStatistics(exposure.getMaskedImage(), afwMath.STDEVCLIP,
203 skyNoise, _ = statObj.getResult(afwMath.STDEVCLIP)
205 statObj = afwMath.makeStatistics(exposure.getMaskedImage().getVariance(),
206 exposure.getMaskedImage().getMask(),
207 afwMath.MEANCLIP, statsCtrl)
208 meanVar, _ = statObj.getResult(afwMath.MEANCLIP)
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']
215 astromOffsetMean = np.nan
216 astromOffsetStd = np.nan
222 summary = afwImage.ExposureSummaryStats(
223 psfSigma=float(psfSigma),
224 psfArea=float(psfArea),
225 psfIxx=float(psfIxx),
226 psfIyy=float(psfIyy),
227 psfIxy=float(psfIxy),
230 zenithDistance=float(zenithDistance),
231 zeroPoint=float(zeroPoint),
233 skyNoise=float(skyNoise),
234 meanVar=float(meanVar),
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
251 def _computePsfStats(self, sources):
252 """Compute psf residual statistics.
254 All residuals are computed using median statistics on the difference
255 between the sources and the models.
259 sources : `lsst.afw.table.SourceCatalog`
260 Source catalog on which to measure the PSF statistics.
264 psfStats : `lsst.pipe.base.Struct`
265 Struct with various psf stats.
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))
280 names = sources.schema.getNames()
281 if self.config.starSelection
not in names
or self.config.starShape +
'_flag' not in names:
285 mask = sources[self.config.starSelection] & (~sources[self.config.starShape +
'_flag'])
286 nPsfStar = mask.sum()
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]
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)
305 psfSize = (psfXX*psfYY - psfXY**2)**0.25
306 psfE1 = (psfXX - psfYY)/(psfXX + psfYY)
307 psfE2 = 2*psfXY/(psfXX + psfYY)
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')
314 psfStarDeltaSizeMedian = np.median(starSize - psfSize)
315 psfStarDeltaSizeScatter =
sigmaMad(starSize - psfSize, scale=
'normal')
316 psfStarScaledDeltaSizeScatter = psfStarDeltaSizeScatter/starSizeMedian**2.
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)
def run(self, exposure, sources, background)
def _computePsfStats(self, sources)