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