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