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