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