Coverage for python/lsst/ip/diffim/dipoleFitTask.py: 11%

442 statements  

« prev     ^ index     » next       coverage.py v6.5.0, created at 2023-03-09 11:56 +0000

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 

23import logging 

24import numpy as np 

25import warnings 

26 

27import lsst.afw.image as afwImage 

28import lsst.meas.base as measBase 

29import lsst.afw.table as afwTable 

30import lsst.afw.detection as afwDet 

31import lsst.geom as geom 

32import lsst.pex.exceptions as pexExcept 

33import lsst.pex.config as pexConfig 

34from lsst.pipe.base import Struct 

35from lsst.utils.timer import timeMethod 

36 

37__all__ = ("DipoleFitTask", "DipoleFitPlugin", "DipoleFitTaskConfig", "DipoleFitPluginConfig", 

38 "DipoleFitAlgorithm") 

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 

45class 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 fitSeparateNegParams = pexConfig.Field( 

77 dtype=bool, default=False, 

78 doc="Include parameters to fit for negative values (flux, gradient) separately from pos.") 

79 

80 # Config params for classification of detected diaSources as dipole or not 

81 minSn = pexConfig.Field( 

82 dtype=float, default=np.sqrt(2) * 5.0, 

83 doc="Minimum quadrature sum of positive+negative lobe S/N to be considered a dipole") 

84 

85 maxFluxRatio = pexConfig.Field( 

86 dtype=float, default=0.65, 

87 doc="Maximum flux ratio in either lobe to be considered a dipole") 

88 

89 maxChi2DoF = pexConfig.Field( 

90 dtype=float, default=0.05, 

91 doc="""Maximum Chi2/DoF significance of fit to be considered a dipole. 

92 Default value means \"Choose a chi2DoF corresponding to a significance level of at most 0.05\" 

93 (note this is actually a significance, not a chi2 value).""") 

94 

95 

96class DipoleFitTaskConfig(measBase.SingleFrameMeasurementConfig): 

97 """Measurement of detected diaSources as dipoles 

98 

99 Currently we keep the "old" DipoleMeasurement algorithms turned on. 

100 """ 

101 

102 def setDefaults(self): 

103 measBase.SingleFrameMeasurementConfig.setDefaults(self) 

104 

105 self.plugins.names = ["base_CircularApertureFlux", 

106 "base_PixelFlags", 

107 "base_SkyCoord", 

108 "base_PsfFlux", 

109 "base_SdssCentroid", 

110 "base_SdssShape", 

111 "base_GaussianFlux", 

112 "base_PeakLikelihoodFlux", 

113 "base_PeakCentroid", 

114 "base_NaiveCentroid", 

115 "ip_diffim_NaiveDipoleCentroid", 

116 "ip_diffim_NaiveDipoleFlux", 

117 "ip_diffim_PsfDipoleFlux", 

118 "ip_diffim_ClassificationDipole", 

119 ] 

120 

121 self.slots.calibFlux = None 

122 self.slots.modelFlux = None 

123 self.slots.gaussianFlux = None 

124 self.slots.shape = "base_SdssShape" 

125 self.slots.centroid = "ip_diffim_NaiveDipoleCentroid" 

126 self.doReplaceWithNoise = False 

127 

128 

129class DipoleFitTask(measBase.SingleFrameMeasurementTask): 

130 """A task that fits a dipole to a difference image, with an optional separate detection image. 

131 

132 Because it subclasses SingleFrameMeasurementTask, and calls 

133 SingleFrameMeasurementTask.run() from its run() method, it still 

134 can be used identically to a standard SingleFrameMeasurementTask. 

135 """ 

136 

137 ConfigClass = DipoleFitTaskConfig 

138 _DefaultName = "ip_diffim_DipoleFit" 

139 

140 def __init__(self, schema, algMetadata=None, **kwargs): 

141 

142 measBase.SingleFrameMeasurementTask.__init__(self, schema, algMetadata, **kwargs) 

143 

144 dpFitPluginConfig = self.config.plugins['ip_diffim_DipoleFit'] 

145 

146 self.dipoleFitter = DipoleFitPlugin(dpFitPluginConfig, name=self._DefaultName, 

147 schema=schema, metadata=algMetadata, 

148 logName=self.log.name) 

149 

150 @timeMethod 

151 def run(self, sources, exposure, posExp=None, negExp=None, **kwargs): 

152 """Run dipole measurement and classification 

153 

154 Parameters 

155 ---------- 

156 sources : `lsst.afw.table.SourceCatalog` 

157 ``diaSources`` that will be measured using dipole measurement 

158 exposure : `lsst.afw.image.Exposure` 

159 The difference exposure on which the ``diaSources`` of the ``sources`` parameter 

160 were detected. If neither ``posExp`` nor ``negExp`` are set, then the dipole is also 

161 fitted directly to this difference image. 

162 posExp : `lsst.afw.image.Exposure`, optional 

163 "Positive" exposure, typically a science exposure, or None if unavailable 

164 When `posExp` is `None`, will compute `posImage = exposure + negExp`. 

165 negExp : `lsst.afw.image.Exposure`, optional 

166 "Negative" exposure, typically a template exposure, or None if unavailable 

167 When `negExp` is `None`, will compute `negImage = posExp - exposure`. 

168 **kwargs 

169 Additional keyword arguments for `lsst.meas.base.sfm.SingleFrameMeasurementTask`. 

170 """ 

171 

172 measBase.SingleFrameMeasurementTask.run(self, sources, exposure, **kwargs) 

173 

174 if not sources: 

175 return 

176 

177 for source in sources: 

178 self.dipoleFitter.measure(source, exposure, posExp, negExp) 

179 

180 

181class DipoleModel: 

182 """Lightweight class containing methods for generating a dipole model for fitting 

183 to sources in diffims, used by DipoleFitAlgorithm. 

184 

185 See also: 

186 `DMTN-007: Dipole characterization for image differencing <https://dmtn-007.lsst.io>`_. 

187 """ 

188 

189 def __init__(self): 

190 import lsstDebug 

191 self.debug = lsstDebug.Info(__name__).debug 

192 self.log = logging.getLogger(__name__) 

193 

194 def makeBackgroundModel(self, in_x, pars=None): 

195 """Generate gradient model (2-d array) with up to 2nd-order polynomial 

196 

197 Parameters 

198 ---------- 

199 in_x : `numpy.array` 

200 (2, w, h)-dimensional `numpy.array`, containing the 

201 input x,y meshgrid providing the coordinates upon which to 

202 compute the gradient. This will typically be generated via 

203 `_generateXYGrid()`. `w` and `h` correspond to the width and 

204 height of the desired grid. 

205 pars : `list` of `float`, optional 

206 Up to 6 floats for up 

207 to 6 2nd-order 2-d polynomial gradient parameters, in the 

208 following order: (intercept, x, y, xy, x**2, y**2). If `pars` 

209 is emtpy or `None`, do nothing and return `None` (for speed). 

210 

211 Returns 

212 ------- 

213 result : `None` or `numpy.array` 

214 return None, or 2-d numpy.array of width/height matching 

215 input bbox, containing computed gradient values. 

216 """ 

217 

218 # Don't fit for other gradient parameters if the intercept is not included. 

219 if (pars is None) or (len(pars) <= 0) or (pars[0] is None): 

220 return 

221 

222 y, x = in_x[0, :], in_x[1, :] 

223 gradient = np.full_like(x, pars[0], dtype='float64') 

224 if len(pars) > 1 and pars[1] is not None: 

225 gradient += pars[1] * x 

226 if len(pars) > 2 and pars[2] is not None: 

227 gradient += pars[2] * y 

228 if len(pars) > 3 and pars[3] is not None: 

229 gradient += pars[3] * (x * y) 

230 if len(pars) > 4 and pars[4] is not None: 

231 gradient += pars[4] * (x * x) 

232 if len(pars) > 5 and pars[5] is not None: 

233 gradient += pars[5] * (y * y) 

234 

235 return gradient 

236 

237 def _generateXYGrid(self, bbox): 

238 """Generate a meshgrid covering the x,y coordinates bounded by bbox 

239 

240 Parameters 

241 ---------- 

242 bbox : `lsst.geom.Box2I` 

243 input Bounding Box defining the coordinate limits 

244 

245 Returns 

246 ------- 

247 in_x : `numpy.array` 

248 (2, w, h)-dimensional numpy array containing the grid indexing over x- and 

249 y- coordinates 

250 """ 

251 

252 x, y = np.mgrid[bbox.getBeginY():bbox.getEndY(), bbox.getBeginX():bbox.getEndX()] 

253 in_x = np.array([y, x]).astype(np.float64) 

254 in_x[0, :] -= np.mean(in_x[0, :]) 

255 in_x[1, :] -= np.mean(in_x[1, :]) 

256 return in_x 

257 

258 def _getHeavyFootprintSubimage(self, fp, badfill=np.nan, grow=0): 

259 """Extract the image from a ``~lsst.afw.detection.HeavyFootprint`` 

260 as an `lsst.afw.image.ImageF`. 

261 

262 Parameters 

263 ---------- 

264 fp : `lsst.afw.detection.HeavyFootprint` 

265 HeavyFootprint to use to generate the subimage 

266 badfill : `float`, optional 

267 Value to fill in pixels in extracted image that are outside the footprint 

268 grow : `int` 

269 Optionally grow the footprint by this amount before extraction 

270 

271 Returns 

272 ------- 

273 subim2 : `lsst.afw.image.ImageF` 

274 An `~lsst.afw.image.ImageF` containing the subimage. 

275 """ 

276 bbox = fp.getBBox() 

277 if grow > 0: 

278 bbox.grow(grow) 

279 

280 subim2 = afwImage.ImageF(bbox, badfill) 

281 fp.getSpans().unflatten(subim2.getArray(), fp.getImageArray(), bbox.getCorners()[0]) 

282 return subim2 

283 

284 def fitFootprintBackground(self, source, posImage, order=1): 

285 """Fit a linear (polynomial) model of given order (max 2) to the background of a footprint. 

286 

287 Only fit the pixels OUTSIDE of the footprint, but within its bounding box. 

288 

289 Parameters 

290 ---------- 

291 source : `lsst.afw.table.SourceRecord` 

292 SourceRecord, the footprint of which is to be fit 

293 posImage : `lsst.afw.image.Exposure` 

294 The exposure from which to extract the footprint subimage 

295 order : `int` 

296 Polynomial order of background gradient to fit. 

297 

298 Returns 

299 ------- 

300 pars : `tuple` of `float` 

301 `tuple` of length (1 if order==0; 3 if order==1; 6 if order == 2), 

302 containing the resulting fit parameters 

303 """ 

304 

305 # TODO look into whether to use afwMath background methods -- see 

306 # http://lsst-web.ncsa.illinois.edu/doxygen/x_masterDoxyDoc/_background_example.html 

307 fp = source.getFootprint() 

308 bbox = fp.getBBox() 

309 bbox.grow(3) 

310 posImg = afwImage.ImageF(posImage.getMaskedImage().getImage(), bbox, afwImage.PARENT) 

311 

312 # This code constructs the footprint image so that we can identify the pixels that are 

313 # outside the footprint (but within the bounding box). These are the pixels used for 

314 # fitting the background. 

315 posHfp = afwDet.HeavyFootprintF(fp, posImage.getMaskedImage()) 

316 posFpImg = self._getHeavyFootprintSubimage(posHfp, grow=3) 

317 

318 isBg = np.isnan(posFpImg.getArray()).ravel() 

319 

320 data = posImg.getArray().ravel() 

321 data = data[isBg] 

322 B = data 

323 

324 x, y = np.mgrid[bbox.getBeginY():bbox.getEndY(), bbox.getBeginX():bbox.getEndX()] 

325 x = x.astype(np.float64).ravel() 

326 x -= np.mean(x) 

327 x = x[isBg] 

328 y = y.astype(np.float64).ravel() 

329 y -= np.mean(y) 

330 y = y[isBg] 

331 b = np.ones_like(x, dtype=np.float64) 

332 

333 M = np.vstack([b]).T # order = 0 

334 if order == 1: 

335 M = np.vstack([b, x, y]).T 

336 elif order == 2: 

337 M = np.vstack([b, x, y, x**2., y**2., x*y]).T 

338 

339 pars = np.linalg.lstsq(M, B, rcond=-1)[0] 

340 return pars 

341 

342 def makeStarModel(self, bbox, psf, xcen, ycen, flux): 

343 """Generate a 2D image model of a single PDF centered at the given coordinates. 

344 

345 Parameters 

346 ---------- 

347 bbox : `lsst.geom.Box` 

348 Bounding box marking pixel coordinates for generated model 

349 psf : TODO: DM-17458 

350 Psf model used to generate the 'star' 

351 xcen : `float` 

352 Desired x-centroid of the 'star' 

353 ycen : `float` 

354 Desired y-centroid of the 'star' 

355 flux : `float` 

356 Desired flux of the 'star' 

357 

358 Returns 

359 ------- 

360 p_Im : `lsst.afw.image.Image` 

361 2-d stellar image of width/height matching input ``bbox``, 

362 containing PSF with given centroid and flux 

363 """ 

364 

365 # Generate the psf image, normalize to flux 

366 psf_img = psf.computeImage(geom.Point2D(xcen, ycen)).convertF() 

367 psf_img_sum = np.nansum(psf_img.getArray()) 

368 psf_img *= (flux/psf_img_sum) 

369 

370 # Clip the PSF image bounding box to fall within the footprint bounding box 

371 psf_box = psf_img.getBBox() 

372 psf_box.clip(bbox) 

373 psf_img = afwImage.ImageF(psf_img, psf_box, afwImage.PARENT) 

374 

375 # Then actually crop the psf image. 

376 # Usually not necessary, but if the dipole is near the edge of the image... 

377 # Would be nice if we could compare original pos_box with clipped pos_box and 

378 # see if it actually was clipped. 

379 p_Im = afwImage.ImageF(bbox) 

380 tmpSubim = afwImage.ImageF(p_Im, psf_box, afwImage.PARENT) 

381 tmpSubim += psf_img 

382 

383 return p_Im 

384 

385 def makeModel(self, x, flux, xcenPos, ycenPos, xcenNeg, ycenNeg, fluxNeg=None, 

386 b=None, x1=None, y1=None, xy=None, x2=None, y2=None, 

387 bNeg=None, x1Neg=None, y1Neg=None, xyNeg=None, x2Neg=None, y2Neg=None, 

388 **kwargs): 

389 """Generate dipole model with given parameters. 

390 

391 This is the function whose sum-of-squared difference from data 

392 is minimized by `lmfit`. 

393 

394 x : TODO: DM-17458 

395 Input independent variable. Used here as the grid on 

396 which to compute the background gradient model. 

397 flux : `float` 

398 Desired flux of the positive lobe of the dipole 

399 xcenPos, ycenPos : `float` 

400 Desired x,y-centroid of the positive lobe of the dipole 

401 xcenNeg, ycenNeg : `float` 

402 Desired x,y-centroid of the negative lobe of the dipole 

403 fluxNeg : `float`, optional 

404 Desired flux of the negative lobe of the dipole, set to 'flux' if None 

405 b, x1, y1, xy, x2, y2 : `float` 

406 Gradient parameters for positive lobe. 

407 bNeg, x1Neg, y1Neg, xyNeg, x2Neg, y2Neg : `float`, optional 

408 Gradient parameters for negative lobe. 

409 They are set to the corresponding positive values if None. 

410 

411 **kwargs : `dict` [`str`] 

412 Keyword arguments passed through ``lmfit`` and 

413 used by this function. These must include: 

414 

415 - ``psf`` Psf model used to generate the 'star' 

416 - ``rel_weight`` Used to signify least-squares weighting of posImage/negImage 

417 relative to diffim. If ``rel_weight == 0`` then posImage/negImage are ignored. 

418 - ``bbox`` Bounding box containing region to be modelled 

419 

420 Returns 

421 ------- 

422 zout : `numpy.array` 

423 Has width and height matching the input bbox, and 

424 contains the dipole model with given centroids and flux(es). If 

425 ``rel_weight`` = 0, this is a 2-d array with dimensions matching 

426 those of bbox; otherwise a stack of three such arrays, 

427 representing the dipole (diffim), positive, and negative images 

428 respectively. 

429 """ 

430 

431 psf = kwargs.get('psf') 

432 rel_weight = kwargs.get('rel_weight') # if > 0, we're including pre-sub. images 

433 fp = kwargs.get('footprint') 

434 bbox = fp.getBBox() 

435 

436 if fluxNeg is None: 

437 fluxNeg = flux 

438 

439 if self.debug: 

440 self.log.debug('%.2f %.2f %.2f %.2f %.2f %.2f', 

441 flux, fluxNeg, xcenPos, ycenPos, xcenNeg, ycenNeg) 

442 if x1 is not None: 

443 self.log.debug(' %.2f %.2f %.2f', b, x1, y1) 

444 if xy is not None: 

445 self.log.debug(' %.2f %.2f %.2f', xy, x2, y2) 

446 

447 posIm = self.makeStarModel(bbox, psf, xcenPos, ycenPos, flux) 

448 negIm = self.makeStarModel(bbox, psf, xcenNeg, ycenNeg, fluxNeg) 

449 

450 in_x = x 

451 if in_x is None: # use the footprint to generate the input grid 

452 y, x = np.mgrid[bbox.getBeginY():bbox.getEndY(), bbox.getBeginX():bbox.getEndX()] 

453 in_x = np.array([x, y]) * 1. 

454 in_x[0, :] -= in_x[0, :].mean() # center it! 

455 in_x[1, :] -= in_x[1, :].mean() 

456 

457 if b is not None: 

458 gradient = self.makeBackgroundModel(in_x, (b, x1, y1, xy, x2, y2)) 

459 

460 # If bNeg is None, then don't fit the negative background separately 

461 if bNeg is not None: 

462 gradientNeg = self.makeBackgroundModel(in_x, (bNeg, x1Neg, y1Neg, xyNeg, x2Neg, y2Neg)) 

463 else: 

464 gradientNeg = gradient 

465 

466 posIm.getArray()[:, :] += gradient 

467 negIm.getArray()[:, :] += gradientNeg 

468 

469 # Generate the diffIm model 

470 diffIm = afwImage.ImageF(bbox) 

471 diffIm += posIm 

472 diffIm -= negIm 

473 

474 zout = diffIm.getArray() 

475 if rel_weight > 0.: 

476 zout = np.append([zout], [posIm.getArray(), negIm.getArray()], axis=0) 

477 

478 return zout 

479 

480 

481class DipoleFitAlgorithm: 

482 """Fit a dipole model using an image difference. 

483 

484 See also: 

485 `DMTN-007: Dipole characterization for image differencing <https://dmtn-007.lsst.io>`_. 

486 """ 

487 

488 # This is just a private version number to sync with the ipython notebooks that I have been 

489 # using for algorithm development. 

490 _private_version_ = '0.0.5' 

491 

492 # Below is a (somewhat incomplete) list of improvements 

493 # that would be worth investigating, given the time: 

494 

495 # todo 1. evaluate necessity for separate parameters for pos- and neg- images 

496 # todo 2. only fit background OUTSIDE footprint (DONE) and dipole params INSIDE footprint (NOT DONE)? 

497 # todo 3. correct normalization of least-squares weights based on variance planes 

498 # todo 4. account for PSFs that vary across the exposures (should be happening by default?) 

499 # todo 5. correctly account for NA/masks (i.e., ignore!) 

500 # todo 6. better exception handling in the plugin 

501 # todo 7. better classification of dipoles (e.g. by comparing chi2 fit vs. monopole?) 

502 # todo 8. (DONE) Initial fast estimate of background gradient(s) params -- perhaps using numpy.lstsq 

503 # todo 9. (NOT NEEDED - see (2)) Initial fast test whether a background gradient needs to be fit 

504 # todo 10. (DONE) better initial estimate for flux when there's a strong gradient 

505 # todo 11. (DONE) requires a new package `lmfit` -- investiate others? (astropy/scipy/iminuit?) 

506 

507 def __init__(self, diffim, posImage=None, negImage=None): 

508 """Algorithm to run dipole measurement on a diaSource 

509 

510 Parameters 

511 ---------- 

512 diffim : `lsst.afw.image.Exposure` 

513 Exposure on which the diaSources were detected 

514 posImage : `lsst.afw.image.Exposure` 

515 "Positive" exposure from which the template was subtracted 

516 negImage : `lsst.afw.image.Exposure` 

517 "Negative" exposure which was subtracted from the posImage 

518 """ 

519 

520 self.diffim = diffim 

521 self.posImage = posImage 

522 self.negImage = negImage 

523 self.psfSigma = None 

524 if diffim is not None: 

525 diffimPsf = diffim.getPsf() 

526 diffimAvgPos = diffimPsf.getAveragePosition() 

527 self.psfSigma = diffimPsf.computeShape(diffimAvgPos).getDeterminantRadius() 

528 

529 self.log = logging.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 

624 # Add the constraints for centroids, fluxes. 

625 # starting constraint - near centroid of footprint 

626 fpCentroid = np.array([fp.getCentroid().getX(), fp.getCentroid().getY()]) 

627 cenNeg = cenPos = fpCentroid 

628 

629 pks = fp.getPeaks() 

630 

631 if len(pks) >= 1: 

632 cenPos = pks[0].getF() # if individual (merged) peaks were detected, use those 

633 if len(pks) >= 2: # peaks are already sorted by centroid flux so take the most negative one 

634 cenNeg = pks[-1].getF() 

635 

636 # For close/faint dipoles the starting locs (min/max) might be way off, let's help them a bit. 

637 # First assume dipole is not separated by more than 5*psfSigma. 

638 maxSep = self.psfSigma * maxSepInSigma 

639 

640 # As an initial guess -- assume the dipole is close to the center of the footprint. 

641 if np.sum(np.sqrt((np.array(cenPos) - fpCentroid)**2.)) > maxSep: 

642 cenPos = fpCentroid 

643 if np.sum(np.sqrt((np.array(cenNeg) - fpCentroid)**2.)) > maxSep: 

644 cenPos = fpCentroid 

645 

646 # parameter hints/constraints: https://lmfit.github.io/lmfit-py/model.html#model-param-hints-section 

647 # might make sense to not use bounds -- see http://lmfit.github.io/lmfit-py/bounds.html 

648 # also see this discussion -- https://github.com/scipy/scipy/issues/3129 

649 gmod.set_param_hint('xcenPos', value=cenPos[0], 

650 min=cenPos[0]-maxSep, max=cenPos[0]+maxSep) 

651 gmod.set_param_hint('ycenPos', value=cenPos[1], 

652 min=cenPos[1]-maxSep, max=cenPos[1]+maxSep) 

653 gmod.set_param_hint('xcenNeg', value=cenNeg[0], 

654 min=cenNeg[0]-maxSep, max=cenNeg[0]+maxSep) 

655 gmod.set_param_hint('ycenNeg', value=cenNeg[1], 

656 min=cenNeg[1]-maxSep, max=cenNeg[1]+maxSep) 

657 

658 # Use the (flux under the dipole)*5 for an estimate. 

659 # Lots of testing showed that having startingFlux be too high was better than too low. 

660 startingFlux = np.nansum(np.abs(diArr) - np.nanmedian(np.abs(diArr))) * 5. 

661 posFlux = negFlux = startingFlux 

662 

663 # TBD: set max. flux limit? 

664 gmod.set_param_hint('flux', value=posFlux, min=0.1) 

665 

666 if separateNegParams: 

667 # TBD: set max negative lobe flux limit? 

668 gmod.set_param_hint('fluxNeg', value=np.abs(negFlux), min=0.1) 

669 

670 # Fixed parameters (don't fit for them if there are no pre-sub images or no gradient fit requested): 

671 # Right now (fitBackground == 1), we fit a linear model to the background and then subtract 

672 # it from the data and then don't fit the background again (this is faster). 

673 # A slower alternative (fitBackground == 2) is to use the estimated background parameters as 

674 # starting points in the integrated model fit. That is currently not performed by default, 

675 # but might be desirable in some cases. 

676 bgParsPos = bgParsNeg = (0., 0., 0.) 

677 if ((rel_weight > 0.) and (fitBackground != 0) and (bgGradientOrder >= 0)): 

678 pbg = 0. 

679 bgFitImage = self.posImage if self.posImage is not None else self.negImage 

680 # Fit the gradient to the background (linear model) 

681 bgParsPos = bgParsNeg = dipoleModel.fitFootprintBackground(source, bgFitImage, 

682 order=bgGradientOrder) 

683 

684 # Generate the gradient and subtract it from the pre-subtraction image data 

685 if fitBackground == 1: 

686 in_x = dipoleModel._generateXYGrid(bbox) 

687 pbg = dipoleModel.makeBackgroundModel(in_x, tuple(bgParsPos)) 

688 z[1, :] -= pbg 

689 z[1, :] -= np.nanmedian(z[1, :]) 

690 posFlux = np.nansum(z[1, :]) 

691 gmod.set_param_hint('flux', value=posFlux*1.5, min=0.1) 

692 

693 if separateNegParams and self.negImage is not None: 

694 bgParsNeg = dipoleModel.fitFootprintBackground(source, self.negImage, 

695 order=bgGradientOrder) 

696 pbg = dipoleModel.makeBackgroundModel(in_x, tuple(bgParsNeg)) 

697 z[2, :] -= pbg 

698 z[2, :] -= np.nanmedian(z[2, :]) 

699 if separateNegParams: 

700 negFlux = np.nansum(z[2, :]) 

701 gmod.set_param_hint('fluxNeg', value=negFlux*1.5, min=0.1) 

702 

703 # Do not subtract the background from the images but include the background parameters in the fit 

704 if fitBackground == 2: 

705 if bgGradientOrder >= 0: 

706 gmod.set_param_hint('b', value=bgParsPos[0]) 

707 if separateNegParams: 

708 gmod.set_param_hint('bNeg', value=bgParsNeg[0]) 

709 if bgGradientOrder >= 1: 

710 gmod.set_param_hint('x1', value=bgParsPos[1]) 

711 gmod.set_param_hint('y1', value=bgParsPos[2]) 

712 if separateNegParams: 

713 gmod.set_param_hint('x1Neg', value=bgParsNeg[1]) 

714 gmod.set_param_hint('y1Neg', value=bgParsNeg[2]) 

715 if bgGradientOrder >= 2: 

716 gmod.set_param_hint('xy', value=bgParsPos[3]) 

717 gmod.set_param_hint('x2', value=bgParsPos[4]) 

718 gmod.set_param_hint('y2', value=bgParsPos[5]) 

719 if separateNegParams: 

720 gmod.set_param_hint('xyNeg', value=bgParsNeg[3]) 

721 gmod.set_param_hint('x2Neg', value=bgParsNeg[4]) 

722 gmod.set_param_hint('y2Neg', value=bgParsNeg[5]) 

723 

724 y, x = np.mgrid[bbox.getBeginY():bbox.getEndY(), bbox.getBeginX():bbox.getEndX()] 

725 in_x = np.array([x, y]).astype(np.float64) 

726 in_x[0, :] -= in_x[0, :].mean() # center it! 

727 in_x[1, :] -= in_x[1, :].mean() 

728 

729 # Instead of explicitly using a mask to ignore flagged pixels, just set the ignored pixels' 

730 # weights to 0 in the fit. TBD: need to inspect mask planes to set this mask. 

731 mask = np.ones_like(z, dtype=bool) # TBD: set mask values to False if the pixels are to be ignored 

732 

733 # I'm not sure about the variance planes in the diffim (or convolved pre-sub. images 

734 # for that matter) so for now, let's just do an un-weighted least-squares fit 

735 # (override weights computed above). 

736 weights = mask.astype(np.float64) 

737 if self.posImage is not None and rel_weight > 0.: 

738 weights = np.array([np.ones_like(diArr), np.ones_like(diArr)*rel_weight, 

739 np.ones_like(diArr)*rel_weight]) 

740 

741 # Set the weights to zero if mask is False 

742 if np.any(~mask): 

743 weights[~mask] = 0. 

744 

745 # Note that although we can, we're not required to set initial values for params here, 

746 # since we set their param_hint's above. 

747 # Can add "method" param to not use 'leastsq' (==levenberg-marquardt), e.g. "method='nelder'" 

748 with warnings.catch_warnings(): 

749 # Ignore lmfit unknown argument warnings: 

750 # "psf, rel_weight, footprint, modelObj" all become pass-through kwargs for makeModel. 

751 warnings.filterwarnings("ignore", "The keyword argument .* does not match", UserWarning) 

752 result = gmod.fit(z, weights=weights, x=in_x, max_nfev=250, 

753 method="leastsq", # TODO: try using `least_squares` here for speed/robustness 

754 verbose=verbose, 

755 # see scipy docs for the meaning of these keywords 

756 fit_kws={'ftol': tol, 'xtol': tol, 'gtol': tol, 

757 # Our model is float32 internally, so we need a larger epsfcn. 

758 'epsfcn': 1e-10}, 

759 psf=self.diffim.getPsf(), # hereon: kwargs that get passed to makeModel() 

760 rel_weight=rel_weight, 

761 footprint=fp, 

762 modelObj=dipoleModel) 

763 

764 if verbose: # the ci_report() seems to fail if neg params are constrained -- TBD why. 

765 # Never wanted in production - this takes a long time (longer than the fit!) 

766 # This is how to get confidence intervals out: 

767 # https://lmfit.github.io/lmfit-py/confidence.html and 

768 # http://cars9.uchicago.edu/software/python/lmfit/model.html 

769 print(result.fit_report(show_correl=False)) 

770 if separateNegParams: 

771 print(result.ci_report()) 

772 

773 return result 

774 

775 def fitDipole(self, source, tol=1e-7, rel_weight=0.1, 

776 fitBackground=1, maxSepInSigma=5., separateNegParams=True, 

777 bgGradientOrder=1, verbose=False, display=False): 

778 """Fit a dipole model to an input ``diaSource`` (wraps `fitDipoleImpl`). 

779 

780 Actually, fits the subimage bounded by the input source's 

781 footprint) and optionally constrain the fit using the 

782 pre-subtraction images self.posImage (science) and 

783 self.negImage (template). Wraps the output into a 

784 `pipeBase.Struct` named tuple after computing additional 

785 statistics such as orientation and SNR. 

786 

787 Parameters 

788 ---------- 

789 source : `lsst.afw.table.SourceRecord` 

790 Record containing the (merged) dipole source footprint detected on the diffim 

791 tol : `float`, optional 

792 Tolerance parameter for scipy.leastsq() optimization 

793 rel_weight : `float`, optional 

794 Weighting of posImage/negImage relative to the diffim in the fit 

795 fitBackground : `int`, {0, 1, 2}, optional 

796 How to fit linear background gradient in posImage/negImage 

797 

798 - 0: do not fit background at all 

799 - 1 (default): pre-fit the background using linear least squares and then do not fit it 

800 as part of the dipole fitting optimization 

801 - 2: pre-fit the background using linear least squares (as in 1), and use the parameter 

802 estimates from that fit as starting parameters for an integrated "re-fit" of the 

803 background as part of the overall dipole fitting optimization. 

804 maxSepInSigma : `float`, optional 

805 Allowed window of centroid parameters relative to peak in input source footprint 

806 separateNegParams : `bool`, optional 

807 Fit separate parameters to the flux and background gradient in 

808 bgGradientOrder : `int`, {0, 1, 2}, optional 

809 Desired polynomial order of background gradient 

810 verbose: `bool`, optional 

811 Be verbose 

812 display 

813 Display input data, best fit model(s) and residuals in a matplotlib window. 

814 

815 Returns 

816 ------- 

817 result : `struct` 

818 `pipeBase.Struct` object containing the fit parameters and other information. 

819 

820 result : `callable` 

821 `lmfit.MinimizerResult` object for debugging and error estimation, etc. 

822 

823 Notes 

824 ----- 

825 Parameter `fitBackground` has three options, thus it is an integer: 

826 

827 """ 

828 

829 fitResult = self.fitDipoleImpl( 

830 source, tol=tol, rel_weight=rel_weight, fitBackground=fitBackground, 

831 maxSepInSigma=maxSepInSigma, separateNegParams=separateNegParams, 

832 bgGradientOrder=bgGradientOrder, verbose=verbose) 

833 

834 # Display images, model fits and residuals (currently uses matplotlib display functions) 

835 if display: 

836 fp = source.getFootprint() 

837 self.displayFitResults(fp, fitResult) 

838 

839 fitParams = fitResult.best_values 

840 if fitParams['flux'] <= 1.: # usually around 0.1 -- the minimum flux allowed -- i.e. bad fit. 

841 out = Struct(posCentroidX=np.nan, posCentroidY=np.nan, 

842 negCentroidX=np.nan, negCentroidY=np.nan, 

843 posFlux=np.nan, negFlux=np.nan, posFluxErr=np.nan, negFluxErr=np.nan, 

844 centroidX=np.nan, centroidY=np.nan, orientation=np.nan, 

845 signalToNoise=np.nan, chi2=np.nan, redChi2=np.nan) 

846 return out, fitResult 

847 

848 centroid = ((fitParams['xcenPos'] + fitParams['xcenNeg']) / 2., 

849 (fitParams['ycenPos'] + fitParams['ycenNeg']) / 2.) 

850 dx, dy = fitParams['xcenPos'] - fitParams['xcenNeg'], fitParams['ycenPos'] - fitParams['ycenNeg'] 

851 angle = np.arctan2(dy, dx) / np.pi * 180. # convert to degrees (should keep as rad?) 

852 

853 # Exctract flux value, compute signalToNoise from flux/variance_within_footprint 

854 # Also extract the stderr of flux estimate. 

855 def computeSumVariance(exposure, footprint): 

856 box = footprint.getBBox() 

857 subim = afwImage.MaskedImageF(exposure.getMaskedImage(), box, origin=afwImage.PARENT) 

858 return np.sqrt(np.nansum(subim.getArrays()[1][:, :])) 

859 

860 fluxVal = fluxVar = fitParams['flux'] 

861 fluxErr = fluxErrNeg = fitResult.params['flux'].stderr 

862 if self.posImage is not None: 

863 fluxVar = computeSumVariance(self.posImage, source.getFootprint()) 

864 else: 

865 fluxVar = computeSumVariance(self.diffim, source.getFootprint()) 

866 

867 fluxValNeg, fluxVarNeg = fluxVal, fluxVar 

868 if separateNegParams: 

869 fluxValNeg = fitParams['fluxNeg'] 

870 fluxErrNeg = fitResult.params['fluxNeg'].stderr 

871 if self.negImage is not None: 

872 fluxVarNeg = computeSumVariance(self.negImage, source.getFootprint()) 

873 

874 try: 

875 signalToNoise = np.sqrt((fluxVal/fluxVar)**2 + (fluxValNeg/fluxVarNeg)**2) 

876 except ZeroDivisionError: # catch divide by zero - should never happen. 

877 signalToNoise = np.nan 

878 

879 out = Struct(posCentroidX=fitParams['xcenPos'], posCentroidY=fitParams['ycenPos'], 

880 negCentroidX=fitParams['xcenNeg'], negCentroidY=fitParams['ycenNeg'], 

881 posFlux=fluxVal, negFlux=-fluxValNeg, posFluxErr=fluxErr, negFluxErr=fluxErrNeg, 

882 centroidX=centroid[0], centroidY=centroid[1], orientation=angle, 

883 signalToNoise=signalToNoise, chi2=fitResult.chisqr, redChi2=fitResult.redchi) 

884 

885 # fitResult may be returned for debugging 

886 return out, fitResult 

887 

888 def displayFitResults(self, footprint, result): 

889 """Display data, model fits and residuals (currently uses matplotlib display functions). 

890 

891 Parameters 

892 ---------- 

893 footprint : TODO: DM-17458 

894 Footprint containing the dipole that was fit 

895 result : `lmfit.MinimizerResult` 

896 `lmfit.MinimizerResult` object returned by `lmfit` optimizer 

897 

898 Returns 

899 ------- 

900 fig : `matplotlib.pyplot.plot` 

901 """ 

902 try: 

903 import matplotlib.pyplot as plt 

904 except ImportError as err: 

905 self.log.warning('Unable to import matplotlib: %s', err) 

906 raise err 

907 

908 def display2dArray(arr, title='Data', extent=None): 

909 """Use `matplotlib.pyplot.imshow` to display a 2-D array with a given coordinate range. 

910 """ 

911 fig = plt.imshow(arr, origin='lower', interpolation='none', cmap='gray', extent=extent) 

912 plt.title(title) 

913 plt.colorbar(fig, cmap='gray') 

914 return fig 

915 

916 z = result.data 

917 fit = result.best_fit 

918 bbox = footprint.getBBox() 

919 extent = (bbox.getBeginX(), bbox.getEndX(), bbox.getBeginY(), bbox.getEndY()) 

920 if z.shape[0] == 3: 

921 fig = plt.figure(figsize=(8, 8)) 

922 for i in range(3): 

923 plt.subplot(3, 3, i*3+1) 

924 display2dArray(z[i, :], 'Data', extent=extent) 

925 plt.subplot(3, 3, i*3+2) 

926 display2dArray(fit[i, :], 'Model', extent=extent) 

927 plt.subplot(3, 3, i*3+3) 

928 display2dArray(z[i, :] - fit[i, :], 'Residual', extent=extent) 

929 return fig 

930 else: 

931 fig = plt.figure(figsize=(8, 2.5)) 

932 plt.subplot(1, 3, 1) 

933 display2dArray(z, 'Data', extent=extent) 

934 plt.subplot(1, 3, 2) 

935 display2dArray(fit, 'Model', extent=extent) 

936 plt.subplot(1, 3, 3) 

937 display2dArray(z - fit, 'Residual', extent=extent) 

938 return fig 

939 

940 plt.show() 

941 

942 

943@measBase.register("ip_diffim_DipoleFit") 

944class DipoleFitPlugin(measBase.SingleFramePlugin): 

945 """A single frame measurement plugin that fits dipoles to all merged (two-peak) ``diaSources``. 

946 

947 This measurement plugin accepts up to three input images in 

948 its `measure` method. If these are provided, it includes data 

949 from the pre-subtraction posImage (science image) and optionally 

950 negImage (template image) to constrain the fit. The meat of the 

951 fitting routines are in the class `~lsst.module.name.DipoleFitAlgorithm`. 

952 

953 Notes 

954 ----- 

955 The motivation behind this plugin and the necessity for including more than 

956 one exposure are documented in DMTN-007 (http://dmtn-007.lsst.io). 

957 

958 This class is named `ip_diffim_DipoleFit` so that it may be used alongside 

959 the existing `ip_diffim_DipoleMeasurement` classes until such a time as those 

960 are deemed to be replaceable by this. 

961 """ 

962 

963 ConfigClass = DipoleFitPluginConfig 

964 DipoleFitAlgorithmClass = DipoleFitAlgorithm # Pointer to the class that performs the fit 

965 

966 FAILURE_EDGE = 1 # too close to the edge 

967 FAILURE_FIT = 2 # failure in the fitting 

968 FAILURE_NOT_DIPOLE = 4 # input source is not a putative dipole to begin with 

969 

970 @classmethod 

971 def getExecutionOrder(cls): 

972 """Set execution order to `FLUX_ORDER`. 

973 

974 This includes algorithms that require both `getShape()` and `getCentroid()`, 

975 in addition to a Footprint and its Peaks. 

976 """ 

977 return cls.FLUX_ORDER 

978 

979 def __init__(self, config, name, schema, metadata, logName=None): 

980 if logName is None: 

981 logName = name 

982 measBase.SingleFramePlugin.__init__(self, config, name, schema, metadata, logName=logName) 

983 

984 self.log = logging.getLogger(logName) 

985 

986 self._setupSchema(config, name, schema, metadata) 

987 

988 def _setupSchema(self, config, name, schema, metadata): 

989 # Get a FunctorKey that can quickly look up the "blessed" centroid value. 

990 self.centroidKey = afwTable.Point2DKey(schema["slot_Centroid"]) 

991 

992 # Add some fields for our outputs, and save their Keys. 

993 # Use setattr() to programmatically set the pos/neg named attributes to values, e.g. 

994 # self.posCentroidKeyX = 'ip_diffim_DipoleFit_pos_centroid_x' 

995 

996 for pos_neg in ['pos', 'neg']: 

997 

998 key = schema.addField( 

999 schema.join(name, pos_neg, "instFlux"), type=float, units="count", 

1000 doc="Dipole {0} lobe flux".format(pos_neg)) 

1001 setattr(self, ''.join((pos_neg, 'FluxKey')), key) 

1002 

1003 key = schema.addField( 

1004 schema.join(name, pos_neg, "instFluxErr"), type=float, units="count", 

1005 doc="1-sigma uncertainty for {0} dipole flux".format(pos_neg)) 

1006 setattr(self, ''.join((pos_neg, 'FluxErrKey')), key) 

1007 

1008 for x_y in ['x', 'y']: 

1009 key = schema.addField( 

1010 schema.join(name, pos_neg, "centroid", x_y), type=float, units="pixel", 

1011 doc="Dipole {0} lobe centroid".format(pos_neg)) 

1012 setattr(self, ''.join((pos_neg, 'CentroidKey', x_y.upper())), key) 

1013 

1014 for x_y in ['x', 'y']: 

1015 key = schema.addField( 

1016 schema.join(name, "centroid", x_y), type=float, units="pixel", 

1017 doc="Dipole centroid") 

1018 setattr(self, ''.join(('centroidKey', x_y.upper())), key) 

1019 

1020 self.fluxKey = schema.addField( 

1021 schema.join(name, "instFlux"), type=float, units="count", 

1022 doc="Dipole overall flux") 

1023 

1024 self.orientationKey = schema.addField( 

1025 schema.join(name, "orientation"), type=float, units="deg", 

1026 doc="Dipole orientation") 

1027 

1028 self.separationKey = schema.addField( 

1029 schema.join(name, "separation"), type=float, units="pixel", 

1030 doc="Pixel separation between positive and negative lobes of dipole") 

1031 

1032 self.chi2dofKey = schema.addField( 

1033 schema.join(name, "chi2dof"), type=float, 

1034 doc="Chi2 per degree of freedom of dipole fit") 

1035 

1036 self.signalToNoiseKey = schema.addField( 

1037 schema.join(name, "signalToNoise"), type=float, 

1038 doc="Estimated signal-to-noise of dipole fit") 

1039 

1040 self.classificationFlagKey = schema.addField( 

1041 schema.join(name, "flag", "classification"), type="Flag", 

1042 doc="Flag indicating diaSource is classified as a dipole") 

1043 

1044 self.classificationAttemptedFlagKey = schema.addField( 

1045 schema.join(name, "flag", "classificationAttempted"), type="Flag", 

1046 doc="Flag indicating diaSource was attempted to be classified as a dipole") 

1047 

1048 self.flagKey = schema.addField( 

1049 schema.join(name, "flag"), type="Flag", 

1050 doc="General failure flag for dipole fit") 

1051 

1052 self.edgeFlagKey = schema.addField( 

1053 schema.join(name, "flag", "edge"), type="Flag", 

1054 doc="Flag set when dipole is too close to edge of image") 

1055 

1056 def measure(self, measRecord, exposure, posExp=None, negExp=None): 

1057 """Perform the non-linear least squares minimization on the putative dipole source. 

1058 

1059 Parameters 

1060 ---------- 

1061 measRecord : `lsst.afw.table.SourceRecord` 

1062 diaSources that will be measured using dipole measurement 

1063 exposure : `lsst.afw.image.Exposure` 

1064 Difference exposure on which the diaSources were detected; `exposure = posExp-negExp` 

1065 If both `posExp` and `negExp` are `None`, will attempt to fit the 

1066 dipole to just the `exposure` with no constraint. 

1067 posExp : `lsst.afw.image.Exposure`, optional 

1068 "Positive" exposure, typically a science exposure, or None if unavailable 

1069 When `posExp` is `None`, will compute `posImage = exposure + negExp`. 

1070 negExp : `lsst.afw.image.Exposure`, optional 

1071 "Negative" exposure, typically a template exposure, or None if unavailable 

1072 When `negExp` is `None`, will compute `negImage = posExp - exposure`. 

1073 

1074 Notes 

1075 ----- 

1076 The main functionality of this routine was placed outside of 

1077 this plugin (into `DipoleFitAlgorithm.fitDipole()`) so that 

1078 `DipoleFitAlgorithm.fitDipole()` can be called separately for 

1079 testing (@see `tests/testDipoleFitter.py`) 

1080 

1081 Returns 

1082 ------- 

1083 result : TODO: DM-17458 

1084 TODO: DM-17458 

1085 """ 

1086 

1087 result = None 

1088 pks = measRecord.getFootprint().getPeaks() 

1089 

1090 # Check if the footprint consists of a putative dipole - else don't fit it. 

1091 if ( 

1092 (len(pks) <= 1) # one peak in the footprint - not a dipole 

1093 or (len(pks) > 1 and (np.sign(pks[0].getPeakValue()) 

1094 == np.sign(pks[-1].getPeakValue()))) # peaks are same sign - not a dipole 

1095 ): 

1096 measRecord.set(self.classificationFlagKey, False) 

1097 measRecord.set(self.classificationAttemptedFlagKey, False) 

1098 self.fail(measRecord, measBase.MeasurementError('not a dipole', self.FAILURE_NOT_DIPOLE)) 

1099 if not self.config.fitAllDiaSources: 

1100 return result 

1101 

1102 try: 

1103 alg = self.DipoleFitAlgorithmClass(exposure, posImage=posExp, negImage=negExp) 

1104 result, _ = alg.fitDipole( 

1105 measRecord, rel_weight=self.config.relWeight, 

1106 tol=self.config.tolerance, 

1107 maxSepInSigma=self.config.maxSeparation, 

1108 fitBackground=self.config.fitBackground, 

1109 separateNegParams=self.config.fitSeparateNegParams, 

1110 verbose=False, display=False) 

1111 except pexExcept.LengthError: 

1112 self.fail(measRecord, measBase.MeasurementError('edge failure', self.FAILURE_EDGE)) 

1113 except Exception as e: 

1114 self.fail(measRecord, measBase.MeasurementError('Exception in dipole fit', self.FAILURE_FIT)) 

1115 self.log.error("Exception in dipole fit. %s: %s", e.__class__.__name__, e) 

1116 

1117 if result is None: 

1118 measRecord.set(self.classificationFlagKey, False) 

1119 measRecord.set(self.classificationAttemptedFlagKey, False) 

1120 return result 

1121 

1122 self.log.debug("Dipole fit result: %d %s", measRecord.getId(), str(result)) 

1123 

1124 if result.posFlux <= 1.: # usually around 0.1 -- the minimum flux allowed -- i.e. bad fit. 

1125 self.fail(measRecord, measBase.MeasurementError('dipole fit failure', self.FAILURE_FIT)) 

1126 

1127 # add chi2, coord/flux uncertainties (TBD), dipole classification 

1128 # Add the relevant values to the measRecord 

1129 measRecord[self.posFluxKey] = result.posFlux 

1130 measRecord[self.posFluxErrKey] = result.signalToNoise # to be changed to actual sigma! 

1131 measRecord[self.posCentroidKeyX] = result.posCentroidX 

1132 measRecord[self.posCentroidKeyY] = result.posCentroidY 

1133 

1134 measRecord[self.negFluxKey] = result.negFlux 

1135 measRecord[self.negFluxErrKey] = result.signalToNoise # to be changed to actual sigma! 

1136 measRecord[self.negCentroidKeyX] = result.negCentroidX 

1137 measRecord[self.negCentroidKeyY] = result.negCentroidY 

1138 

1139 # Dia source flux: average of pos+neg 

1140 measRecord[self.fluxKey] = (abs(result.posFlux) + abs(result.negFlux))/2. 

1141 measRecord[self.orientationKey] = result.orientation 

1142 measRecord[self.separationKey] = np.sqrt((result.posCentroidX - result.negCentroidX)**2. 

1143 + (result.posCentroidY - result.negCentroidY)**2.) 

1144 measRecord[self.centroidKeyX] = result.centroidX 

1145 measRecord[self.centroidKeyY] = result.centroidY 

1146 

1147 measRecord[self.signalToNoiseKey] = result.signalToNoise 

1148 measRecord[self.chi2dofKey] = result.redChi2 

1149 

1150 self.doClassify(measRecord, result.chi2) 

1151 

1152 def doClassify(self, measRecord, chi2val): 

1153 """Classify a source as a dipole. 

1154 

1155 Parameters 

1156 ---------- 

1157 measRecord : TODO: DM-17458 

1158 TODO: DM-17458 

1159 chi2val : TODO: DM-17458 

1160 TODO: DM-17458 

1161 

1162 Notes 

1163 ----- 

1164 Sources are classified as dipoles, or not, according to three criteria: 

1165 

1166 1. Does the total signal-to-noise surpass the ``minSn``? 

1167 2. Are the pos/neg fluxes greater than 1.0 and no more than 0.65 (``maxFluxRatio``) 

1168 of the total flux? By default this will never happen since ``posFlux == negFlux``. 

1169 3. Is it a good fit (``chi2dof`` < 1)? (Currently not used.) 

1170 """ 

1171 

1172 # First, does the total signal-to-noise surpass the minSn? 

1173 passesSn = measRecord[self.signalToNoiseKey] > self.config.minSn 

1174 

1175 # Second, are the pos/neg fluxes greater than 1.0 and no more than 0.65 (param maxFluxRatio) 

1176 # of the total flux? By default this will never happen since posFlux = negFlux. 

1177 passesFluxPos = (abs(measRecord[self.posFluxKey]) 

1178 / (measRecord[self.fluxKey]*2.)) < self.config.maxFluxRatio 

1179 passesFluxPos &= (abs(measRecord[self.posFluxKey]) >= 1.0) 

1180 passesFluxNeg = (abs(measRecord[self.negFluxKey]) 

1181 / (measRecord[self.fluxKey]*2.)) < self.config.maxFluxRatio 

1182 passesFluxNeg &= (abs(measRecord[self.negFluxKey]) >= 1.0) 

1183 allPass = (passesSn and passesFluxPos and passesFluxNeg) # and passesChi2) 

1184 

1185 # Third, is it a good fit (chi2dof < 1)? 

1186 # Use scipy's chi2 cumulative distrib to estimate significance 

1187 # This doesn't really work since I don't trust the values in the variance plane (which 

1188 # affects the least-sq weights, which affects the resulting chi2). 

1189 # But I'm going to keep this here for future use. 

1190 if False: 

1191 from scipy.stats import chi2 

1192 ndof = chi2val / measRecord[self.chi2dofKey] 

1193 significance = chi2.cdf(chi2val, ndof) 

1194 passesChi2 = significance < self.config.maxChi2DoF 

1195 allPass = allPass and passesChi2 

1196 

1197 measRecord.set(self.classificationAttemptedFlagKey, True) 

1198 

1199 if allPass: # Note cannot pass `allPass` into the `measRecord.set()` call below...? 

1200 measRecord.set(self.classificationFlagKey, True) 

1201 else: 

1202 measRecord.set(self.classificationFlagKey, False) 

1203 

1204 def fail(self, measRecord, error=None): 

1205 """Catch failures and set the correct flags. 

1206 """ 

1207 

1208 measRecord.set(self.flagKey, True) 

1209 if error is not None: 

1210 if error.getFlagBit() == self.FAILURE_EDGE: 

1211 self.log.warning('DipoleFitPlugin not run on record %d: %s', measRecord.getId(), str(error)) 

1212 measRecord.set(self.edgeFlagKey, True) 

1213 if error.getFlagBit() == self.FAILURE_FIT: 

1214 self.log.warning('DipoleFitPlugin failed on record %d: %s', measRecord.getId(), str(error)) 

1215 measRecord.set(self.flagKey, True) 

1216 if error.getFlagBit() == self.FAILURE_NOT_DIPOLE: 

1217 self.log.debug('DipoleFitPlugin not run on record %d: %s', 

1218 measRecord.getId(), str(error)) 

1219 measRecord.set(self.classificationAttemptedFlagKey, False) 

1220 measRecord.set(self.flagKey, True) 

1221 else: 

1222 self.log.warning('DipoleFitPlugin failed on record %d', measRecord.getId())