lsst.ip.diffim  13.0-17-gc2cefa3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Groups Pages
dipoleFitTask.py
Go to the documentation of this file.
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 #
22 from __future__ import absolute_import, division, print_function
23 
24 import numpy as np
25 import warnings
26 
27 from builtins import str
28 from builtins import range
29 from builtins import object
30 
31 import lsst.afw.geom as afwGeom
32 import lsst.afw.image as afwImage
33 import lsst.meas.base as measBase
34 import lsst.afw.table as afwTable
35 import lsst.afw.detection as afwDet
36 from lsst.log import Log
37 import lsst.pex.exceptions as pexExcept
38 import lsst.pex.config as pexConfig
39 from lsst.pipe.base import Struct
40 
41 __all__ = ("DipoleFitTask", "DipoleFitPlugin", "DipoleFitTaskConfig", "DipoleFitPluginConfig")
42 
43 
44 # Create a new measurement task (`DipoleFitTask`) that can handle all other SFM tasks but can
45 # pass a separate pos- and neg- exposure/image to the `DipoleFitPlugin`s `run()` method.
46 
47 
48 class DipoleFitPluginConfig(measBase.SingleFramePluginConfig):
49  """!Configuration for DipoleFitPlugin
50  """
51 
52  fitAllDiaSources = pexConfig.Field(
53  dtype=float, default=False,
54  doc="""Attempte dipole fit of all diaSources (otherwise just the ones consisting of overlapping
55  positive and negative footprints)""")
56 
57  maxSeparation = pexConfig.Field(
58  dtype=float, default=5.,
59  doc="Assume dipole is not separated by more than maxSeparation * psfSigma")
60 
61  relWeight = pexConfig.Field(
62  dtype=float, default=0.5,
63  doc="""Relative weighting of pre-subtraction images (higher -> greater influence of pre-sub.
64  images on fit)""")
65 
66  tolerance = pexConfig.Field(
67  dtype=float, default=1e-7,
68  doc="Fit tolerance")
69 
70  fitBackground = pexConfig.Field(
71  dtype=int, default=1,
72  doc="""Set whether and how to fit for linear gradient in pre-sub. images. Possible values:
73  0: do not fit background at all
74  1 (default): pre-fit the background using linear least squares and then do not fit it as part
75  of the dipole fitting optimization
76  2: pre-fit the background using linear least squares (as in 1), and use the parameter
77  estimates from that fit as starting parameters for an integrated "re-fit" of the background
78  """)
79 
80  fitSeparateNegParams = pexConfig.Field(
81  dtype=bool, default=False,
82  doc="Include parameters to fit for negative values (flux, gradient) separately from pos.")
83 
84  # Config params for classification of detected diaSources as dipole or not
85  minSn = pexConfig.Field(
86  dtype=float, default=np.sqrt(2) * 5.0,
87  doc="Minimum quadrature sum of positive+negative lobe S/N to be considered a dipole")
88 
89  maxFluxRatio = pexConfig.Field(
90  dtype=float, default=0.65,
91  doc="Maximum flux ratio in either lobe to be considered a dipole")
92 
93  maxChi2DoF = pexConfig.Field(
94  dtype=float, default=0.05,
95  doc="""Maximum Chi2/DoF significance of fit to be considered a dipole.
96  Default value means \"Choose a chi2DoF corresponding to a significance level of at most 0.05\"
97  (note this is actually a significance, not a chi2 value).""")
98 
99 
100 class DipoleFitTaskConfig(measBase.SingleFrameMeasurementConfig):
101  """!Measurement of detected diaSources as dipoles
102 
103  Currently we keep the "old" DipoleMeasurement algorithms turned on.
104  """
105 
106  def setDefaults(self):
107  measBase.SingleFrameMeasurementConfig.setDefaults(self)
108 
109  self.plugins.names = ["base_CircularApertureFlux",
110  "base_PixelFlags",
111  "base_SkyCoord",
112  "base_PsfFlux",
113  "base_SdssCentroid",
114  "base_SdssShape",
115  "base_GaussianCentroid",
116  "base_GaussianFlux",
117  "base_PeakLikelihoodFlux",
118  "base_PeakCentroid",
119  "base_NaiveCentroid",
120  "ip_diffim_NaiveDipoleCentroid",
121  "ip_diffim_NaiveDipoleFlux",
122  "ip_diffim_PsfDipoleFlux",
123  "ip_diffim_ClassificationDipole",
124  ]
125 
126  self.slots.calibFlux = None
127  self.slots.modelFlux = None
128  self.slots.instFlux = None
129  self.slots.shape = "base_SdssShape"
130  self.slots.centroid = "ip_diffim_NaiveDipoleCentroid"
131  self.doReplaceWithNoise = False
132 
133 
134 class DipoleFitTask(measBase.SingleFrameMeasurementTask):
135  """!Subclass of SingleFrameMeasurementTask which accepts up to three input images in its run() method.
136 
137  Because it subclasses SingleFrameMeasurementTask, and calls
138  SingleFrameMeasurementTask.run() from its run() method, it still
139  can be used identically to a standard SingleFrameMeasurementTask.
140  """
141 
142  ConfigClass = DipoleFitTaskConfig
143  _DefaultName = "ip_diffim_DipoleFit"
144 
145  def __init__(self, schema, algMetadata=None, **kwds):
146 
147  measBase.SingleFrameMeasurementTask.__init__(self, schema, algMetadata, **kwds)
148 
149  dpFitPluginConfig = self.config.plugins['ip_diffim_DipoleFit']
150 
151  self.dipoleFitter = DipoleFitPlugin(dpFitPluginConfig, name=self._DefaultName,
152  schema=schema, metadata=algMetadata)
153 
154  def run(self, sources, exposure, posExp=None, negExp=None, **kwds):
155  """!Run dipole measurement and classification
156 
157  @param sources diaSources that will be measured using dipole measurement
158  @param exposure Difference exposure on which the diaSources were detected; exposure = posExp - negExp
159  @param posExp "Positive" exposure, typically a science exposure, or None if unavailable
160  @param negExp "Negative" exposure, typically a template exposure, or None if unavailable
161  @param **kwds Sent to SingleFrameMeasurementTask
162 
163  @note When `posExp` is `None`, will compute `posImage = exposure + negExp`.
164  Likewise, when `negExp` is `None`, will compute `negImage = posExp - exposure`.
165  If both `posExp` and `negExp` are `None`, will attempt to fit the dipole to just the `exposure`
166  with no constraint.
167  """
168 
169  measBase.SingleFrameMeasurementTask.run(self, sources, exposure, **kwds)
170 
171  if not sources:
172  return
173 
174  for source in sources:
175  self.dipoleFitter.measure(source, exposure, posExp, negExp)
176 
177 
178 class DipoleModel(object):
179  """!Lightweight class containing methods for generating a dipole model for fitting
180  to sources in diffims, used by DipoleFitAlgorithm.
181 
182  This code is documented in DMTN-007 (http://dmtn-007.lsst.io).
183  """
184 
185  def __init__(self):
186  import lsstDebug
187  self.debug = lsstDebug.Info(__name__).debug
188  self.log = Log.getLogger(__name__)
189 
190  def makeBackgroundModel(self, in_x, pars=None):
191  """!Generate gradient model (2-d array) with up to 2nd-order polynomial
192 
193  @param in_x (2, w, h)-dimensional `numpy.array`, containing the
194  input x,y meshgrid providing the coordinates upon which to
195  compute the gradient. This will typically be generated via
196  `_generateXYGrid()`. `w` and `h` correspond to the width and
197  height of the desired grid.
198  @param pars Up to 6 floats for up
199  to 6 2nd-order 2-d polynomial gradient parameters, in the
200  following order: (intercept, x, y, xy, x**2, y**2). If `pars`
201  is emtpy or `None`, do nothing and return `None` (for speed).
202 
203  @return None, or 2-d numpy.array of width/height matching
204  input bbox, containing computed gradient values.
205  """
206 
207  # Don't fit for other gradient parameters if the intercept is not included.
208  if (pars is None) or (len(pars) <= 0) or (pars[0] is None):
209  return
210 
211  y, x = in_x[0, :], in_x[1, :]
212  gradient = np.full_like(x, pars[0], dtype='float64')
213  if len(pars) > 1 and pars[1] is not None:
214  gradient += pars[1] * x
215  if len(pars) > 2 and pars[2] is not None:
216  gradient += pars[2] * y
217  if len(pars) > 3 and pars[3] is not None:
218  gradient += pars[3] * (x * y)
219  if len(pars) > 4 and pars[4] is not None:
220  gradient += pars[4] * (x * x)
221  if len(pars) > 5 and pars[5] is not None:
222  gradient += pars[5] * (y * y)
223 
224  return gradient
225 
226  def _generateXYGrid(self, bbox):
227  """!Generate a meshgrid covering the x,y coordinates bounded by bbox
228 
229  @param bbox input BoundingBox defining the coordinate limits
230  @return in_x (2, w, h)-dimensional numpy array containing the grid indexing over x- and
231  y- coordinates
232 
233  @see makeBackgroundModel
234  """
235 
236  x, y = np.mgrid[bbox.getBeginY():bbox.getEndY(), bbox.getBeginX():bbox.getEndX()]
237  in_x = np.array([y, x]).astype(np.float64)
238  in_x[0, :] -= np.mean(in_x[0, :])
239  in_x[1, :] -= np.mean(in_x[1, :])
240  return in_x
241 
242  def _getHeavyFootprintSubimage(self, fp, badfill=np.nan, grow=0):
243  """!Extract the image from a `HeavyFootprint` as an `afwImage.ImageF`.
244 
245  @param fp HeavyFootprint to use to generate the subimage
246  @param badfill Value to fill in pixels in extracted image that are outside the footprint
247  @param grow Optionally grow the footprint by this amount before extraction
248 
249  @return an `afwImage.ImageF` containing the subimage
250  """
251  bbox = fp.getBBox()
252  if grow > 0:
253  bbox.grow(grow)
254 
255  subim2 = afwImage.ImageF(bbox, badfill)
256  fp.getSpans().unflatten(subim2.getArray(), fp.getImageArray(), bbox.getCorners()[0])
257  return subim2
258 
259  def fitFootprintBackground(self, source, posImage, order=1):
260  """!Fit a linear (polynomial) model of given order (max 2) to the background of a footprint.
261 
262  Only fit the pixels OUTSIDE of the footprint, but within its bounding box.
263 
264  @param source SourceRecord, the footprint of which is to be fit
265  @param posImage The exposure from which to extract the footprint subimage
266  @param order Polynomial order of background gradient to fit.
267 
268  @return pars `tuple` of length (1 if order==0; 3 if order==1; 6 if order == 2),
269  containing the resulting fit parameters
270 
271  @todo look into whether to use afwMath background methods -- see
272  http://lsst-web.ncsa.illinois.edu/doxygen/x_masterDoxyDoc/_background_example.html
273  """
274 
275  fp = source.getFootprint()
276  bbox = fp.getBBox()
277  bbox.grow(3)
278  posImg = afwImage.ImageF(posImage.getMaskedImage().getImage(), bbox, afwImage.PARENT)
279 
280  # This code constructs the footprint image so that we can identify the pixels that are
281  # outside the footprint (but within the bounding box). These are the pixels used for
282  # fitting the background.
283  posHfp = afwDet.HeavyFootprintF(fp, posImage.getMaskedImage())
284  posFpImg = self._getHeavyFootprintSubimage(posHfp, grow=3)
285 
286  isBg = np.isnan(posFpImg.getArray()).ravel()
287 
288  data = posImg.getArray().ravel()
289  data = data[isBg]
290  B = data
291 
292  x, y = np.mgrid[bbox.getBeginY():bbox.getEndY(), bbox.getBeginX():bbox.getEndX()]
293  x = x.astype(np.float64).ravel()
294  x -= np.mean(x)
295  x = x[isBg]
296  y = y.astype(np.float64).ravel()
297  y -= np.mean(y)
298  y = y[isBg]
299  b = np.ones_like(x, dtype=np.float64)
300 
301  M = np.vstack([b]).T # order = 0
302  if order == 1:
303  M = np.vstack([b, x, y]).T
304  elif order == 2:
305  M = np.vstack([b, x, y, x**2., y**2., x*y]).T
306 
307  pars = np.linalg.lstsq(M, B)[0]
308  return pars
309 
310  def makeStarModel(self, bbox, psf, xcen, ycen, flux):
311  """!Generate model (2-d Image) of a 'star' (single PSF) centered at given coordinates
312 
313  @param bbox Bounding box marking pixel coordinates for generated model
314  @param psf Psf model used to generate the 'star'
315  @param xcen Desired x-centroid of the 'star'
316  @param ycen Desired y-centroid of the 'star'
317  @param flux Desired flux of the 'star'
318 
319  @return 2-d stellar `afwImage.Image` of width/height matching input `bbox`,
320  containing PSF with given centroid and flux
321  """
322 
323  # Generate the psf image, normalize to flux
324  psf_img = psf.computeImage(afwGeom.Point2D(xcen, ycen)).convertF()
325  psf_img_sum = np.nansum(psf_img.getArray())
326  psf_img *= (flux/psf_img_sum)
327 
328  # Clip the PSF image bounding box to fall within the footprint bounding box
329  psf_box = psf_img.getBBox()
330  psf_box.clip(bbox)
331  psf_img = afwImage.ImageF(psf_img, psf_box, afwImage.PARENT)
332 
333  # Then actually crop the psf image.
334  # Usually not necessary, but if the dipole is near the edge of the image...
335  # Would be nice if we could compare original pos_box with clipped pos_box and
336  # see if it actually was clipped.
337  p_Im = afwImage.ImageF(bbox)
338  tmpSubim = afwImage.ImageF(p_Im, psf_box, afwImage.PARENT)
339  tmpSubim += psf_img
340 
341  return p_Im
342 
343  def makeModel(self, x, flux, xcenPos, ycenPos, xcenNeg, ycenNeg, fluxNeg=None,
344  b=None, x1=None, y1=None, xy=None, x2=None, y2=None,
345  bNeg=None, x1Neg=None, y1Neg=None, xyNeg=None, x2Neg=None, y2Neg=None,
346  **kwargs):
347  """!Generate dipole model with given parameters.
348 
349  This is the function whose sum-of-squared difference from data
350  is minimized by `lmfit`.
351 
352  @param x Input independent variable. Used here as the grid on
353  which to compute the background gradient model.
354  @param flux Desired flux of the positive lobe of the dipole
355  @param xcenPos Desired x-centroid of the positive lobe of the dipole
356  @param ycenPos Desired y-centroid of the positive lobe of the dipole
357  @param xcenNeg Desired x-centroid of the negative lobe of the dipole
358  @param ycenNeg Desired y-centroid of the negative lobe of the dipole
359  @param fluxNeg Desired flux of the negative lobe of the dipole, set to 'flux' if None
360  @param b, x1, y1, xy, x2, y2 Gradient parameters for positive lobe.
361  @param bNeg, x1Neg, y1Neg, xyNeg, x2Neg, y2Neg Gradient parameters for negative lobe.
362  They are set to the corresponding positive values if None.
363 
364  @param **kwargs Keyword arguments passed through `lmfit` and
365  used by this function. These must include:
366  - `psf` Psf model used to generate the 'star'
367  - `rel_weight` Used to signify least-squares weighting of posImage/negImage
368  relative to diffim. If `rel_weight == 0` then posImage/negImage are ignored.
369  - `bbox` Bounding box containing region to be modelled
370 
371  @see `makeBackgroundModel` for further parameter descriptions.
372 
373  @return `numpy.array` of width/height matching input bbox,
374  containing dipole model with given centroids and flux(es). If
375  `rel_weight` = 0, this is a 2-d array with dimensions matching
376  those of bbox; otherwise a stack of three such arrays,
377  representing the dipole (diffim), positive and negative images
378  respectively.
379  """
380 
381  psf = kwargs.get('psf')
382  rel_weight = kwargs.get('rel_weight') # if > 0, we're including pre-sub. images
383  fp = kwargs.get('footprint')
384  bbox = fp.getBBox()
385 
386  if fluxNeg is None:
387  fluxNeg = flux
388 
389  if self.debug:
390  self.log.debug('%.2f %.2f %.2f %.2f %.2f %.2f',
391  flux, fluxNeg, xcenPos, ycenPos, xcenNeg, ycenNeg)
392  if x1 is not None:
393  self.log.debug(' %.2f %.2f %.2f', b, x1, y1)
394  if xy is not None:
395  self.log.debug(' %.2f %.2f %.2f', xy, x2, y2)
396 
397  posIm = self.makeStarModel(bbox, psf, xcenPos, ycenPos, flux)
398  negIm = self.makeStarModel(bbox, psf, xcenNeg, ycenNeg, fluxNeg)
399 
400  in_x = x
401  if in_x is None: # use the footprint to generate the input grid
402  y, x = np.mgrid[bbox.getBeginY():bbox.getEndY(), bbox.getBeginX():bbox.getEndX()]
403  in_x = np.array([x, y]) * 1.
404  in_x[0, :] -= in_x[0, :].mean() # center it!
405  in_x[1, :] -= in_x[1, :].mean()
406 
407  if b is not None:
408  gradient = self.makeBackgroundModel(in_x, (b, x1, y1, xy, x2, y2))
409 
410  # If bNeg is None, then don't fit the negative background separately
411  if bNeg is not None:
412  gradientNeg = self.makeBackgroundModel(in_x, (bNeg, x1Neg, y1Neg, xyNeg, x2Neg, y2Neg))
413  else:
414  gradientNeg = gradient
415 
416  posIm.getArray()[:, :] += gradient
417  negIm.getArray()[:, :] += gradientNeg
418 
419  # Generate the diffIm model
420  diffIm = afwImage.ImageF(bbox)
421  diffIm += posIm
422  diffIm -= negIm
423 
424  zout = diffIm.getArray()
425  if rel_weight > 0.:
426  zout = np.append([zout], [posIm.getArray(), negIm.getArray()], axis=0)
427 
428  return zout
429 
430 
431 class DipoleFitAlgorithm(object):
432  """!Lightweight class containing methods for fitting a dipole model in
433  a diffim, used by DipoleFitPlugin.
434 
435  This code is documented in DMTN-007 (http://dmtn-007.lsst.io).
436 
437  Below is a (somewhat incomplete) list of improvements
438  that would be worth investigating, given the time:
439 
440  @todo 1. evaluate necessity for separate parameters for pos- and neg- images
441  @todo 2. only fit background OUTSIDE footprint (DONE) and dipole params INSIDE footprint (NOT DONE)?
442  @todo 3. correct normalization of least-squares weights based on variance planes
443  @todo 4. account for PSFs that vary across the exposures (should be happening by default?)
444  @todo 5. correctly account for NA/masks (i.e., ignore!)
445  @todo 6. better exception handling in the plugin
446  @todo 7. better classification of dipoles (e.g. by comparing chi2 fit vs. monopole?)
447  @todo 8. (DONE) Initial fast estimate of background gradient(s) params -- perhaps using numpy.lstsq
448  @todo 9. (NOT NEEDED - see (2)) Initial fast test whether a background gradient needs to be fit
449  @todo 10. (DONE) better initial estimate for flux when there's a strong gradient
450  @todo 11. (DONE) requires a new package `lmfit` -- investiate others? (astropy/scipy/iminuit?)
451  """
452 
453  # This is just a private version number to sync with the ipython notebooks that I have been
454  # using for algorithm development.
455  _private_version_ = '0.0.5'
456 
457  def __init__(self, diffim, posImage=None, negImage=None):
458  """!Algorithm to run dipole measurement on a diaSource
459 
460  @param diffim Exposure on which the diaSources were detected
461  @param posImage "Positive" exposure from which the template was subtracted
462  @param negImage "Negative" exposure which was subtracted from the posImage
463  """
464 
465  self.diffim = diffim
466  self.posImage = posImage
467  self.negImage = negImage
468  self.psfSigma = None
469  if diffim is not None:
470  self.psfSigma = diffim.getPsf().computeShape().getDeterminantRadius()
471 
472  self.log = Log.getLogger(__name__)
473 
474  import lsstDebug
475  self.debug = lsstDebug.Info(__name__).debug
476 
477  def fitDipoleImpl(self, source, tol=1e-7, rel_weight=0.5,
478  fitBackground=1, bgGradientOrder=1, maxSepInSigma=5.,
479  separateNegParams=True, verbose=False):
480  """!Fit a dipole model to an input difference image.
481 
482  Actually, fits the subimage bounded by the input source's
483  footprint) and optionally constrain the fit using the
484  pre-subtraction images posImage and negImage.
485 
486  @return `lmfit.MinimizerResult` object containing the fit
487  parameters and other information.
488 
489  @see `fitDipole()`
490  """
491 
492  # Only import lmfit if someone wants to use the new DipoleFitAlgorithm.
493  import lmfit
494 
495  fp = source.getFootprint()
496  bbox = fp.getBBox()
497  subim = afwImage.MaskedImageF(self.diffim.getMaskedImage(), bbox=bbox, origin=afwImage.PARENT)
498 
499  z = diArr = subim.getArrays()[0]
500  weights = 1. / subim.getArrays()[2] # get the weights (=1/variance)
501 
502  if rel_weight > 0. and ((self.posImage is not None) or (self.negImage is not None)):
503  if self.negImage is not None:
504  negSubim = afwImage.MaskedImageF(self.negImage.getMaskedImage(), bbox, origin=afwImage.PARENT)
505  if self.posImage is not None:
506  posSubim = afwImage.MaskedImageF(self.posImage.getMaskedImage(), bbox, origin=afwImage.PARENT)
507  if self.posImage is None: # no science image provided; generate it from diffim + negImage
508  posSubim = subim.clone()
509  posSubim += negSubim
510  if self.negImage is None: # no template provided; generate it from the posImage - diffim
511  negSubim = posSubim.clone()
512  negSubim -= subim
513 
514  z = np.append([z], [posSubim.getArrays()[0],
515  negSubim.getArrays()[0]], axis=0)
516  # Weight the pos/neg images by rel_weight relative to the diffim
517  weights = np.append([weights], [1. / posSubim.getArrays()[2] * rel_weight,
518  1. / negSubim.getArrays()[2] * rel_weight], axis=0)
519  else:
520  rel_weight = 0. # a short-cut for "don't include the pre-subtraction data"
521 
522  # It seems that `lmfit` requires a static functor as its optimized method, which eliminates
523  # the ability to pass a bound method or other class method. Here we write a wrapper which
524  # makes this possible.
525  def dipoleModelFunctor(x, flux, xcenPos, ycenPos, xcenNeg, ycenNeg, fluxNeg=None,
526  b=None, x1=None, y1=None, xy=None, x2=None, y2=None,
527  bNeg=None, x1Neg=None, y1Neg=None, xyNeg=None, x2Neg=None, y2Neg=None,
528  **kwargs):
529  """!Generate dipole model with given parameters.
530 
531  It simply defers to `modelObj.makeModel()`, where `modelObj` comes
532  out of `kwargs['modelObj']`.
533 
534  @see DipoleModel.makeModel
535  """
536  modelObj = kwargs.pop('modelObj')
537  return modelObj.makeModel(x, flux, xcenPos, ycenPos, xcenNeg, ycenNeg, fluxNeg=fluxNeg,
538  b=b, x1=x1, y1=y1, xy=xy, x2=x2, y2=y2,
539  bNeg=bNeg, x1Neg=x1Neg, y1Neg=y1Neg, xyNeg=xyNeg,
540  x2Neg=x2Neg, y2Neg=y2Neg, **kwargs)
541 
542  dipoleModel = DipoleModel()
543 
544  modelFunctor = dipoleModelFunctor # dipoleModel.makeModel does not work for now.
545  # Create the lmfit model (lmfit uses scipy 'leastsq' option by default - Levenberg-Marquardt)
546  # Note we can also tell it to drop missing values from the data.
547  gmod = lmfit.Model(modelFunctor, verbose=verbose, missing='drop')
548  # independent_vars=independent_vars) #, param_names=param_names)
549 
550  # Add the constraints for centroids, fluxes.
551  # starting constraint - near centroid of footprint
552  fpCentroid = np.array([fp.getCentroid().getX(), fp.getCentroid().getY()])
553  cenNeg = cenPos = fpCentroid
554 
555  pks = fp.getPeaks()
556 
557  if len(pks) >= 1:
558  cenPos = pks[0].getF() # if individual (merged) peaks were detected, use those
559  if len(pks) >= 2: # peaks are already sorted by centroid flux so take the most negative one
560  cenNeg = pks[-1].getF()
561 
562  # For close/faint dipoles the starting locs (min/max) might be way off, let's help them a bit.
563  # First assume dipole is not separated by more than 5*psfSigma.
564  maxSep = self.psfSigma * maxSepInSigma
565 
566  # As an initial guess -- assume the dipole is close to the center of the footprint.
567  if np.sum(np.sqrt((np.array(cenPos) - fpCentroid)**2.)) > maxSep:
568  cenPos = fpCentroid
569  if np.sum(np.sqrt((np.array(cenNeg) - fpCentroid)**2.)) > maxSep:
570  cenPos = fpCentroid
571 
572  # parameter hints/constraints: https://lmfit.github.io/lmfit-py/model.html#model-param-hints-section
573  # might make sense to not use bounds -- see http://lmfit.github.io/lmfit-py/bounds.html
574  # also see this discussion -- https://github.com/scipy/scipy/issues/3129
575  gmod.set_param_hint('xcenPos', value=cenPos[0],
576  min=cenPos[0]-maxSep, max=cenPos[0]+maxSep)
577  gmod.set_param_hint('ycenPos', value=cenPos[1],
578  min=cenPos[1]-maxSep, max=cenPos[1]+maxSep)
579  gmod.set_param_hint('xcenNeg', value=cenNeg[0],
580  min=cenNeg[0]-maxSep, max=cenNeg[0]+maxSep)
581  gmod.set_param_hint('ycenNeg', value=cenNeg[1],
582  min=cenNeg[1]-maxSep, max=cenNeg[1]+maxSep)
583 
584  # Use the (flux under the dipole)*5 for an estimate.
585  # Lots of testing showed that having startingFlux be too high was better than too low.
586  startingFlux = np.nansum(np.abs(diArr) - np.nanmedian(np.abs(diArr))) * 5.
587  posFlux = negFlux = startingFlux
588 
589  # TBD: set max. flux limit?
590  gmod.set_param_hint('flux', value=posFlux, min=0.1)
591 
592  if separateNegParams:
593  # TBD: set max negative lobe flux limit?
594  gmod.set_param_hint('fluxNeg', value=np.abs(negFlux), min=0.1)
595 
596  # Fixed parameters (don't fit for them if there are no pre-sub images or no gradient fit requested):
597  # Right now (fitBackground == 1), we fit a linear model to the background and then subtract
598  # it from the data and then don't fit the background again (this is faster).
599  # A slower alternative (fitBackground == 2) is to use the estimated background parameters as
600  # starting points in the integrated model fit. That is currently not performed by default,
601  # but might be desirable in some cases.
602  bgParsPos = bgParsNeg = (0., 0., 0.)
603  if ((rel_weight > 0.) and (fitBackground != 0) and (bgGradientOrder >= 0)):
604  pbg = 0.
605  bgFitImage = self.posImage if self.posImage is not None else self.negImage
606  # Fit the gradient to the background (linear model)
607  bgParsPos = bgParsNeg = dipoleModel.fitFootprintBackground(source, bgFitImage,
608  order=bgGradientOrder)
609 
610  # Generate the gradient and subtract it from the pre-subtraction image data
611  if fitBackground == 1:
612  in_x = dipoleModel._generateXYGrid(bbox)
613  pbg = dipoleModel.makeBackgroundModel(in_x, tuple(bgParsPos))
614  z[1, :] -= pbg
615  z[1, :] -= np.nanmedian(z[1, :])
616  posFlux = np.nansum(z[1, :])
617  gmod.set_param_hint('flux', value=posFlux*1.5, min=0.1)
618 
619  if separateNegParams and self.negImage is not None:
620  bgParsNeg = dipoleModel.fitFootprintBackground(source, self.negImage,
621  order=bgGradientOrder)
622  pbg = dipoleModel.makeBackgroundModel(in_x, tuple(bgParsNeg))
623  z[2, :] -= pbg
624  z[2, :] -= np.nanmedian(z[2, :])
625  if separateNegParams:
626  negFlux = np.nansum(z[2, :])
627  gmod.set_param_hint('fluxNeg', value=negFlux*1.5, min=0.1)
628 
629  # Do not subtract the background from the images but include the background parameters in the fit
630  if fitBackground == 2:
631  if bgGradientOrder >= 0:
632  gmod.set_param_hint('b', value=bgParsPos[0])
633  if separateNegParams:
634  gmod.set_param_hint('bNeg', value=bgParsNeg[0])
635  if bgGradientOrder >= 1:
636  gmod.set_param_hint('x1', value=bgParsPos[1])
637  gmod.set_param_hint('y1', value=bgParsPos[2])
638  if separateNegParams:
639  gmod.set_param_hint('x1Neg', value=bgParsNeg[1])
640  gmod.set_param_hint('y1Neg', value=bgParsNeg[2])
641  if bgGradientOrder >= 2:
642  gmod.set_param_hint('xy', value=bgParsPos[3])
643  gmod.set_param_hint('x2', value=bgParsPos[4])
644  gmod.set_param_hint('y2', value=bgParsPos[5])
645  if separateNegParams:
646  gmod.set_param_hint('xyNeg', value=bgParsNeg[3])
647  gmod.set_param_hint('x2Neg', value=bgParsNeg[4])
648  gmod.set_param_hint('y2Neg', value=bgParsNeg[5])
649 
650  y, x = np.mgrid[bbox.getBeginY():bbox.getEndY(), bbox.getBeginX():bbox.getEndX()]
651  in_x = np.array([x, y]).astype(np.float)
652  in_x[0, :] -= in_x[0, :].mean() # center it!
653  in_x[1, :] -= in_x[1, :].mean()
654 
655  # Instead of explicitly using a mask to ignore flagged pixels, just set the ignored pixels'
656  # weights to 0 in the fit. TBD: need to inspect mask planes to set this mask.
657  mask = np.ones_like(z, dtype=bool) # TBD: set mask values to False if the pixels are to be ignored
658 
659  # I'm not sure about the variance planes in the diffim (or convolved pre-sub. images
660  # for that matter) so for now, let's just do an un-weighted least-squares fit
661  # (override weights computed above).
662  weights = mask.astype(np.float64)
663  if self.posImage is not None and rel_weight > 0.:
664  weights = np.array([np.ones_like(diArr), np.ones_like(diArr)*rel_weight,
665  np.ones_like(diArr)*rel_weight])
666 
667  # Set the weights to zero if mask is False
668  if np.any(~mask):
669  weights[~mask] = 0.
670 
671  # Note that although we can, we're not required to set initial values for params here,
672  # since we set their param_hint's above.
673  # Can add "method" param to not use 'leastsq' (==levenberg-marquardt), e.g. "method='nelder'"
674  with warnings.catch_warnings():
675  warnings.simplefilter("ignore") # temporarily turn off silly lmfit warnings
676  result = gmod.fit(z, weights=weights, x=in_x,
677  verbose=verbose,
678  fit_kws={'ftol': tol, 'xtol': tol, 'gtol': tol,
679  'maxfev': 250}, # see scipy docs
680  psf=self.diffim.getPsf(), # hereon: kwargs that get passed to genDipoleModel()
681  rel_weight=rel_weight,
682  footprint=fp,
683  modelObj=dipoleModel)
684 
685  if verbose: # the ci_report() seems to fail if neg params are constrained -- TBD why.
686  # Never wanted in production - this takes a long time (longer than the fit!)
687  # This is how to get confidence intervals out:
688  # https://lmfit.github.io/lmfit-py/confidence.html and
689  # http://cars9.uchicago.edu/software/python/lmfit/model.html
690  print(result.fit_report(show_correl=False))
691  if separateNegParams:
692  print(result.ci_report())
693 
694  return result
695 
696  def fitDipole(self, source, tol=1e-7, rel_weight=0.1,
697  fitBackground=1, maxSepInSigma=5., separateNegParams=True,
698  bgGradientOrder=1, verbose=False, display=False):
699  """!Wrapper around `fitDipoleImpl()` which performs the fit of a dipole
700  model to an input diaSource.
701 
702  Actually, fits the subimage bounded by the input source's
703  footprint) and optionally constrain the fit using the
704  pre-subtraction images self.posImage (science) and
705  self.negImage (template). Wraps the output into a
706  `pipeBase.Struct` named tuple after computing additional
707  statistics such as orientation and SNR.
708 
709  @param source Record containing the (merged) dipole source footprint detected on the diffim
710  @param tol Tolerance parameter for scipy.leastsq() optimization
711  @param rel_weight Weighting of posImage/negImage relative to the diffim in the fit
712  @param fitBackground How to fit linear background gradient in posImage/negImage (see notes)
713  @param bgGradientOrder Desired polynomial order of background gradient (allowed are [0,1,2])
714  @param maxSepInSigma Allowed window of centroid parameters relative to peak in input source footprint
715  @param separateNegParams Fit separate parameters to the flux and background gradient in
716  the negative images? If true, this adds a separate parameter for the negative flux, and [1, 3, or 6]
717  additional parameters to fit for the background gradient in the negImage. Otherwise, the flux and
718  gradient parameters are constrained to be exactly equal in the fit.
719  @param verbose Be verbose
720  @param display Display input data, best fit model(s) and residuals in a matplotlib window.
721 
722  @return `pipeBase.Struct` object containing the fit parameters and other information.
723  @return `lmfit.MinimizerResult` object for debugging and error estimation, etc.
724 
725  @note Parameter `fitBackground` has three options, thus it is an integer:
726  - 0: do not fit background at all
727  - 1 (default): pre-fit the background using linear least squares and then do not fit it as part
728  of the dipole fitting optimization
729  - 2: pre-fit the background using linear least squares (as in 1), and use the parameter
730  estimates from that fit as starting parameters for an integrated "re-fit" of the background
731  as part of the overall dipole fitting optimization.
732  """
733 
734  fitResult = self.fitDipoleImpl(
735  source, tol=tol, rel_weight=rel_weight, fitBackground=fitBackground,
736  maxSepInSigma=maxSepInSigma, separateNegParams=separateNegParams,
737  bgGradientOrder=bgGradientOrder, verbose=verbose)
738 
739  # Display images, model fits and residuals (currently uses matplotlib display functions)
740  if display:
741  fp = source.getFootprint()
742  self.displayFitResults(fp, fitResult)
743 
744  fitParams = fitResult.best_values
745  if fitParams['flux'] <= 1.: # usually around 0.1 -- the minimum flux allowed -- i.e. bad fit.
746  out = Struct(posCentroidX=np.nan, posCentroidY=np.nan,
747  negCentroidX=np.nan, negCentroidY=np.nan,
748  posFlux=np.nan, negFlux=np.nan, posFluxSigma=np.nan, negFluxSigma=np.nan,
749  centroidX=np.nan, centroidY=np.nan, orientation=np.nan,
750  signalToNoise=np.nan, chi2=np.nan, redChi2=np.nan)
751  return out, fitResult
752 
753  centroid = ((fitParams['xcenPos'] + fitParams['xcenNeg']) / 2.,
754  (fitParams['ycenPos'] + fitParams['ycenNeg']) / 2.)
755  dx, dy = fitParams['xcenPos'] - fitParams['xcenNeg'], fitParams['ycenPos'] - fitParams['ycenNeg']
756  angle = np.arctan2(dy, dx) / np.pi * 180. # convert to degrees (should keep as rad?)
757 
758  # Exctract flux value, compute signalToNoise from flux/variance_within_footprint
759  # Also extract the stderr of flux estimate.
760  def computeSumVariance(exposure, footprint):
761  box = footprint.getBBox()
762  subim = afwImage.MaskedImageF(exposure.getMaskedImage(), box, origin=afwImage.PARENT)
763  return np.sqrt(np.nansum(subim.getArrays()[1][:, :]))
764 
765  fluxVal = fluxVar = fitParams['flux']
766  fluxErr = fluxErrNeg = fitResult.params['flux'].stderr
767  if self.posImage is not None:
768  fluxVar = computeSumVariance(self.posImage, source.getFootprint())
769  else:
770  fluxVar = computeSumVariance(self.diffim, source.getFootprint())
771 
772  fluxValNeg, fluxVarNeg = fluxVal, fluxVar
773  if separateNegParams:
774  fluxValNeg = fitParams['fluxNeg']
775  fluxErrNeg = fitResult.params['fluxNeg'].stderr
776  if self.negImage is not None:
777  fluxVarNeg = computeSumVariance(self.negImage, source.getFootprint())
778 
779  try:
780  signalToNoise = np.sqrt((fluxVal/fluxVar)**2 + (fluxValNeg/fluxVarNeg)**2)
781  except: # catch divide by zero - should never happen.
782  signalToNoise = np.nan
783 
784  out = Struct(posCentroidX=fitParams['xcenPos'], posCentroidY=fitParams['ycenPos'],
785  negCentroidX=fitParams['xcenNeg'], negCentroidY=fitParams['ycenNeg'],
786  posFlux=fluxVal, negFlux=-fluxValNeg, posFluxSigma=fluxErr, negFluxSigma=fluxErrNeg,
787  centroidX=centroid[0], centroidY=centroid[1], orientation=angle,
788  signalToNoise=signalToNoise, chi2=fitResult.chisqr, redChi2=fitResult.redchi)
789 
790  # fitResult may be returned for debugging
791  return out, fitResult
792 
793  def displayFitResults(self, footprint, result):
794  """!Display data, model fits and residuals (currently uses matplotlib display functions).
795 
796  @param footprint Footprint containing the dipole that was fit
797  @param result `lmfit.MinimizerResult` object returned by `lmfit` optimizer
798  """
799  try:
800  import matplotlib.pyplot as plt
801  except ImportError as err:
802  self.log.warn('Unable to import matplotlib: %s', err)
803  raise err
804 
805  def display2dArray(arr, title='Data', extent=None):
806  """!Use `matplotlib.pyplot.imshow` to display a 2-D array with a given coordinate range.
807  """
808  fig = plt.imshow(arr, origin='lower', interpolation='none', cmap='gray', extent=extent)
809  plt.title(title)
810  plt.colorbar(fig, cmap='gray')
811  return fig
812 
813  z = result.data
814  fit = result.best_fit
815  bbox = footprint.getBBox()
816  extent = (bbox.getBeginX(), bbox.getEndX(), bbox.getBeginY(), bbox.getEndY())
817  if z.shape[0] == 3:
818  fig = plt.figure(figsize=(8, 8))
819  for i in range(3):
820  plt.subplot(3, 3, i*3+1)
821  display2dArray(z[i, :], 'Data', extent=extent)
822  plt.subplot(3, 3, i*3+2)
823  display2dArray(fit[i, :], 'Model', extent=extent)
824  plt.subplot(3, 3, i*3+3)
825  display2dArray(z[i, :] - fit[i, :], 'Residual', extent=extent)
826  return fig
827  else:
828  fig = plt.figure(figsize=(8, 2.5))
829  plt.subplot(1, 3, 1)
830  display2dArray(z, 'Data', extent=extent)
831  plt.subplot(1, 3, 2)
832  display2dArray(fit, 'Model', extent=extent)
833  plt.subplot(1, 3, 3)
834  display2dArray(z - fit, 'Residual', extent=extent)
835  return fig
836 
837  plt.show()
838 
839 
840 @measBase.register("ip_diffim_DipoleFit")
841 class DipoleFitPlugin(measBase.SingleFramePlugin):
842  """!Subclass of SingleFramePlugin which fits dipoles to all merged (two-peak) diaSources
843 
844  Accepts up to three input images in its `measure` method. If these are
845  provided, it includes data from the pre-subtraction posImage
846  (science image) and optionally negImage (template image) to
847  constrain the fit. The meat of the fitting routines are in the
848  class DipoleFitAlgorithm.
849 
850  The motivation behind this plugin and the necessity for including more than
851  one exposure are documented in DMTN-007 (http://dmtn-007.lsst.io).
852 
853  This class is named `ip_diffim_DipoleFit` so that it may be used alongside
854  the existing `ip_diffim_DipoleMeasurement` classes until such a time as those
855  are deemed to be replaceable by this.
856  """
857 
858  ConfigClass = DipoleFitPluginConfig
859  DipoleFitAlgorithmClass = DipoleFitAlgorithm # Pointer to the class that performs the fit
860 
861  FAILURE_EDGE = 1 # too close to the edge
862  FAILURE_FIT = 2 # failure in the fitting
863  FAILURE_NOT_DIPOLE = 4 # input source is not a putative dipole to begin with
864 
865  @classmethod
867  """!Set execution order to `FLUX_ORDER`.
868 
869  This includes algorithms that require both `getShape()` and `getCentroid()`,
870  in addition to a Footprint and its Peaks.
871  """
872  return cls.FLUX_ORDER
873 
874  def __init__(self, config, name, schema, metadata):
875  measBase.SingleFramePlugin.__init__(self, config, name, schema, metadata)
876 
877  self.log = Log.getLogger(name)
878 
879  self._setupSchema(config, name, schema, metadata)
880 
881  def _setupSchema(self, config, name, schema, metadata):
882  # Get a FunctorKey that can quickly look up the "blessed" centroid value.
883  self.centroidKey = afwTable.Point2DKey(schema["slot_Centroid"])
884 
885  # Add some fields for our outputs, and save their Keys.
886  # Use setattr() to programmatically set the pos/neg named attributes to values, e.g.
887  # self.posCentroidKeyX = 'ip_diffim_DipoleFit_pos_centroid_x'
888 
889  for pos_neg in ['pos', 'neg']:
890 
891  key = schema.addField(
892  schema.join(name, pos_neg, "flux"), type=float, units="count",
893  doc="Dipole {0} lobe flux".format(pos_neg))
894  setattr(self, ''.join((pos_neg, 'FluxKey')), key)
895 
896  key = schema.addField(
897  schema.join(name, pos_neg, "fluxSigma"), type=float, units="pixel",
898  doc="1-sigma uncertainty for {0} dipole flux".format(pos_neg))
899  setattr(self, ''.join((pos_neg, 'FluxSigmaKey')), key)
900 
901  for x_y in ['x', 'y']:
902  key = schema.addField(
903  schema.join(name, pos_neg, "centroid", x_y), type=float, units="pixel",
904  doc="Dipole {0} lobe centroid".format(pos_neg))
905  setattr(self, ''.join((pos_neg, 'CentroidKey', x_y.upper())), key)
906 
907  for x_y in ['x', 'y']:
908  key = schema.addField(
909  schema.join(name, "centroid", x_y), type=float, units="pixel",
910  doc="Dipole centroid")
911  setattr(self, ''.join(('centroidKey', x_y.upper())), key)
912 
913  self.fluxKey = schema.addField(
914  schema.join(name, "flux"), type=float, units="count",
915  doc="Dipole overall flux")
916 
917  self.orientationKey = schema.addField(
918  schema.join(name, "orientation"), type=float, units="deg",
919  doc="Dipole orientation")
920 
921  self.separationKey = schema.addField(
922  schema.join(name, "separation"), type=float, units="pixel",
923  doc="Pixel separation between positive and negative lobes of dipole")
924 
925  self.chi2dofKey = schema.addField(
926  schema.join(name, "chi2dof"), type=float,
927  doc="Chi2 per degree of freedom of dipole fit")
928 
929  self.signalToNoiseKey = schema.addField(
930  schema.join(name, "signalToNoise"), type=float,
931  doc="Estimated signal-to-noise of dipole fit")
932 
933  self.classificationFlagKey = schema.addField(
934  schema.join(name, "flag", "classification"), type="Flag",
935  doc="Flag indicating diaSource is classified as a dipole")
936 
937  self.classificationAttemptedFlagKey = schema.addField(
938  schema.join(name, "flag", "classificationAttempted"), type="Flag",
939  doc="Flag indicating diaSource was attempted to be classified as a dipole")
940 
941  self.flagKey = schema.addField(
942  schema.join(name, "flag"), type="Flag",
943  doc="General failure flag for dipole fit")
944 
945  self.edgeFlagKey = schema.addField(
946  schema.join(name, "flag", "edge"), type="Flag",
947  doc="Flag set when dipole is too close to edge of image")
948 
949  def measure(self, measRecord, exposure, posExp=None, negExp=None):
950  """!Perform the non-linear least squares minimization on the putative dipole source.
951 
952  @param measRecord diaSources that will be measured using dipole measurement
953  @param exposure Difference exposure on which the diaSources were detected; `exposure = posExp-negExp`
954  @param posExp "Positive" exposure, typically a science exposure, or None if unavailable
955  @param negExp "Negative" exposure, typically a template exposure, or None if unavailable
956 
957  @note When `posExp` is `None`, will compute `posImage = exposure + negExp`.
958  Likewise, when `negExp` is `None`, will compute `negImage = posExp - exposure`.
959  If both `posExp` and `negExp` are `None`, will attempt to fit the dipole to just the `exposure`
960  with no constraint.
961 
962  The main functionality of this routine was placed outside of
963  this plugin (into `DipoleFitAlgorithm.fitDipole()`) so that
964  `DipoleFitAlgorithm.fitDipole()` can be called separately for
965  testing (@see `tests/testDipoleFitter.py`)
966  """
967 
968  result = None
969  pks = measRecord.getFootprint().getPeaks()
970 
971  # Check if the footprint consists of a putative dipole - else don't fit it.
972  if (
973  (len(pks) <= 1) or # one peak in the footprint - not a dipole
974  (len(pks) > 1 and (np.sign(pks[0].getPeakValue()) ==
975  np.sign(pks[-1].getPeakValue()))) # peaks are same sign - not a dipole
976  ):
977  measRecord.set(self.classificationFlagKey, False)
978  measRecord.set(self.classificationAttemptedFlagKey, False)
979  self.fail(measRecord, measBase.MeasurementError('not a dipole', self.FAILURE_NOT_DIPOLE))
980  if not self.config.fitAllDiaSources:
981  return result
982 
983  try:
984  alg = self.DipoleFitAlgorithmClass(exposure, posImage=posExp, negImage=negExp)
985  result, _ = alg.fitDipole(
986  measRecord, rel_weight=self.config.relWeight,
987  tol=self.config.tolerance,
988  maxSepInSigma=self.config.maxSeparation,
989  fitBackground=self.config.fitBackground,
990  separateNegParams=self.config.fitSeparateNegParams,
991  verbose=False, display=False)
992  except pexExcept.LengthError:
993  self.fail(measRecord, measBase.MeasurementError('edge failure', self.FAILURE_EDGE))
994  except Exception:
995  self.fail(measRecord, measBase.MeasurementError('dipole fit failure', self.FAILURE_FIT))
996 
997  if result is None:
998  measRecord.set(self.classificationFlagKey, False)
999  measRecord.set(self.classificationAttemptedFlagKey, False)
1000  return result
1001 
1002  self.log.debug("Dipole fit result: %d %s", measRecord.getId(), str(result))
1003 
1004  if result.posFlux <= 1.: # usually around 0.1 -- the minimum flux allowed -- i.e. bad fit.
1005  self.fail(measRecord, measBase.MeasurementError('dipole fit failure', self.FAILURE_FIT))
1006 
1007  # add chi2, coord/flux uncertainties (TBD), dipole classification
1008  # Add the relevant values to the measRecord
1009  measRecord[self.posFluxKey] = result.posFlux
1010  measRecord[self.posFluxSigmaKey] = result.signalToNoise # to be changed to actual sigma!
1011  measRecord[self.posCentroidKeyX] = result.posCentroidX
1012  measRecord[self.posCentroidKeyY] = result.posCentroidY
1013 
1014  measRecord[self.negFluxKey] = result.negFlux
1015  measRecord[self.negFluxSigmaKey] = result.signalToNoise # to be changed to actual sigma!
1016  measRecord[self.negCentroidKeyX] = result.negCentroidX
1017  measRecord[self.negCentroidKeyY] = result.negCentroidY
1018 
1019  # Dia source flux: average of pos+neg
1020  measRecord[self.fluxKey] = (abs(result.posFlux) + abs(result.negFlux))/2.
1021  measRecord[self.orientationKey] = result.orientation
1022  measRecord[self.separationKey] = np.sqrt((result.posCentroidX - result.negCentroidX)**2. +
1023  (result.posCentroidY - result.negCentroidY)**2.)
1024  measRecord[self.centroidKeyX] = result.centroidX
1025  measRecord[self.centroidKeyY] = result.centroidY
1026 
1027  measRecord[self.signalToNoiseKey] = result.signalToNoise
1028  measRecord[self.chi2dofKey] = result.redChi2
1029 
1030  self.doClassify(measRecord, result.chi2)
1031 
1032  def doClassify(self, measRecord, chi2val):
1033  """!Determine if source is classified as dipole via three criteria:
1034  - does the total signal-to-noise surpass the minSn?
1035  - are the pos/neg fluxes greater than 1.0 and no more than 0.65 (param `maxFluxRatio`)
1036  of the total flux? By default this will never happen since `posFlux == negFlux`.
1037  - is it a good fit (`chi2dof` < 1)? (Currently not used.)
1038  """
1039 
1040  # First, does the total signal-to-noise surpass the minSn?
1041  passesSn = measRecord[self.signalToNoiseKey] > self.config.minSn
1042 
1043  # Second, are the pos/neg fluxes greater than 1.0 and no more than 0.65 (param maxFluxRatio)
1044  # of the total flux? By default this will never happen since posFlux = negFlux.
1045  passesFluxPos = (abs(measRecord[self.posFluxKey]) /
1046  (measRecord[self.fluxKey]*2.)) < self.config.maxFluxRatio
1047  passesFluxPos &= (abs(measRecord[self.posFluxKey]) >= 1.0)
1048  passesFluxNeg = (abs(measRecord[self.negFluxKey]) /
1049  (measRecord[self.fluxKey]*2.)) < self.config.maxFluxRatio
1050  passesFluxNeg &= (abs(measRecord[self.negFluxKey]) >= 1.0)
1051  allPass = (passesSn and passesFluxPos and passesFluxNeg) # and passesChi2)
1052 
1053  # Third, is it a good fit (chi2dof < 1)?
1054  # Use scipy's chi2 cumulative distrib to estimate significance
1055  # This doesn't really work since I don't trust the values in the variance plane (which
1056  # affects the least-sq weights, which affects the resulting chi2).
1057  # But I'm going to keep this here for future use.
1058  if False:
1059  from scipy.stats import chi2
1060  ndof = chi2val / measRecord[self.chi2dofKey]
1061  significance = chi2.cdf(chi2val, ndof)
1062  passesChi2 = significance < self.config.maxChi2DoF
1063  allPass = allPass and passesChi2
1064 
1065  measRecord.set(self.classificationAttemptedFlagKey, True)
1066 
1067  if allPass: # Note cannot pass `allPass` into the `measRecord.set()` call below...?
1068  measRecord.set(self.classificationFlagKey, True)
1069  else:
1070  measRecord.set(self.classificationFlagKey, False)
1071 
1072  def fail(self, measRecord, error=None):
1073  """!Catch failures and set the correct flags.
1074  """
1075 
1076  measRecord.set(self.flagKey, True)
1077  if error is not None:
1078  if error.getFlagBit() == self.FAILURE_EDGE:
1079  self.log.warn('DipoleFitPlugin not run on record %d: %s', measRecord.getId(), str(error))
1080  measRecord.set(self.edgeFlagKey, True)
1081  if error.getFlagBit() == self.FAILURE_FIT:
1082  self.log.warn('DipoleFitPlugin failed on record %d: %s', measRecord.getId(), str(error))
1083  measRecord.set(self.flagKey, True)
1084  if error.getFlagBit() == self.FAILURE_NOT_DIPOLE:
1085  self.log.debug('DipoleFitPlugin not run on record %d: %s',
1086  measRecord.getId(), str(error))
1087  measRecord.set(self.classificationAttemptedFlagKey, False)
1088  measRecord.set(self.flagKey, True)
1089  else:
1090  self.log.warn('DipoleFitPlugin failed on record %d', measRecord.getId())
Subclass of SingleFramePlugin which fits dipoles to all merged (two-peak) diaSources.
Lightweight class containing methods for generating a dipole model for fitting to sources in diffims...
Lightweight class containing methods for fitting a dipole model in a diffim, used by DipoleFitPlugin...
def doClassify
Determine if source is classified as dipole via three criteria:
def displayFitResults
Display data, model fits and residuals (currently uses matplotlib display functions).
def _getHeavyFootprintSubimage
Extract the image from a HeavyFootprint as an afwImage.ImageF.
def fitDipole
Wrapper around fitDipoleImpl() which performs the fit of a dipole model to an input diaSource...
def _generateXYGrid
Generate a meshgrid covering the x,y coordinates bounded by bbox.
def makeBackgroundModel
Generate gradient model (2-d array) with up to 2nd-order polynomial.
def fitFootprintBackground
Fit a linear (polynomial) model of given order (max 2) to the background of a footprint.
def measure
Perform the non-linear least squares minimization on the putative dipole source.
def makeStarModel
Generate model (2-d Image) of a &#39;star&#39; (single PSF) centered at given coordinates.
def getExecutionOrder
Set execution order to FLUX_ORDER.
Measurement of detected diaSources as dipoles.
def makeModel
Generate dipole model with given parameters.
Subclass of SingleFrameMeasurementTask which accepts up to three input images in its run() method...
def fitDipoleImpl
Fit a dipole model to an input difference image.
def run
Run dipole measurement and classification.
def fail
Catch failures and set the correct flags.
def __init__
Algorithm to run dipole measurement on a diaSource.