LSST Applications g0f08755f38+82efc23009,g12f32b3c4e+e7bdf1200e,g1653933729+a8ce1bb630,g1a0ca8cf93+50eff2b06f,g28da252d5a+52db39f6a5,g2bbee38e9b+37c5a29d61,g2bc492864f+37c5a29d61,g2cdde0e794+c05ff076ad,g3156d2b45e+41e33cbcdc,g347aa1857d+37c5a29d61,g35bb328faa+a8ce1bb630,g3a166c0a6a+37c5a29d61,g3e281a1b8c+fb992f5633,g414038480c+7f03dfc1b0,g41af890bb2+11b950c980,g5fbc88fb19+17cd334064,g6b1c1869cb+12dd639c9a,g781aacb6e4+a8ce1bb630,g80478fca09+72e9651da0,g82479be7b0+04c31367b4,g858d7b2824+82efc23009,g9125e01d80+a8ce1bb630,g9726552aa6+8047e3811d,ga5288a1d22+e532dc0a0b,gae0086650b+a8ce1bb630,gb58c049af0+d64f4d3760,gc28159a63d+37c5a29d61,gcf0d15dbbd+2acd6d4d48,gd7358e8bfb+778a810b6e,gda3e153d99+82efc23009,gda6a2b7d83+2acd6d4d48,gdaeeff99f8+1711a396fd,ge2409df99d+6b12de1076,ge79ae78c31+37c5a29d61,gf0baf85859+d0a5978c5a,gf3967379c6+4954f8c433,gfb92a5be7c+82efc23009,gfec2e1e490+2aaed99252,w.2024.46
LSST Data Management Base Package
Loading...
Searching...
No Matches
dipoleFitTask.py
Go to the documentation of this file.
2# LSST Data Management System
3# Copyright 2008-2016 AURA/LSST.
4#
5# This product includes software developed by the
6# LSST Project (http://www.lsst.org/).
7#
8# This program is free software: you can redistribute it and/or modify
9# it under the terms of the GNU General Public License as published by
10# the Free Software Foundation, either version 3 of the License, or
11# (at your option) any later version.
12#
13# This program is distributed in the hope that it will be useful,
14# but WITHOUT ANY WARRANTY; without even the implied warranty of
15# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16# GNU General Public License for more details.
17#
18# You should have received a copy of the LSST License Statement and
19# the GNU General Public License along with this program. If not,
20# see <https://www.lsstcorp.org/LegalNotices/>.
21#
22
23import logging
24import numpy as np
25import warnings
26
27import lsst.afw.image as afwImage
28import lsst.meas.base as measBase
29import lsst.afw.detection as afwDet
30import lsst.geom as geom
31import lsst.pex.exceptions as pexExcept
32import lsst.pex.config as pexConfig
33from lsst.pipe.base import Struct
34from lsst.utils.timer import timeMethod
35
36__all__ = ("DipoleFitTask", "DipoleFitPlugin", "DipoleFitTaskConfig", "DipoleFitPluginConfig",
37 "DipoleFitAlgorithm")
38
39
40# Create a new measurement task (`DipoleFitTask`) that can handle all other SFM tasks but can
41# pass a separate pos- and neg- exposure/image to the `DipoleFitPlugin`s `run()` method.
42
43
44class DipoleFitPluginConfig(measBase.SingleFramePluginConfig):
45 """Configuration for DipoleFitPlugin
46 """
47
48 fitAllDiaSources = pexConfig.Field(
49 dtype=bool, default=False,
50 doc="""Attempte dipole fit of all diaSources (otherwise just the ones consisting of overlapping
51 positive and negative footprints)""")
52
53 maxSeparation = pexConfig.Field(
54 dtype=float, default=5.,
55 doc="Assume dipole is not separated by more than maxSeparation * psfSigma")
56
57 relWeight = pexConfig.Field(
58 dtype=float, default=0.5,
59 doc="""Relative weighting of pre-subtraction images (higher -> greater influence of pre-sub.
60 images on fit)""")
61
62 tolerance = pexConfig.Field(
63 dtype=float, default=1e-7,
64 doc="Fit tolerance")
65
66 fitBackground = pexConfig.Field(
67 dtype=int, default=1,
68 doc="Set whether and how to fit for linear gradient in pre-sub. images. Possible values:"
69 "0: do not fit background at all"
70 "1 (default): pre-fit the background using linear least squares and then do not fit it as part"
71 "of the dipole fitting optimization"
72 "2: pre-fit the background using linear least squares (as in 1), and use the parameter"
73 "estimates from that fit as starting parameters for an integrated re-fit of the background")
74
75 fitSeparateNegParams = pexConfig.Field(
76 dtype=bool, default=False,
77 doc="Include parameters to fit for negative values (flux, gradient) separately from pos.")
78
79 # Config params for classification of detected diaSources as dipole or not
80 minSn = pexConfig.Field(
81 dtype=float, default=np.sqrt(2) * 5.0,
82 doc="Minimum quadrature sum of positive+negative lobe S/N to be considered a dipole")
83
84 maxFluxRatio = pexConfig.Field(
85 dtype=float, default=0.65,
86 doc="Maximum flux ratio in either lobe to be considered a dipole")
87
88 maxChi2DoF = pexConfig.Field(
89 dtype=float, default=0.05,
90 doc="""Maximum Chi2/DoF significance of fit to be considered a dipole.
91 Default value means \"Choose a chi2DoF corresponding to a significance level of at most 0.05\"
92 (note this is actually a significance, not a chi2 value).""")
93
94 maxFootprintArea = pexConfig.Field(
95 dtype=int, default=1_200,
96 doc=("Maximum area for footprints before they are ignored as large; "
97 "non-positive means no threshold applied"
98 "Threshold chosen for HSC and DECam data, see DM-38741 for details."))
99
100
101class DipoleFitTaskConfig(measBase.SingleFrameMeasurementConfig):
102
103 def setDefaults(self):
104 measBase.SingleFrameMeasurementConfig.setDefaults(self)
105
106 self.plugins.names = ["base_SdssCentroid",
107 "ip_diffim_DipoleFit",
108 "base_CircularApertureFlux",
109 "base_PixelFlags",
110 "base_SkyCoord",
111 "base_PsfFlux",
112 ]
113 # Only measure the apertures we need to report in the alert stream.
114 self.plugins["base_CircularApertureFlux"].radii = [12.0]
115
116 self.slots.calibFlux = None
117 self.slots.modelFlux = None
118 self.slots.gaussianFlux = None
119 self.slots.shape = None
120 # This will be switched to "ip_diffim_DipoleFit" as this task runs.
121 self.slots.centroid = "base_SdssCentroid"
123
124
125class DipoleFitTask(measBase.SingleFrameMeasurementTask):
126 """A task that fits a dipole to a difference image, with an optional
127 separate detection image.
128
129 Because it subclasses SingleFrameMeasurementTask, and calls
130 SingleFrameMeasurementTask.run() from its run() method, it still
131 can be used identically to a standard SingleFrameMeasurementTask.
132 """
133
134 ConfigClass = DipoleFitTaskConfig
135 _DefaultName = "dipoleFit"
136
137 def __init__(self, schema, algMetadata=None, **kwargs):
138 super().__init__(schema, algMetadata, **kwargs)
139
140 # Enforce a specific plugin order, so that DipoleFit can fall back on
141 # SdssCentroid for non-dipoles
142 self.plugins_pre = self.plugins.copy()
143 self.plugins_post = self.plugins.copy()
144 self.plugins_pre.clear()
145 self.plugins_pre["base_SdssCentroid"] = self.plugins["base_SdssCentroid"]
146 self.plugins_post.pop("base_SdssCentroid")
147 self.dipoleFit = self.plugins_post.pop("ip_diffim_DipoleFit")
148 del self.plugins
149
150 @timeMethod
151 def run(self, sources, exposure, posExp=None, negExp=None, **kwargs):
152 """Run dipole measurement and classification.
153
154 Run SdssCentroid first, then switch the centroid slot, then DipoleFit
155 then the rest; DipoleFit will fall back on SdssCentroid for sources
156 not containing positive+negative peaks.
157
158 Parameters
159 ----------
160 sources : `lsst.afw.table.SourceCatalog`
161 ``diaSources`` that will be measured using dipole measurement.
162 exposure : `lsst.afw.image.Exposure`
163 The difference exposure on which the ``sources`` were detected.
164 If neither ``posExp`` nor ``negExp`` are set, then the dipole is also
165 fitted directly to this difference image.
166 posExp : `lsst.afw.image.Exposure`, optional
167 "Positive" exposure, typically a science exposure, or None if unavailable
168 When `posExp` is `None`, will compute `posImage = exposure + negExp`.
169 negExp : `lsst.afw.image.Exposure`, optional
170 "Negative" exposure, typically a template exposure, or None if unavailable
171 When `negExp` is `None`, will compute `negImage = posExp - exposure`.
172 **kwargs
173 Additional keyword arguments for `lsst.meas.base.sfm.SingleFrameMeasurementTask`.
174 """
175 # Run plugins in a very specific order, so DipoleFitPlugin has a
176 # centroid to fall back on.
178 super().run(sources, exposure, **kwargs)
179
180 for source in sources:
181 self.dipoleFit.measureDipoles(source, exposure, posExp, negExp)
182 # Use the new DipoleFit outputs for subsequent measurements, now that
183 # non-dipoles have been filled in with the earlier centroid values.
184 sources.schema.getAliasMap().set("slot_Centroid", "ip_diffim_DipoleFit")
185
186 self.plugins = self.plugins_post
187 super().run(sources, exposure, **kwargs)
188
189
191 """Lightweight class containing methods for generating a dipole model for fitting
192 to sources in diffims, used by DipoleFitAlgorithm.
193
194 See also:
195 `DMTN-007: Dipole characterization for image differencing <https://dmtn-007.lsst.io>`_.
196 """
197
198 def __init__(self):
199 import lsstDebug
200 self.debug = lsstDebug.Info(__name__).debug
201 self.log = logging.getLogger(__name__)
202
203 def makeBackgroundModel(self, in_x, pars=None):
204 """Generate gradient model (2-d array) with up to 2nd-order polynomial
205
206 Parameters
207 ----------
208 in_x : `numpy.array`
209 (2, w, h)-dimensional `numpy.array`, containing the
210 input x,y meshgrid providing the coordinates upon which to
211 compute the gradient. This will typically be generated via
212 `_generateXYGrid()`. `w` and `h` correspond to the width and
213 height of the desired grid.
214 pars : `list` of `float`, optional
215 Up to 6 floats for up
216 to 6 2nd-order 2-d polynomial gradient parameters, in the
217 following order: (intercept, x, y, xy, x**2, y**2). If `pars`
218 is emtpy or `None`, do nothing and return `None` (for speed).
219
220 Returns
221 -------
222 result : `None` or `numpy.array`
223 return None, or 2-d numpy.array of width/height matching
224 input bbox, containing computed gradient values.
225 """
226
227 # Don't fit for other gradient parameters if the intercept is not included.
228 if (pars is None) or (len(pars) <= 0) or (pars[0] is None):
229 return
230
231 y, x = in_x[0, :], in_x[1, :]
232 gradient = np.full_like(x, pars[0], dtype='float64')
233 if len(pars) > 1 and pars[1] is not None:
234 gradient += pars[1] * x
235 if len(pars) > 2 and pars[2] is not None:
236 gradient += pars[2] * y
237 if len(pars) > 3 and pars[3] is not None:
238 gradient += pars[3] * (x * y)
239 if len(pars) > 4 and pars[4] is not None:
240 gradient += pars[4] * (x * x)
241 if len(pars) > 5 and pars[5] is not None:
242 gradient += pars[5] * (y * y)
243
244 return gradient
245
246 def _generateXYGrid(self, bbox):
247 """Generate a meshgrid covering the x,y coordinates bounded by bbox
248
249 Parameters
250 ----------
251 bbox : `lsst.geom.Box2I`
252 input Bounding Box defining the coordinate limits
253
254 Returns
255 -------
256 in_x : `numpy.array`
257 (2, w, h)-dimensional numpy array containing the grid indexing over x- and
258 y- coordinates
259 """
260
261 x, y = np.mgrid[bbox.getBeginY():bbox.getEndY(), bbox.getBeginX():bbox.getEndX()]
262 in_x = np.array([y, x]).astype(np.float64)
263 in_x[0, :] -= np.mean(in_x[0, :])
264 in_x[1, :] -= np.mean(in_x[1, :])
265 return in_x
266
267 def _getHeavyFootprintSubimage(self, fp, badfill=np.nan, grow=0):
268 """Extract the image from a ``~lsst.afw.detection.HeavyFootprint``
269 as an `lsst.afw.image.ImageF`.
270
271 Parameters
272 ----------
273 fp : `lsst.afw.detection.HeavyFootprint`
274 HeavyFootprint to use to generate the subimage
275 badfill : `float`, optional
276 Value to fill in pixels in extracted image that are outside the footprint
277 grow : `int`
278 Optionally grow the footprint by this amount before extraction
279
280 Returns
281 -------
282 subim2 : `lsst.afw.image.ImageF`
283 An `~lsst.afw.image.ImageF` containing the subimage.
284 """
285 bbox = fp.getBBox()
286 if grow > 0:
287 bbox.grow(grow)
288
289 subim2 = afwImage.ImageF(bbox, badfill)
290 fp.getSpans().unflatten(subim2.array, fp.getImageArray(), bbox.getCorners()[0])
291 return subim2
292
293 def fitFootprintBackground(self, source, posImage, order=1):
294 """Fit a linear (polynomial) model of given order (max 2) to the background of a footprint.
295
296 Only fit the pixels OUTSIDE of the footprint, but within its bounding box.
297
298 Parameters
299 ----------
300 source : `lsst.afw.table.SourceRecord`
301 SourceRecord, the footprint of which is to be fit
302 posImage : `lsst.afw.image.Exposure`
303 The exposure from which to extract the footprint subimage
304 order : `int`
305 Polynomial order of background gradient to fit.
306
307 Returns
308 -------
309 pars : `tuple` of `float`
310 `tuple` of length (1 if order==0; 3 if order==1; 6 if order == 2),
311 containing the resulting fit parameters
312 """
313
314 # TODO look into whether to use afwMath background methods -- see
315 # http://lsst-web.ncsa.illinois.edu/doxygen/x_masterDoxyDoc/_background_example.html
316 fp = source.getFootprint()
317 bbox = fp.getBBox()
318 bbox.grow(3)
319 posImg = afwImage.ImageF(posImage.image, bbox, afwImage.PARENT)
320
321 # This code constructs the footprint image so that we can identify the pixels that are
322 # outside the footprint (but within the bounding box). These are the pixels used for
323 # fitting the background.
324 posHfp = afwDet.HeavyFootprintF(fp, posImage.getMaskedImage())
325 posFpImg = self._getHeavyFootprintSubimage(posHfp, grow=3)
326
327 isBg = np.isnan(posFpImg.array).ravel()
328
329 data = posImg.array.ravel()
330 data = data[isBg]
331 B = data
332
333 x, y = np.mgrid[bbox.getBeginY():bbox.getEndY(), bbox.getBeginX():bbox.getEndX()]
334 x = x.astype(np.float64).ravel()
335 x -= np.mean(x)
336 x = x[isBg]
337 y = y.astype(np.float64).ravel()
338 y -= np.mean(y)
339 y = y[isBg]
340 b = np.ones_like(x, dtype=np.float64)
341
342 M = np.vstack([b]).T # order = 0
343 if order == 1:
344 M = np.vstack([b, x, y]).T
345 elif order == 2:
346 M = np.vstack([b, x, y, x**2., y**2., x*y]).T
347
348 pars = np.linalg.lstsq(M, B, rcond=-1)[0]
349 return pars
350
351 def makeStarModel(self, bbox, psf, xcen, ycen, flux):
352 """Generate a 2D image model of a single PDF centered at the given coordinates.
353
354 Parameters
355 ----------
356 bbox : `lsst.geom.Box`
357 Bounding box marking pixel coordinates for generated model
358 psf : TODO: DM-17458
359 Psf model used to generate the 'star'
360 xcen : `float`
361 Desired x-centroid of the 'star'
362 ycen : `float`
363 Desired y-centroid of the 'star'
364 flux : `float`
365 Desired flux of the 'star'
366
367 Returns
368 -------
369 p_Im : `lsst.afw.image.Image`
370 2-d stellar image of width/height matching input ``bbox``,
371 containing PSF with given centroid and flux
372 """
373
374 # Generate the psf image, normalize to flux
375 psf_img = psf.computeImage(geom.Point2D(xcen, ycen)).convertF()
376 psf_img_sum = np.nansum(psf_img.array)
377 psf_img *= (flux/psf_img_sum)
378
379 # Clip the PSF image bounding box to fall within the footprint bounding box
380 psf_box = psf_img.getBBox()
381 psf_box.clip(bbox)
382 psf_img = afwImage.ImageF(psf_img, psf_box, afwImage.PARENT)
383
384 # Then actually crop the psf image.
385 # Usually not necessary, but if the dipole is near the edge of the image...
386 # Would be nice if we could compare original pos_box with clipped pos_box and
387 # see if it actually was clipped.
388 p_Im = afwImage.ImageF(bbox)
389 tmpSubim = afwImage.ImageF(p_Im, psf_box, afwImage.PARENT)
390 tmpSubim += psf_img
391
392 return p_Im
393
394 def makeModel(self, x, flux, xcenPos, ycenPos, xcenNeg, ycenNeg, fluxNeg=None,
395 b=None, x1=None, y1=None, xy=None, x2=None, y2=None,
396 bNeg=None, x1Neg=None, y1Neg=None, xyNeg=None, x2Neg=None, y2Neg=None,
397 **kwargs):
398 """Generate dipole model with given parameters.
399
400 This is the function whose sum-of-squared difference from data
401 is minimized by `lmfit`.
402
403 x : TODO: DM-17458
404 Input independent variable. Used here as the grid on
405 which to compute the background gradient model.
406 flux : `float`
407 Desired flux of the positive lobe of the dipole
408 xcenPos, ycenPos : `float`
409 Desired x,y-centroid of the positive lobe of the dipole
410 xcenNeg, ycenNeg : `float`
411 Desired x,y-centroid of the negative lobe of the dipole
412 fluxNeg : `float`, optional
413 Desired flux of the negative lobe of the dipole, set to 'flux' if None
414 b, x1, y1, xy, x2, y2 : `float`
415 Gradient parameters for positive lobe.
416 bNeg, x1Neg, y1Neg, xyNeg, x2Neg, y2Neg : `float`, optional
417 Gradient parameters for negative lobe.
418 They are set to the corresponding positive values if None.
419
420 **kwargs : `dict` [`str`]
421 Keyword arguments passed through ``lmfit`` and
422 used by this function. These must include:
423
424 - ``psf`` Psf model used to generate the 'star'
425 - ``rel_weight`` Used to signify least-squares weighting of posImage/negImage
426 relative to diffim. If ``rel_weight == 0`` then posImage/negImage are ignored.
427 - ``bbox`` Bounding box containing region to be modelled
428
429 Returns
430 -------
431 zout : `numpy.array`
432 Has width and height matching the input bbox, and
433 contains the dipole model with given centroids and flux(es). If
434 ``rel_weight`` = 0, this is a 2-d array with dimensions matching
435 those of bbox; otherwise a stack of three such arrays,
436 representing the dipole (diffim), positive, and negative images
437 respectively.
438 """
439
440 psf = kwargs.get('psf')
441 rel_weight = kwargs.get('rel_weight') # if > 0, we're including pre-sub. images
442 fp = kwargs.get('footprint')
443 bbox = fp.getBBox()
444
445 if fluxNeg is None:
446 fluxNeg = flux
447
448 self.log.debug('flux: %.2f fluxNeg: %.2f x+: %.2f x-: %.2f y+: %.2f y-: %.2f ',
449 flux, fluxNeg, xcenPos, xcenNeg, ycenPos, ycenNeg)
450 if x1 is not None:
451 self.log.debug(' b: %.2f x1: %.2f y1: %.2f', b, x1, y1)
452 if xy is not None:
453 self.log.debug(' xy: %.2f x2: %.2f y2: %.2f', xy, x2, y2)
454
455 posIm = self.makeStarModel(bbox, psf, xcenPos, ycenPos, flux)
456 negIm = self.makeStarModel(bbox, psf, xcenNeg, ycenNeg, fluxNeg)
457
458 in_x = x
459 if in_x is None: # use the footprint to generate the input grid
460 y, x = np.mgrid[bbox.getBeginY():bbox.getEndY(), bbox.getBeginX():bbox.getEndX()]
461 in_x = np.array([x, y]) * 1.
462 in_x[0, :] -= in_x[0, :].mean() # center it!
463 in_x[1, :] -= in_x[1, :].mean()
464
465 if b is not None:
466 gradient = self.makeBackgroundModel(in_x, (b, x1, y1, xy, x2, y2))
467
468 # If bNeg is None, then don't fit the negative background separately
469 if bNeg is not None:
470 gradientNeg = self.makeBackgroundModel(in_x, (bNeg, x1Neg, y1Neg, xyNeg, x2Neg, y2Neg))
471 else:
472 gradientNeg = gradient
473
474 posIm.array[:, :] += gradient
475 negIm.array[:, :] += gradientNeg
476
477 # Generate the diffIm model
478 diffIm = afwImage.ImageF(bbox)
479 diffIm += posIm
480 diffIm -= negIm
481
482 zout = diffIm.array
483 if rel_weight > 0.:
484 zout = np.append([zout], [posIm.array, negIm.array], axis=0)
485
486 return zout
487
488
490 """Fit a dipole model using an image difference.
491
492 See also:
493 `DMTN-007: Dipole characterization for image differencing <https://dmtn-007.lsst.io>`_.
494 """
495
496 # This is just a private version number to sync with the ipython notebooks that I have been
497 # using for algorithm development.
498 _private_version_ = '0.0.5'
499
500 # Below is a (somewhat incomplete) list of improvements
501 # that would be worth investigating, given the time:
502
503 # todo 1. evaluate necessity for separate parameters for pos- and neg- images
504 # todo 2. only fit background OUTSIDE footprint (DONE) and dipole params INSIDE footprint (NOT DONE)?
505 # todo 3. correct normalization of least-squares weights based on variance planes
506 # todo 4. account for PSFs that vary across the exposures (should be happening by default?)
507 # todo 5. correctly account for NA/masks (i.e., ignore!)
508 # todo 6. better exception handling in the plugin
509 # todo 7. better classification of dipoles (e.g. by comparing chi2 fit vs. monopole?)
510 # todo 8. (DONE) Initial fast estimate of background gradient(s) params -- perhaps using numpy.lstsq
511 # todo 9. (NOT NEEDED - see (2)) Initial fast test whether a background gradient needs to be fit
512 # todo 10. (DONE) better initial estimate for flux when there's a strong gradient
513 # todo 11. (DONE) requires a new package `lmfit` -- investiate others? (astropy/scipy/iminuit?)
514
515 def __init__(self, diffim, posImage=None, negImage=None):
516 """Algorithm to run dipole measurement on a diaSource
517
518 Parameters
519 ----------
520 diffim : `lsst.afw.image.Exposure`
521 Exposure on which the diaSources were detected
522 posImage : `lsst.afw.image.Exposure`
523 "Positive" exposure from which the template was subtracted
524 negImage : `lsst.afw.image.Exposure`
525 "Negative" exposure which was subtracted from the posImage
526 """
527
528 self.diffim = diffim
529 self.posImage = posImage
530 self.negImage = negImage
531 self.psfSigma = None
532 if diffim is not None:
533 diffimPsf = diffim.getPsf()
534 diffimAvgPos = diffimPsf.getAveragePosition()
535 self.psfSigma = diffimPsf.computeShape(diffimAvgPos).getDeterminantRadius()
536
537 self.log = logging.getLogger(__name__)
538
539 import lsstDebug
540 self.debug = lsstDebug.Info(__name__).debug
541
542 def fitDipoleImpl(self, source, tol=1e-7, rel_weight=0.5,
543 fitBackground=1, bgGradientOrder=1, maxSepInSigma=5.,
544 separateNegParams=True, verbose=False):
545 """Fit a dipole model to an input difference image.
546
547 Actually, fits the subimage bounded by the input source's
548 footprint) and optionally constrain the fit using the
549 pre-subtraction images posImage and negImage.
550
551 Parameters
552 ----------
553 source : TODO: DM-17458
554 TODO: DM-17458
555 tol : float, optional
556 TODO: DM-17458
557 rel_weight : `float`, optional
558 TODO: DM-17458
559 fitBackground : `int`, optional
560 TODO: DM-17458
561 bgGradientOrder : `int`, optional
562 TODO: DM-17458
563 maxSepInSigma : `float`, optional
564 TODO: DM-17458
565 separateNegParams : `bool`, optional
566 TODO: DM-17458
567 verbose : `bool`, optional
568 TODO: DM-17458
569
570 Returns
571 -------
572 result : `lmfit.MinimizerResult`
573 return `lmfit.MinimizerResult` object containing the fit
574 parameters and other information.
575 """
576
577 # Only import lmfit if someone wants to use the new DipoleFitAlgorithm.
578 import lmfit
579
580 fp = source.getFootprint()
581 bbox = fp.getBBox()
582 subim = afwImage.MaskedImageF(self.diffim.getMaskedImage(), bbox=bbox, origin=afwImage.PARENT)
583
584 z = diArr = subim.image.array
585 # Make sure we don't overwrite buffers.
586 z = z.copy()
587 weights = 1. / subim.variance.array # get the weights (=1/variance)
588
589 if rel_weight > 0. and ((self.posImage is not None) or (self.negImage is not None)):
590 if self.negImage is not None:
591 negSubim = afwImage.MaskedImageF(self.negImage.getMaskedImage(), bbox, origin=afwImage.PARENT)
592 if self.posImage is not None:
593 posSubim = afwImage.MaskedImageF(self.posImage.getMaskedImage(), bbox, origin=afwImage.PARENT)
594 if self.posImage is None: # no science image provided; generate it from diffim + negImage
595 posSubim = subim.clone()
596 posSubim += negSubim
597 if self.negImage is None: # no template provided; generate it from the posImage - diffim
598 negSubim = posSubim.clone()
599 negSubim -= subim
600
601 z = np.append([z], [posSubim.image.array,
602 negSubim.image.array], axis=0)
603 # Weight the pos/neg images by rel_weight relative to the diffim
604 weights = np.append([weights], [1. / posSubim.variance.array * rel_weight,
605 1. / negSubim.variance.array * rel_weight], axis=0)
606 else:
607 rel_weight = 0. # a short-cut for "don't include the pre-subtraction data"
608
609 # It seems that `lmfit` requires a static functor as its optimized method, which eliminates
610 # the ability to pass a bound method or other class method. Here we write a wrapper which
611 # makes this possible.
612 def dipoleModelFunctor(x, flux, xcenPos, ycenPos, xcenNeg, ycenNeg, fluxNeg=None,
613 b=None, x1=None, y1=None, xy=None, x2=None, y2=None,
614 bNeg=None, x1Neg=None, y1Neg=None, xyNeg=None, x2Neg=None, y2Neg=None,
615 **kwargs):
616 """Generate dipole model with given parameters.
617
618 It simply defers to `modelObj.makeModel()`, where `modelObj` comes
619 out of `kwargs['modelObj']`.
620 """
621 modelObj = kwargs.pop('modelObj')
622 return modelObj.makeModel(x, flux, xcenPos, ycenPos, xcenNeg, ycenNeg, fluxNeg=fluxNeg,
623 b=b, x1=x1, y1=y1, xy=xy, x2=x2, y2=y2,
624 bNeg=bNeg, x1Neg=x1Neg, y1Neg=y1Neg, xyNeg=xyNeg,
625 x2Neg=x2Neg, y2Neg=y2Neg, **kwargs)
626
627 dipoleModel = DipoleModel()
628
629 modelFunctor = dipoleModelFunctor # dipoleModel.makeModel does not work for now.
630 # Create the lmfit model (lmfit uses scipy 'leastsq' option by default - Levenberg-Marquardt)
631 # We have to (later) filter out the nans by hand in our input to gmod.fit().
632 # The only independent variable in the model is "x"; lmfit tries to
633 # introspect variables and parameters from the function signature, but
634 # gets it wrong for the model signature above.
635 gmod = lmfit.Model(modelFunctor, independent_vars=["x"], verbose=verbose)
636
637 # Add the constraints for centroids, fluxes.
638 # starting constraint - near centroid of footprint
639 fpCentroid = np.array([fp.getCentroid().getX(), fp.getCentroid().getY()])
640 cenNeg = cenPos = fpCentroid
641
642 pks = fp.getPeaks()
643
644 if len(pks) >= 1:
645 cenPos = pks[0].getF() # if individual (merged) peaks were detected, use those
646 if len(pks) >= 2: # peaks are already sorted by centroid flux so take the most negative one
647 cenNeg = pks[-1].getF()
648
649 # For close/faint dipoles the starting locs (min/max) might be way off, let's help them a bit.
650 # First assume dipole is not separated by more than 5*psfSigma.
651 maxSep = self.psfSigma * maxSepInSigma
652
653 # As an initial guess -- assume the dipole is close to the center of the footprint.
654 if np.sum(np.sqrt((np.array(cenPos) - fpCentroid)**2.)) > maxSep:
655 cenPos = fpCentroid
656 if np.sum(np.sqrt((np.array(cenNeg) - fpCentroid)**2.)) > maxSep:
657 cenPos = fpCentroid
658
659 # parameter hints/constraints: https://lmfit.github.io/lmfit-py/model.html#model-param-hints-section
660 # might make sense to not use bounds -- see http://lmfit.github.io/lmfit-py/bounds.html
661 # also see this discussion -- https://github.com/scipy/scipy/issues/3129
662 gmod.set_param_hint('xcenPos', value=cenPos[0],
663 min=cenPos[0]-maxSep, max=cenPos[0]+maxSep)
664 gmod.set_param_hint('ycenPos', value=cenPos[1],
665 min=cenPos[1]-maxSep, max=cenPos[1]+maxSep)
666 gmod.set_param_hint('xcenNeg', value=cenNeg[0],
667 min=cenNeg[0]-maxSep, max=cenNeg[0]+maxSep)
668 gmod.set_param_hint('ycenNeg', value=cenNeg[1],
669 min=cenNeg[1]-maxSep, max=cenNeg[1]+maxSep)
670
671 # Use the (flux under the dipole)*5 for an estimate.
672 # Lots of testing showed that having startingFlux be too high was better than too low.
673 startingFlux = np.nansum(np.abs(diArr) - np.nanmedian(np.abs(diArr))) * 5.
674 posFlux = negFlux = startingFlux
675
676 # TBD: set max. flux limit?
677 gmod.set_param_hint('flux', value=posFlux, min=0.1)
678
679 if separateNegParams:
680 # TBD: set max negative lobe flux limit?
681 gmod.set_param_hint('fluxNeg', value=np.abs(negFlux), min=0.1)
682
683 # Fixed parameters (don't fit for them if there are no pre-sub images or no gradient fit requested):
684 # Right now (fitBackground == 1), we fit a linear model to the background and then subtract
685 # it from the data and then don't fit the background again (this is faster).
686 # A slower alternative (fitBackground == 2) is to use the estimated background parameters as
687 # starting points in the integrated model fit. That is currently not performed by default,
688 # but might be desirable in some cases.
689 bgParsPos = bgParsNeg = (0., 0., 0.)
690 if ((rel_weight > 0.) and (fitBackground != 0) and (bgGradientOrder >= 0)):
691 pbg = 0.
692 bgFitImage = self.posImage if self.posImage is not None else self.negImage
693 # Fit the gradient to the background (linear model)
694 bgParsPos = bgParsNeg = dipoleModel.fitFootprintBackground(source, bgFitImage,
695 order=bgGradientOrder)
696
697 # Generate the gradient and subtract it from the pre-subtraction image data
698 if fitBackground == 1:
699 in_x = dipoleModel._generateXYGrid(bbox)
700 pbg = dipoleModel.makeBackgroundModel(in_x, tuple(bgParsPos))
701 z[1, :] -= pbg
702 z[1, :] -= np.nanmedian(z[1, :])
703 posFlux = np.nansum(z[1, :])
704 gmod.set_param_hint('flux', value=posFlux*1.5, min=0.1)
705
706 if separateNegParams and self.negImage is not None:
707 bgParsNeg = dipoleModel.fitFootprintBackground(source, self.negImage,
708 order=bgGradientOrder)
709 pbg = dipoleModel.makeBackgroundModel(in_x, tuple(bgParsNeg))
710 z[2, :] -= pbg
711 z[2, :] -= np.nanmedian(z[2, :])
712 if separateNegParams:
713 negFlux = np.nansum(z[2, :])
714 gmod.set_param_hint('fluxNeg', value=negFlux*1.5, min=0.1)
715
716 # Do not subtract the background from the images but include the background parameters in the fit
717 if fitBackground == 2:
718 if bgGradientOrder >= 0:
719 gmod.set_param_hint('b', value=bgParsPos[0])
720 if separateNegParams:
721 gmod.set_param_hint('bNeg', value=bgParsNeg[0])
722 if bgGradientOrder >= 1:
723 gmod.set_param_hint('x1', value=bgParsPos[1])
724 gmod.set_param_hint('y1', value=bgParsPos[2])
725 if separateNegParams:
726 gmod.set_param_hint('x1Neg', value=bgParsNeg[1])
727 gmod.set_param_hint('y1Neg', value=bgParsNeg[2])
728 if bgGradientOrder >= 2:
729 gmod.set_param_hint('xy', value=bgParsPos[3])
730 gmod.set_param_hint('x2', value=bgParsPos[4])
731 gmod.set_param_hint('y2', value=bgParsPos[5])
732 if separateNegParams:
733 gmod.set_param_hint('xyNeg', value=bgParsNeg[3])
734 gmod.set_param_hint('x2Neg', value=bgParsNeg[4])
735 gmod.set_param_hint('y2Neg', value=bgParsNeg[5])
736
737 y, x = np.mgrid[bbox.getBeginY():bbox.getEndY(), bbox.getBeginX():bbox.getEndX()]
738 in_x = np.array([x, y]).astype(np.float64)
739 in_x[0, :] -= in_x[0, :].mean() # center it!
740 in_x[1, :] -= in_x[1, :].mean()
741
742 # Instead of explicitly using a mask to ignore flagged pixels, just set the ignored pixels'
743 # weights to 0 in the fit. TBD: need to inspect mask planes to set this mask.
744 mask = np.ones_like(z, dtype=bool) # TBD: set mask values to False if the pixels are to be ignored
745
746 # I'm not sure about the variance planes in the diffim (or convolved pre-sub. images
747 # for that matter) so for now, let's just do an un-weighted least-squares fit
748 # (override weights computed above).
749 weights = mask.astype(np.float64)
750 if self.posImage is not None and rel_weight > 0.:
751 weights = np.array([np.ones_like(diArr), np.ones_like(diArr)*rel_weight,
752 np.ones_like(diArr)*rel_weight])
753
754 # Set the weights to zero if mask is False
755 if np.any(~mask):
756 weights[~mask] = 0.
757
758 # Filter out any nans, and make the weights 0.
759 nans = (np.isnan(z) | np.isnan(weights))
760 nNans = nans.sum()
761 if nNans > 0:
762 if nNans < len(z):
763 z[nans] = np.nanmedian(z)
764 else:
765 z[nans] = 0
766 weights[nans] = 0
767
768 # Note that although we can, we're not required to set initial values for params here,
769 # since we set their param_hint's above.
770 # Can add "method" param to not use 'leastsq' (==levenberg-marquardt), e.g. "method='nelder'"
771 with warnings.catch_warnings():
772 # Ignore lmfit unknown argument warnings:
773 # "psf, rel_weight, footprint, modelObj" all become pass-through kwargs for makeModel.
774 warnings.filterwarnings("ignore", "The keyword argument .* does not match", UserWarning)
775 result = gmod.fit(z, weights=weights, x=in_x, max_nfev=250,
776 method="leastsq", # TODO: try using `least_squares` here for speed/robustness
777 verbose=verbose,
778 # see scipy docs for the meaning of these keywords
779 fit_kws={'ftol': tol, 'xtol': tol, 'gtol': tol,
780 # Our model is float32 internally, so we need a larger epsfcn.
781 'epsfcn': 1e-8},
782 psf=self.diffim.getPsf(), # hereon: kwargs that get passed to makeModel()
783 rel_weight=rel_weight,
784 footprint=fp,
785 modelObj=dipoleModel)
786
787 if verbose: # the ci_report() seems to fail if neg params are constrained -- TBD why.
788 # Never wanted in production - this takes a long time (longer than the fit!)
789 # This is how to get confidence intervals out:
790 # https://lmfit.github.io/lmfit-py/confidence.html and
791 # http://cars9.uchicago.edu/software/python/lmfit/model.html
792 print(result.fit_report(show_correl=False))
793 if separateNegParams:
794 print(result.ci_report())
795
796 return result
797
798 def fitDipole(self, source, tol=1e-7, rel_weight=0.1,
799 fitBackground=1, maxSepInSigma=5., separateNegParams=True,
800 bgGradientOrder=1, verbose=False, display=False):
801 """Fit a dipole model to an input ``diaSource`` (wraps `fitDipoleImpl`).
802
803 Actually, fits the subimage bounded by the input source's
804 footprint) and optionally constrain the fit using the
805 pre-subtraction images self.posImage (science) and
806 self.negImage (template). Wraps the output into a
807 `pipeBase.Struct` named tuple after computing additional
808 statistics such as orientation and SNR.
809
810 Parameters
811 ----------
812 source : `lsst.afw.table.SourceRecord`
813 Record containing the (merged) dipole source footprint detected on the diffim
814 tol : `float`, optional
815 Tolerance parameter for scipy.leastsq() optimization
816 rel_weight : `float`, optional
817 Weighting of posImage/negImage relative to the diffim in the fit
818 fitBackground : `int`, {0, 1, 2}, optional
819 How to fit linear background gradient in posImage/negImage
820
821 - 0: do not fit background at all
822 - 1 (default): pre-fit the background using linear least squares and then do not fit it
823 as part of the dipole fitting optimization
824 - 2: pre-fit the background using linear least squares (as in 1), and use the parameter
825 estimates from that fit as starting parameters for an integrated "re-fit" of the
826 background as part of the overall dipole fitting optimization.
827 maxSepInSigma : `float`, optional
828 Allowed window of centroid parameters relative to peak in input source footprint
829 separateNegParams : `bool`, optional
830 Fit separate parameters to the flux and background gradient in
831 bgGradientOrder : `int`, {0, 1, 2}, optional
832 Desired polynomial order of background gradient
833 verbose: `bool`, optional
834 Be verbose
835 display
836 Display input data, best fit model(s) and residuals in a matplotlib window.
837
838 Returns
839 -------
840 result : `struct`
841 `pipeBase.Struct` object containing the fit parameters and other information.
842
843 result : `callable`
844 `lmfit.MinimizerResult` object for debugging and error estimation, etc.
845
846 Notes
847 -----
848 Parameter `fitBackground` has three options, thus it is an integer:
849
850 """
851
852 fitResult = self.fitDipoleImpl(
853 source, tol=tol, rel_weight=rel_weight, fitBackground=fitBackground,
854 maxSepInSigma=maxSepInSigma, separateNegParams=separateNegParams,
855 bgGradientOrder=bgGradientOrder, verbose=verbose)
856
857 # Display images, model fits and residuals (currently uses matplotlib display functions)
858 if display:
859 fp = source.getFootprint()
860 self.displayFitResults(fp, fitResult)
861
862 fitParams = fitResult.best_values
863 if fitParams['flux'] <= 1.: # usually around 0.1 -- the minimum flux allowed -- i.e. bad fit.
864 return None, fitResult
865
866 centroid = ((fitParams['xcenPos'] + fitParams['xcenNeg']) / 2.,
867 (fitParams['ycenPos'] + fitParams['ycenNeg']) / 2.)
868 dx, dy = fitParams['xcenPos'] - fitParams['xcenNeg'], fitParams['ycenPos'] - fitParams['ycenNeg']
869 angle = np.arctan2(dy, dx)
870
871 # Exctract flux value, compute signalToNoise from flux/variance_within_footprint
872 # Also extract the stderr of flux estimate.
873 def computeSumVariance(exposure, footprint):
874 return np.sqrt(np.nansum(exposure[footprint.getBBox(), afwImage.PARENT].variance.array))
875
876 fluxVal = fluxVar = fitParams['flux']
877 fluxErr = fluxErrNeg = fitResult.params['flux'].stderr
878 if self.posImage is not None:
879 fluxVar = computeSumVariance(self.posImage, source.getFootprint())
880 else:
881 fluxVar = computeSumVariance(self.diffim, source.getFootprint())
882
883 fluxValNeg, fluxVarNeg = fluxVal, fluxVar
884 if separateNegParams:
885 fluxValNeg = fitParams['fluxNeg']
886 fluxErrNeg = fitResult.params['fluxNeg'].stderr
887 if self.negImage is not None:
888 fluxVarNeg = computeSumVariance(self.negImage, source.getFootprint())
889
890 try:
891 signalToNoise = np.sqrt((fluxVal/fluxVar)**2 + (fluxValNeg/fluxVarNeg)**2)
892 except ZeroDivisionError: # catch divide by zero - should never happen.
893 signalToNoise = np.nan
894
895 out = Struct(posCentroidX=fitParams['xcenPos'], posCentroidY=fitParams['ycenPos'],
896 negCentroidX=fitParams['xcenNeg'], negCentroidY=fitParams['ycenNeg'],
897 posFlux=fluxVal, negFlux=-fluxValNeg, posFluxErr=fluxErr, negFluxErr=fluxErrNeg,
898 centroidX=centroid[0], centroidY=centroid[1], orientation=angle,
899 signalToNoise=signalToNoise, chi2=fitResult.chisqr, redChi2=fitResult.redchi,
900 nData=fitResult.ndata)
901
902 # fitResult may be returned for debugging
903 return out, fitResult
904
905 def displayFitResults(self, footprint, result):
906 """Display data, model fits and residuals (currently uses matplotlib display functions).
907
908 Parameters
909 ----------
910 footprint : `lsst.afw.detection.Footprint`
911 Footprint containing the dipole that was fit
912 result : `lmfit.MinimizerResult`
913 `lmfit.MinimizerResult` object returned by `lmfit` optimizer
914 """
915 try:
916 import matplotlib.pyplot as plt
917 except ImportError as err:
918 self.log.warning('Unable to import matplotlib: %s', err)
919 raise err
920
921 def display2dArray(arr, title='Data', extent=None):
922 """Use `matplotlib.pyplot.imshow` to display a 2-D array with a given coordinate range.
923 """
924 fig = plt.imshow(arr, origin='lower', interpolation='none', cmap='gray', extent=extent)
925 plt.title(title)
926 plt.colorbar(fig, cmap='gray')
927 return fig
928
929 z = result.data
930 fit = result.best_fit
931 bbox = footprint.getBBox()
932 extent = (bbox.getBeginX(), bbox.getEndX(), bbox.getBeginY(), bbox.getEndY())
933 if z.shape[0] == 3:
934 plt.figure(figsize=(8, 8))
935 for i in range(3):
936 plt.subplot(3, 3, i*3+1)
937 display2dArray(z[i, :], 'Data', extent=extent)
938 plt.subplot(3, 3, i*3+2)
939 display2dArray(fit[i, :], 'Model', extent=extent)
940 plt.subplot(3, 3, i*3+3)
941 display2dArray(z[i, :] - fit[i, :], 'Residual', extent=extent)
942 else:
943 plt.figure(figsize=(8, 2.5))
944 plt.subplot(1, 3, 1)
945 display2dArray(z, 'Data', extent=extent)
946 plt.subplot(1, 3, 2)
947 display2dArray(fit, 'Model', extent=extent)
948 plt.subplot(1, 3, 3)
949 display2dArray(z - fit, 'Residual', extent=extent)
950
951 plt.show()
952
953
954@measBase.register("ip_diffim_DipoleFit")
955class DipoleFitPlugin(measBase.SingleFramePlugin):
956 """A single frame measurement plugin that fits dipoles to all merged (two-peak) ``diaSources``.
957
958 This measurement plugin accepts up to three input images in
959 its `measure` method. If these are provided, it includes data
960 from the pre-subtraction posImage (science image) and optionally
961 negImage (template image) to constrain the fit. The meat of the
962 fitting routines are in the class `~lsst.module.name.DipoleFitAlgorithm`.
963
964 Notes
965 -----
966 The motivation behind this plugin and the necessity for including more than
967 one exposure are documented in DMTN-007 (http://dmtn-007.lsst.io).
968
969 This class is named `ip_diffim_DipoleFit` so that it may be used alongside
970 the existing `ip_diffim_DipoleMeasurement` classes until such a time as those
971 are deemed to be replaceable by this.
972 """
973
974 ConfigClass = DipoleFitPluginConfig
975 DipoleFitAlgorithmClass = DipoleFitAlgorithm # Pointer to the class that performs the fit
976
977 FAILURE_EDGE = 1 # too close to the edge
978 FAILURE_FIT = 2 # failure in the fitting
979 FAILURE_NOT_DIPOLE = 4 # input source is not a putative dipole to begin with
980 FAILURE_TOO_LARGE = 8 # input source is too large to be fit
981
982 @classmethod
984 """This algorithm simultaneously fits the centroid and flux, and does
985 not require any previous centroid fit.
986 """
987 return cls.CENTROID_ORDER
988
989 def __init__(self, config, name, schema, metadata, logName=None):
990 if logName is None:
991 logName = name
992 measBase.SingleFramePlugin.__init__(self, config, name, schema, metadata, logName=logName)
993
994 self.log = logging.getLogger(logName)
995
996 self._setupSchema(config, name, schema, metadata)
997
998 def _setupSchema(self, config, name, schema, metadata):
999 """Add fields for the outputs, and save the keys for fast assignment.
1000 """
1001 self.posFluxKey = measBase.FluxResultKey.addFields(schema,
1002 schema.join(name, "pos"),
1003 "Dipole positive lobe instrumental flux.")
1004 self.negFluxKey = measBase.FluxResultKey.addFields(schema,
1005 schema.join(name, "neg"),
1006 "Dipole negative lobe instrumental flux.")
1007 doc = "Dipole overall instrumental flux (mean of absolute value of positive and negative lobes)."
1008 self.fluxKey = measBase.FluxResultKey.addFields(schema, name, doc)
1009
1010 self.posCentroidKey = measBase.CentroidResultKey.addFields(schema,
1011 schema.join(name, "pos"),
1012 "Dipole positive lobe centroid position.",
1013 measBase.UncertaintyEnum.NO_UNCERTAINTY)
1014 self.negCentroidKey = measBase.CentroidResultKey.addFields(schema,
1015 schema.join(name, "neg"),
1016 "Dipole negative lobe centroid position.",
1017 measBase.UncertaintyEnum.NO_UNCERTAINTY)
1018 self.centroidKey = measBase.CentroidResultKey.addFields(schema,
1019 name,
1020 "Dipole centroid position.",
1021 measBase.UncertaintyEnum.SIGMA_ONLY)
1022
1023 self.orientationKey = schema.addField(
1024 schema.join(name, "orientation"), type=float, units="rad",
1025 doc="Dipole orientation. Convention is CCW from +x on image.")
1026
1027 self.separationKey = schema.addField(
1028 schema.join(name, "separation"), type=float, units="pixel",
1029 doc="Pixel separation between positive and negative lobes of dipole")
1030
1031 self.chi2dofKey = schema.addField(
1032 schema.join(name, "chi2dof"), type=float,
1033 doc="Chi2 per degree of freedom (chi2/(nData-nVariables)) of dipole fit")
1034
1035 self.nDataKey = schema.addField(
1036 schema.join(name, "nData"), type=np.int64,
1037 doc="Number of data points in the dipole fit")
1038
1039 self.signalToNoiseKey = schema.addField(
1040 schema.join(name, "signalToNoise"), type=float,
1041 doc="Estimated signal-to-noise of dipole fit")
1042
1043 self.classificationFlagKey = schema.addField(
1044 schema.join(name, "flag", "classification"), type="Flag",
1045 doc="Flag indicating diaSource is classified as a dipole")
1046
1047 self.classificationAttemptedFlagKey = schema.addField(
1048 schema.join(name, "flag", "classificationAttempted"), type="Flag",
1049 doc="Flag indicating diaSource was attempted to be classified as a dipole")
1050
1051 self.flagKey = schema.addField(
1052 schema.join(name, "flag"), type="Flag",
1053 doc="General failure flag for dipole fit")
1054
1055 self.edgeFlagKey = schema.addField(
1056 schema.join(name, "flag", "edge"), type="Flag",
1057 doc="Flag set when dipole is too close to edge of image")
1058
1059 def measureDipoles(self, measRecord, exposure, posExp=None, negExp=None):
1060 """Perform the non-linear least squares minimization on the putative dipole source.
1061
1062 Parameters
1063 ----------
1064 measRecord : `lsst.afw.table.SourceRecord`
1065 diaSources that will be measured using dipole measurement
1066 exposure : `lsst.afw.image.Exposure`
1067 Difference exposure on which the diaSources were detected; `exposure = posExp-negExp`
1068 If both `posExp` and `negExp` are `None`, will attempt to fit the
1069 dipole to just the `exposure` with no constraint.
1070 posExp : `lsst.afw.image.Exposure`, optional
1071 "Positive" exposure, typically a science exposure, or None if unavailable
1072 When `posExp` is `None`, will compute `posImage = exposure + negExp`.
1073 negExp : `lsst.afw.image.Exposure`, optional
1074 "Negative" exposure, typically a template exposure, or None if unavailable
1075 When `negExp` is `None`, will compute `negImage = posExp - exposure`.
1076
1077 Notes
1078 -----
1079 The main functionality of this routine was placed outside of
1080 this plugin (into `DipoleFitAlgorithm.fitDipole()`) so that
1081 `DipoleFitAlgorithm.fitDipole()` can be called separately for
1082 testing (@see `tests/testDipoleFitter.py`)
1083 """
1084 result = None
1085 pks = measRecord.getFootprint().getPeaks()
1086
1087 # Check if the footprint consists of a putative dipole - else don't fit it.
1088 if (
1089 # One peak in the footprint (not a dipole)
1090 ((nPeaks := len(pks)) <= 1)
1091 # Peaks are the same sign (not a dipole); peaks are ordered
1092 # from highest to lowest.
1093 or (nPeaks > 1 and (np.sign(pks[0].getPeakValue())
1094 == np.sign(pks[-1].getPeakValue())))
1095 ):
1096 if not self.config.fitAllDiaSources:
1097 # Non-dipoles fall back on the centroid slot for positions,
1098 # errors, and the failure flag, if we're not fitting them.
1099 measRecord[self.centroidKey.getX()] = measRecord.getX()
1100 measRecord[self.centroidKey.getY()] = measRecord.getY()
1101 self.centroidKey.getCentroidErr().set(measRecord, measRecord.getCentroidErr())
1102 measRecord[self.flagKey] = measRecord.getCentroidFlag()
1103 return
1104
1105 # Footprint is too large (not a dipole).
1106 if ((area := measRecord.getFootprint().getArea()) > self.config.maxFootprintArea):
1107 self.fail(measRecord, measBase.MeasurementError(f"{area} > {self.config.maxFootprintArea}",
1109 return
1110
1111 try:
1112 alg = self.DipoleFitAlgorithmClass(exposure, posImage=posExp, negImage=negExp)
1113 result, _ = alg.fitDipole(
1114 measRecord, rel_weight=self.config.relWeight,
1115 tol=self.config.tolerance,
1116 maxSepInSigma=self.config.maxSeparation,
1117 fitBackground=self.config.fitBackground,
1118 separateNegParams=self.config.fitSeparateNegParams,
1119 verbose=False, display=False)
1120 except pexExcept.LengthError:
1121 self.fail(measRecord, measBase.MeasurementError('edge failure', self.FAILURE_EDGEFAILURE_EDGE))
1122 except Exception as e:
1123 errorMessage = f"Exception in dipole fit. {e.__class__.__name__}: {e}"
1124 self.fail(measRecord, measBase.MeasurementError(errorMessage, self.FAILURE_FITFAILURE_FIT))
1125
1126 self.log.debug("Dipole fit result: %d %s", measRecord.getId(), str(result))
1127
1128 if result is None:
1129 self.fail(measRecord, measBase.MeasurementError("bad dipole fit", self.FAILURE_FITFAILURE_FIT))
1130 return
1131
1132 # add chi2, coord/flux uncertainties (TBD), dipole classification
1133 # Add the relevant values to the measRecord
1134 measRecord[self.posFluxKey.getInstFlux()] = result.posFlux
1135 measRecord[self.posFluxKey.getInstFluxErr()] = result.signalToNoise # to be changed to actual sigma!
1136 measRecord[self.posCentroidKey.getX()] = result.posCentroidX
1137 measRecord[self.posCentroidKey.getY()] = result.posCentroidY
1138
1139 measRecord[self.negFluxKey.getInstFlux()] = result.negFlux
1140 measRecord[self.negFluxKey.getInstFluxErr()] = result.signalToNoise # to be changed to actual sigma!
1141 measRecord[self.negCentroidKey.getX()] = result.negCentroidX
1142 measRecord[self.negCentroidKey.getY()] = result.negCentroidY
1143
1144 # Dia source flux: average of pos+neg
1145 measRecord[self.fluxKey.getInstFlux()] = (abs(result.posFlux) + abs(result.negFlux))/2.
1146 measRecord[self.orientationKey] = result.orientation
1147 measRecord[self.separationKey] = np.sqrt((result.posCentroidX - result.negCentroidX)**2.
1148 + (result.posCentroidY - result.negCentroidY)**2.)
1149 measRecord[self.centroidKey.getX()] = result.centroidX
1150 measRecord[self.centroidKey.getY()] = result.centroidY
1151
1152 measRecord[self.signalToNoiseKey] = result.signalToNoise
1153 measRecord[self.chi2dofKey] = result.redChi2
1154
1155 if result.nData >= 1:
1156 measRecord[self.nDataKey] = result.nData
1157 else:
1158 measRecord[self.nDataKey] = 0
1159
1160 self.doClassify(measRecord, result.chi2)
1161
1162 def doClassify(self, measRecord, chi2val):
1163 """Classify a source as a dipole.
1164
1165 Parameters
1166 ----------
1167 measRecord : TODO: DM-17458
1168 TODO: DM-17458
1169 chi2val : TODO: DM-17458
1170 TODO: DM-17458
1171
1172 Notes
1173 -----
1174 Sources are classified as dipoles, or not, according to three criteria:
1175
1176 1. Does the total signal-to-noise surpass the ``minSn``?
1177 2. Are the pos/neg fluxes greater than 1.0 and no more than 0.65 (``maxFluxRatio``)
1178 of the total flux? By default this will never happen since ``posFlux == negFlux``.
1179 3. Is it a good fit (``chi2dof`` < 1)? (Currently not used.)
1180 """
1181
1182 # First, does the total signal-to-noise surpass the minSn?
1183 passesSn = measRecord[self.signalToNoiseKey] > self.config.minSn
1184
1185 # Second, are the pos/neg fluxes greater than 1.0 and no more than 0.65 (param maxFluxRatio)
1186 # of the total flux? By default this will never happen since posFlux = negFlux.
1187 passesFluxPos = (abs(measRecord[self.posFluxKey.getInstFlux()])
1188 / (measRecord[self.fluxKey.getInstFlux()]*2.)) < self.config.maxFluxRatio
1189 passesFluxPos &= (abs(measRecord[self.posFluxKey.getInstFlux()]) >= 1.0)
1190 passesFluxNeg = (abs(measRecord[self.negFluxKey.getInstFlux()])
1191 / (measRecord[self.fluxKey.getInstFlux()]*2.)) < self.config.maxFluxRatio
1192 passesFluxNeg &= (abs(measRecord[self.negFluxKey.getInstFlux()]) >= 1.0)
1193 allPass = (passesSn and passesFluxPos and passesFluxNeg) # and passesChi2)
1194
1195 # Third, is it a good fit (chi2dof < 1)?
1196 # Use scipy's chi2 cumulative distrib to estimate significance
1197 # This doesn't really work since I don't trust the values in the variance plane (which
1198 # affects the least-sq weights, which affects the resulting chi2).
1199 # But I'm going to keep this here for future use.
1200 if False:
1201 from scipy.stats import chi2
1202 ndof = chi2val / measRecord[self.chi2dofKey]
1203 significance = chi2.cdf(chi2val, ndof)
1204 passesChi2 = significance < self.config.maxChi2DoF
1205 allPass = allPass and passesChi2
1206
1207 measRecord.set(self.classificationAttemptedFlagKey, True)
1208
1209 if allPass: # Note cannot pass `allPass` into the `measRecord.set()` call below...?
1210 measRecord.set(self.classificationFlagKey, True)
1211 else:
1212 measRecord.set(self.classificationFlagKey, False)
1213
1214 def fail(self, measRecord, error=None):
1215 """Catch failures and set the correct flags.
1216
1217 Fallback on the current slot centroid positions, but set the dipole
1218 failure flag, since we attempted to fit the source.
1219 """
1220 measRecord[self.centroidKey.getX()] = measRecord.getX()
1221 measRecord[self.centroidKey.getY()] = measRecord.getY()
1222 self.centroidKey.getCentroidErr().set(measRecord, measRecord.getCentroidErr())
1223
1224 measRecord.set(self.flagKey, True)
1225 if error is not None:
1226 if error.getFlagBit() == self.FAILURE_EDGEFAILURE_EDGE:
1227 self.log.debug('DipoleFitPlugin not run on record %d: %s', measRecord.getId(), str(error))
1228 measRecord.set(self.edgeFlagKey, True)
1229 if error.getFlagBit() == self.FAILURE_FITFAILURE_FIT:
1230 self.log.warning('DipoleFitPlugin failed on record %d: %s', measRecord.getId(), str(error))
1231 if error.getFlagBit() == self.FAILURE_TOO_LARGEFAILURE_TOO_LARGE:
1232 self.log.debug('DipoleFitPlugin not run on record with too large footprint %d: %s',
1233 measRecord.getId(), str(error))
1234 else:
1235 self.log.warning('DipoleFitPlugin failed on record %d', measRecord.getId())
fitDipole(self, source, tol=1e-7, rel_weight=0.1, fitBackground=1, maxSepInSigma=5., separateNegParams=True, bgGradientOrder=1, verbose=False, display=False)
fitDipoleImpl(self, source, tol=1e-7, rel_weight=0.5, fitBackground=1, bgGradientOrder=1, maxSepInSigma=5., separateNegParams=True, verbose=False)
__init__(self, diffim, posImage=None, negImage=None)
measureDipoles(self, measRecord, exposure, posExp=None, negExp=None)
_setupSchema(self, config, name, schema, metadata)
__init__(self, config, name, schema, metadata, logName=None)
__init__(self, schema, algMetadata=None, **kwargs)
run(self, sources, exposure, posExp=None, negExp=None, **kwargs)
makeStarModel(self, bbox, psf, xcen, ycen, flux)
_getHeavyFootprintSubimage(self, fp, badfill=np.nan, grow=0)
makeModel(self, x, flux, xcenPos, ycenPos, xcenNeg, ycenNeg, fluxNeg=None, b=None, x1=None, y1=None, xy=None, x2=None, y2=None, bNeg=None, x1Neg=None, y1Neg=None, xyNeg=None, x2Neg=None, y2Neg=None, **kwargs)
fitFootprintBackground(self, source, posImage, order=1)
Reports attempts to exceed implementation-defined length limits for some classes.
Definition Runtime.h:76