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