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