LSST Applications g0b6bd0c080+a72a5dd7e6,g1182afd7b4+2a019aa3bb,g17e5ecfddb+2b8207f7de,g1d67935e3f+06cf436103,g38293774b4+ac198e9f13,g396055baef+6a2097e274,g3b44f30a73+6611e0205b,g480783c3b1+98f8679e14,g48ccf36440+89c08d0516,g4b93dc025c+98f8679e14,g5c4744a4d9+a302e8c7f0,g613e996a0d+e1c447f2e0,g6c8d09e9e7+25247a063c,g7271f0639c+98f8679e14,g7a9cd813b8+124095ede6,g9d27549199+a302e8c7f0,ga1cf026fa3+ac198e9f13,ga32aa97882+7403ac30ac,ga786bb30fb+7a139211af,gaa63f70f4e+9994eb9896,gabf319e997+ade567573c,gba47b54d5d+94dc90c3ea,gbec6a3398f+06cf436103,gc6308e37c7+07dd123edb,gc655b1545f+ade567573c,gcc9029db3c+ab229f5caf,gd01420fc67+06cf436103,gd877ba84e5+06cf436103,gdb4cecd868+6f279b5b48,ge2d134c3d5+cc4dbb2e3f,ge448b5faa6+86d1ceac1d,gecc7e12556+98f8679e14,gf3ee170dca+25247a063c,gf4ac96e456+ade567573c,gf9f5ea5b4d+ac198e9f13,gff490e6085+8c2580be5c,w.2022.27
LSST Data Management Base Package
quickFrameMeasurement.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
22import numpy as np
23import scipy.ndimage as ndImage
24
25import lsst.afw.detection as afwDetect
26import lsst.afw.table as afwTable
27import lsst.geom as geom
28import lsst.meas.base as measBase
29import lsst.daf.base as dafBase
30import lsst.pipe.base as pipeBase
31import lsst.pex.config as pexConfig
32from lsst.meas.base import MeasurementError
33from lsst.meas.algorithms.installGaussianPsf import InstallGaussianPsfTask
34
35
36class QuickFrameMeasurementTaskConfig(pexConfig.Config):
37 """Config class for the QuickFrameMeasurementTask.
38 """
39 installPsf = pexConfig.ConfigurableField(
40 target=InstallGaussianPsfTask,
41 doc="Task for installing an initial PSF",
42 )
43 maxNonRoundness = pexConfig.Field(
44 dtype=float,
45 doc="Ratio of xx to yy (or vice versa) above which to cut, in order to exclude spectra",
46 default=5.,
47 )
48 maxExtendedness = pexConfig.Field(
49 dtype=float,
50 doc="Max absolute value of xx and yy above which to cut, in order to exclude large/things",
51 default=100,
52 )
53 doExtendednessCut = pexConfig.Field(
54 dtype=bool,
55 doc="Apply the extendeness cut, as definted by maxExtendedness",
56 default=False,
57 )
58 centroidPixelPercentile = pexConfig.Field(
59 dtype=float,
60 doc="The image's percentile value which the centroid must be greater than to pass the final peak"
61 " check. Ignored if doCheckCentroidPixelValue is False",
62 default=90,
63 )
64 doCheckCentroidPixelValue = pexConfig.Field(
65 dtype=bool,
66 doc="Check that the centroid found is actually in the centroidPixelPercentile percentile of the"
67 " image? Set to False for donut images.",
68 default=True,
69 )
70 initialPsfWidth = pexConfig.Field(
71 dtype=float,
72 doc="Guess at the initial PSF FWHM in pixels.",
73 default=10,
74 )
75 nSigmaDetection = pexConfig.Field(
76 dtype=float,
77 doc="Number of sigma for the detection limit.",
78 default=20,
79 )
80 nPixMinDetection = pexConfig.Field(
81 dtype=int,
82 doc="Minimum number of pixels in a detected source.",
83 default=10,
84 )
85 donutDiameter = pexConfig.Field(
86 dtype=int,
87 doc="The expected diameter of donuts in a donut image, in pixels.",
88 default=400,
89 )
90
91 def setDefaults(self):
92 super().setDefaults()
93 self.installPsfinstallPsf.fwhm = self.initialPsfWidthinitialPsfWidth
94
95
96class QuickFrameMeasurementTask(pipeBase.Task):
97 """WARNING: An experimental new task with changable API! Do not rely on yet!
98
99 This task finds the centroid of the brightest source in a given CCD-image
100 and returns its centroid and a rough estimate of the seeing/PSF.
101
102 It is designed for speed, such that it can be used in observing scripts
103 to provide pointing offsets, allowing subsequent pointings to place
104 a source at an exact pixel position.
105
106 The approach taken here is deliberately sub-optimal in the detection and
107 measurement sense, with all optimisation being done for speed and robustness
108 of the result.
109
110 A small set of unit tests exist for this task, which run automatically
111 if afwdata is setup. These, however, are stricky unit tests, and will not
112 catch algorithmic regressions. TODO: DM-29038 exists to merge a regression
113 real test which runs against 1,000 LATISS images, but is therefore slow
114 and requires access to the data.
115
116 Parameters
117 ----------
119 Configuration class for the QuickFrameMeasurementTask.
120 display : `lsst.afw.display.Display`, optional
121 The display to use for showing the images, detections and centroids.
122
123 Returns
124 -------
125 result : `lsst.pipe.base.Struct`
126 Return strucure containing whether the task was successful, the main
127 source's centroid, its the aperture fluxes, the ixx and iyy of the
128 source, and the median ixx, iyy of the detections in the exposure.
129 See run() method for further details.
130
131 Raises
132 ------
133 This task should *never* raise, as the run() method is enclosed in an
134 except Exception block, so that it will never fail during observing.
135 Failure modes should be limited to returning a return Struct() with the same
136 structure as the success case, with all value set to np.nan but with
137 result.success=False.
138 """
139 ConfigClass = QuickFrameMeasurementTaskConfig
140 _DefaultName = 'quickFrameMeasurementTask'
141
142 def __init__(self, config, *, display=None, **kwargs):
143 super().__init__(config=config, **kwargs)
144 self.makeSubtask("installPsf")
145
146 self.displaydisplay = None
147 if display:
148 self.displaydisplay = display
149
150 self.centroidNamecentroidName = "base_SdssCentroid"
151 self.shapeNameshapeName = "base_SdssShape"
152 self.schemaschema = afwTable.SourceTable.makeMinimalSchema()
153 self.schemaschema.getAliasMap().set("slot_Centroid", self.centroidNamecentroidName)
154 self.schemaschema.getAliasMap().set("slot_Shape", self.shapeNameshapeName)
155 self.controlcontrol = measBase.SdssCentroidControl()
156 self.centroidercentroider = measBase.SdssCentroidAlgorithm(self.controlcontrol, self.centroidNamecentroidName, self.schemaschema)
157 self.sdssShapesdssShape = measBase.SdssShapeControl()
158 self.shapershaper = measBase.SdssShapeAlgorithm(self.sdssShapesdssShape, self.shapeNameshapeName, self.schemaschema)
159 self.apFluxControlapFluxControl = measBase.ApertureFluxControl()
161 self.apFluxerapFluxer = measBase.CircularApertureFluxAlgorithm(self.apFluxControlapFluxControl, "aperFlux",
162 self.schemaschema, md)
163
164 self.tabletable = afwTable.SourceTable.make(self.schemaschema) # make sure to call this last!
165
166 @staticmethod
167 def detectObjectsInExp(exp, nSigma, nPixMin, grow=0):
168 """Run a very basic but fast threshold-based object detection on an exposure
169 Return the footPrintSet for the objects in a postISR exposure.
170
171 Parameters
172 ----------
174 Image in which to detect objects.
175 nSigma : `float`
176 nSigma above image's stddev at which to set the detection threshold.
177 nPixMin : `int`
178 Minimum number of pixels for detection.
179 grow : `int`
180 Grow the detected footprint by this many pixels.
181
182 Returns
183 -------
184 footPrintSet : `lsst.afw.detection.FootprintSet`
185 FootprintSet containing the detections.
186 """
187 threshold = afwDetect.Threshold(nSigma, afwDetect.Threshold.STDEV)
188 footPrintSet = afwDetect.FootprintSet(exp.getMaskedImage(), threshold, "DETECTED", nPixMin)
189 if grow > 0:
190 isotropic = True
191 footPrintSet = afwDetect.FootprintSet(footPrintSet, grow, isotropic)
192 return footPrintSet
193
194 @staticmethod
195 def checkResult(exp, centroid, srcNum, percentile):
196 """Perform a final check that centroid location is actually bright.
197
198 Parameters
199 ----------
201 The exposure on which to operate
202 centroid : `tuple` of `float`
203 Location of the centroid in pixel coordinates
204 scrNum : `int`
205 Number of the source in the source catalog. Only used if the check
206 is failed, for debug purposes.
207 percentile : `float`
208 Image's percentile above which the pixel containing the centroid
209 must be in order to pass the check.
210
211 Raises
212 ------
213 ValueError
214 Raised if the centroid's pixel is not above the percentile threshold
215 """
216 threshold = np.percentile(exp.image.array, percentile)
217 pixelValue = exp.image[centroid]
218 if pixelValue < threshold:
219 msg = (f"Final centroid pixel value check failed: srcNum {srcNum} at {centroid}"
220 f" has central pixel = {pixelValue:3f} <"
221 f" {percentile} percentile of image = {threshold:3f}")
222 raise ValueError(msg)
223 return
224
225 @staticmethod
226 def _calcMedianXxYy(objData):
227 """Return the median ixx and iyy for object in the image.
228 """
229 medianXx = np.nanmedian([element['xx'] for element in objData.values()])
230 medianYy = np.nanmedian([element['yy'] for element in objData.values()])
231 return medianXx, medianYy
232
233 @staticmethod
234 def _getCenterOfMass(exp, nominalCentroid, boxSize):
235 """Get the centre of mass around a point in the image.
236
237 Parameters
238 ----------
240 The exposure in question.
241 nominalCentroid : `tuple` of `float`
242 Nominal location of the centroid in pixel coordinates.
243 boxSize : `int`
244 The size of the box around the nominalCentroid in which to measure
245 the centre of mass.
246
247 Returns
248 -------
249 com : `tuple` of `float`
250 The locaiton of the centre of mass of the brightest source in pixel
251 coordinates.
252 """
253 centroidPoint = geom.Point2I(nominalCentroid)
254 extent = geom.Extent2I(1, 1)
255 bbox = geom.Box2I(centroidPoint, extent)
256 bbox = bbox.dilatedBy(int(boxSize//2))
257 bbox = bbox.clippedTo(exp.getBBox())
258 data = exp[bbox].image.array
259 xy0 = exp[bbox].getXY0()
260
261 peak = ndImage.center_of_mass(data)
262 peak = (peak[1], peak[0]) # numpy coords returned
263 com = geom.Point2D(xy0)
264 com.shift(geom.Extent2D(*peak))
265 return (com[0], com[1])
266
267 def _calcBrightestObjSrcNum(self, objData):
268 """Find the brightest source which passes the cuts among the sources.
269
270 Parameters
271 ----------
272 objData : `dict` of `dict`
273 Dictionary, keyed by source number, containing the measurements.
274
275 Returns
276 -------
277 srcNum : `int`
278 The source number of the brightest source which passes the cuts.
279 """
280 max70, max70srcNum = -1, -1
281 max25, max25srcNum = -1, -1
282
283 for srcNum in sorted(objData.keys()): # srcNum not contiguous so don't use a list comp
284 # skip flag used rather than continue statements so we have all the
285 # metrics computed for debug purposes as this task is whack-a-mole
286 skip = False
287 xx = objData[srcNum]['xx']
288 yy = objData[srcNum]['yy']
289
290 xx = max(xx, 1e-9) # need to protect against division by zero
291 yy = max(yy, 1e-9) # because we don't `continue` on zero moments
292
293 if self.config.doExtendednessCut:
294 if xx > self.config.maxExtendedness or yy > self.config.maxExtendedness:
295 skip = True
296
297 nonRoundness = xx/yy
298 nonRoundness = max(nonRoundness, 1/nonRoundness)
299 if nonRoundness > self.config.maxNonRoundness:
300 skip = True
301
302 if self.log.isEnabledFor(self.log.DEBUG):
303 text = f"src {srcNum}: {objData[srcNum]['xCentroid']:.0f}, {objData[srcNum]['yCentroid']:.0f}"
304 text += f" - xx={xx:.1f}, yy={yy:.1f}, nonRound={nonRoundness:.1f}"
305 text += f" - ap70={objData[srcNum]['apFlux70']:,.0f}"
306 text += f" - ap25={objData[srcNum]['apFlux25']:,.0f}"
307 text += f" - skip={skip}"
308 self.log.debug(text)
309
310 if skip:
311 continue
312
313 ap70 = objData[srcNum]['apFlux70']
314 ap25 = objData[srcNum]['apFlux25']
315 if ap70 > max70:
316 max70 = ap70
317 max70srcNum = srcNum
318 if ap25 > max25:
319 max25 = ap25
320 max25srcNum = srcNum
321 if max70srcNum != max25srcNum:
322 self.log.warning("WARNING! Max apFlux70 for different object than with max apFlux25")
323
324 if max70srcNum >= 0: # starts as -1, return None if nothing is acceptable
325 return max70srcNum
326 return None
327
328 def _measureFp(self, fp, exp):
329 """Run the measurements on a footprint.
330
331 Parameters
332 ----------
334 The footprint to measure.
336 The footprint's parent exposure.
337
338 Returns
339 -------
341 The source record containing the measurements.
342 """
343 src = self.tabletable.makeRecord()
344 src.setFootprint(fp)
345 self.centroidercentroider.measure(src, exp)
346 self.shapershaper.measure(src, exp)
347 self.apFluxerapFluxer.measure(src, exp)
348 return src
349
350 def _getDataFromSrcRecord(self, src):
351 """Extract the shapes and centroids from a source record.
352
353 Parameters
354 ----------
356 The source record from which to extract the measurements.
357
358 Returns
359 -------
360 srcData : `lsst.pipe.base.Struct`
361 The struct containing the extracted measurements.
362 """
363 pScale = self.plateScaleplateScale
364 xx = np.sqrt(src['base_SdssShape_xx'])*2.355*pScale # 2.355 for FWHM, pScale for platescale from exp
365 yy = np.sqrt(src['base_SdssShape_yy'])*2.355*pScale
366 xCentroid = src['base_SdssCentroid_x']
367 yCentroid = src['base_SdssCentroid_y']
368 # apFluxes available: 70, 50, 35, 25, 17, 12 9, 6, 4.5, 3
369 apFlux70 = src['aperFlux_70_0_instFlux']
370 apFlux25 = src['aperFlux_25_0_instFlux']
371 return pipeBase.Struct(xx=xx,
372 yy=yy,
373 xCentroid=xCentroid,
374 yCentroid=yCentroid,
375 apFlux70=apFlux70,
376 apFlux25=apFlux25)
377
378 @staticmethod
379 def _getDataFromFootprintOnly(fp, exp):
380 """Get the shape, centroid and flux from a footprint.
381
382 Parameters
383 ----------
385 The footprint to measure.
387 The footprint's parent exposure.
388
389 Returns
390 -------
391 srcData : `lsst.pipe.base.Struct`
392 The struct containing the extracted measurements.
393 """
394 xx = fp.getShape().getIxx()
395 yy = fp.getShape().getIyy()
396 xCentroid, yCentroid = fp.getCentroid()
397 apFlux70 = np.sum(exp[fp.getBBox()].image.array)
398 apFlux25 = np.sum(exp[fp.getBBox()].image.array)
399 return pipeBase.Struct(xx=xx,
400 yy=yy,
401 xCentroid=xCentroid,
402 yCentroid=yCentroid,
403 apFlux70=apFlux70,
404 apFlux25=apFlux25)
405
406 @staticmethod
407 def _measurementResultToDict(measurementResult):
408 """Convenience function to repackage measurement results to a dict.
409
410 Parameters
411 ----------
412 measurementResult : `lsst.afw.table.SourceRecord`
413 The source record to convert to a dict.
414
415 Returns
416 -------
417 objData : `dict`
418 The dict containing the extracted data.
419 """
420 objData = {}
421 objData['xx'] = measurementResult.xx
422 objData['yy'] = measurementResult.yy
423 objData['xCentroid'] = measurementResult.xCentroid
424 objData['yCentroid'] = measurementResult.yCentroid
425 objData['apFlux70'] = measurementResult.apFlux70
426 objData['apFlux25'] = measurementResult.apFlux25
427 return objData
428
429 @staticmethod
430 def _makeEmptyReturnStruct():
431 """Make the default/template return struct, with defaults to False/nan.
432
433 Returns
434 -------
435 objData : `lsst.pipe.base.Struct`
436 The default template return structure.
437 """
438 result = pipeBase.Struct()
439 result.success = False
440 result.brightestObjCentroid = (np.nan, np.nan)
441 result.brightestObjCentroidCofM = None
442 result.brightestObj_xXyY = (np.nan, np.nan)
443 result.brightestObjApFlux70 = np.nan
444 result.brightestObjApFlux25 = np.nan
445 result.medianXxYy = (np.nan, np.nan)
446 return result
447
448 def run(self, exp, *, donutDiameter=None, doDisplay=False):
449 """Calculate position, flux and shape of the brightest star in an image.
450
451 Given an an assembled (and at least minimally ISRed exposure),
452 quickly and robustly calculate the centroid of the
453 brightest star in the image.
454
455 Parameters
456 ----------
458 The exposure in which to find and measure the brightest star.
459 donutDiameter : `int` or `float`, optional
460 The expected diameter of donuts in pixels for use in the centre of
461 mass centroid measurement. If None is provided, the config option
462 is used.
463 doDisplay : `bool`
464 Display the image and found sources. A diplay object must have
465 been passed to the task constructor.
466
467 Returns
468 -------
469 result : `lsst.pipe.base.Struct`
470 Struct containing:
471 Whether the task ran successfully and found the object (bool)
472 The object's centroid (float, float)
473 The object's ixx, iyy (float, float)
474 The object's 70 pixel aperture flux (float)
475 The object's 25 pixel aperture flux (float)
476 The images's median ixx, iyy (float, float)
477 If unsuccessful, the success field is False and all other results
478 are np.nan of the expected shape.
479
480 Notes
481 -----
482 Because of this task's involvement in observing scripts, the run method
483 should *never* raise. Failure modes are noted by returning a Struct with
484 the same structure as the success case, with all value set to np.nan and
485 result.success=False.
486 """
487 try:
488 result = self._run_run(exp=exp, donutDiameter=donutDiameter, doDisplay=doDisplay)
489 return result
490 except Exception as e:
491 self.log.warning("Failed to find main source centroid %s", e)
492 result = self._makeEmptyReturnStruct_makeEmptyReturnStruct()
493 return result
494
495 def _run(self, exp, *, donutDiameter=None, doDisplay=False):
496 """The actual run method, called by run()
497
498 Behaviour is documented in detail in the main run().
499 """
500 if donutDiameter is None:
501 donutDiameter = self.config.donutDiameter
502
503 self.plateScaleplateScale = exp.getWcs().getPixelScale().asArcseconds()
504 median = np.nanmedian(exp.image.array)
505 exp.image -= median # is put back later
506 self.installPsf.run(exp)
507 sources = self.detectObjectsInExpdetectObjectsInExp(exp, nSigma=self.config.nSigmaDetection,
508 nPixMin=self.config.nPixMinDetection)
509
510 if doDisplay:
511 if self.displaydisplay is None:
512 raise RuntimeError("Display failed as no display provided during init()")
513 self.displaydisplay.mtv(exp)
514
515 fpSet = sources.getFootprints()
516 self.log.info("Found %d sources in exposure", len(fpSet))
517
518 objData = {}
519 nMeasured = 0
520
521 for srcNum, fp in enumerate(fpSet):
522 try:
523 src = self._measureFp_measureFp(fp, exp)
524 result = self._getDataFromSrcRecord_getDataFromSrcRecord(src)
525 except MeasurementError:
526 try:
527 # gets shape and centroid from footprint
528 result = self._getDataFromFootprintOnly_getDataFromFootprintOnly(fp, exp)
529 except MeasurementError as e:
530 self.log.info("Skipped measuring source %s: %s", srcNum, e)
531 continue
532 objData[srcNum] = self._measurementResultToDict_measurementResultToDict(result)
533 nMeasured += 1
534
535 self.log.info("Measured %d of %d sources in exposure", nMeasured, len(fpSet))
536
537 medianXxYy = self._calcMedianXxYy_calcMedianXxYy(objData)
538
539 brightestObjSrcNum = self._calcBrightestObjSrcNum_calcBrightestObjSrcNum(objData)
540 if brightestObjSrcNum is None:
541 raise RuntimeError("No sources in image passed cuts")
542
543 x = objData[brightestObjSrcNum]['xCentroid']
544 y = objData[brightestObjSrcNum]['yCentroid']
545 brightestObjCentroid = (x, y)
546 xx = objData[brightestObjSrcNum]['xx']
547 yy = objData[brightestObjSrcNum]['yy']
548 brightestObjApFlux70 = objData[brightestObjSrcNum]['apFlux70']
549 brightestObjApFlux25 = objData[brightestObjSrcNum]['apFlux25']
550
551 exp.image += median # put background back in
552 if self.config.doCheckCentroidPixelValue:
553 self.checkResultcheckResult(exp, brightestObjCentroid, brightestObjSrcNum,
554 self.config.centroidPixelPercentile)
555
556 boxSize = donutDiameter * 1.3 # allow some slack, as cutting off side of donut is very bad
557 centreOfMass = self._getCenterOfMass_getCenterOfMass(exp, brightestObjCentroid, boxSize)
558
559 result = self._makeEmptyReturnStruct_makeEmptyReturnStruct()
560 result.success = True
561 result.brightestObjCentroid = brightestObjCentroid
562 result.brightestObj_xXyY = (xx, yy)
563 result.brightestObjApFlux70 = brightestObjApFlux70
564 result.brightestObjApFlux25 = brightestObjApFlux25
565 result.medianXxYy = medianXxYy
566 result.brightestObjCentroidCofM = centreOfMass
567
568 return result
int max
Class to describe the properties of a detected object from an image.
Definition: Footprint.h:63
A set of Footprints, associated with a MaskedImage.
Definition: FootprintSet.h:55
A Threshold is used to pass a threshold value to detection algorithms.
Definition: Threshold.h:43
A class to contain the data, WCS, and other information needed to describe an image of the sky.
Definition: Exposure.h:72
Record class that contains measurements made on a single exposure.
Definition: Source.h:78
Class for storing generic metadata.
Definition: PropertySet.h:66
An integer coordinate rectangle.
Definition: Box.h:55
def run(self, exp, *donutDiameter=None, doDisplay=False)
def _run(self, exp, *donutDiameter=None, doDisplay=False)
daf::base::PropertySet * set
Definition: fits.cc:912
def mtv(data, frame=None, title="", wcs=None, *args, **kwargs)
Definition: ds9.py:92
def measure(mi, x, y, size, statistic, stats)
Definition: fringe.py:519
def isEnabledFor(loggername, level)