Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1# 

2# LSST Data Management System 

3# Copyright 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 numpy as np 

24 

25import lsst.afw.image as afwImage 

26import lsst.afw.math as afwMath 

27import lsst.geom as geom 

28import lsst.log 

29import lsst.meas.algorithms as measAlg 

30import lsst.pex.config as pexConfig 

31import lsst.pipe.base as pipeBase 

32 

33 

34from .imageMapReduce import (ImageMapReduceConfig, ImageMapReduceTask, 

35 ImageMapper) 

36 

37__all__ = ("DecorrelateALKernelTask", "DecorrelateALKernelConfig", 

38 "DecorrelateALKernelMapper", "DecorrelateALKernelMapReduceConfig", 

39 "DecorrelateALKernelSpatialConfig", "DecorrelateALKernelSpatialTask") 

40 

41 

42class DecorrelateALKernelConfig(pexConfig.Config): 

43 """Configuration parameters for the DecorrelateALKernelTask 

44 """ 

45 

46 ignoreMaskPlanes = pexConfig.ListField( 

47 dtype=str, 

48 doc="""Mask planes to ignore for sigma-clipped statistics""", 

49 default=("INTRP", "EDGE", "DETECTED", "SAT", "CR", "BAD", "NO_DATA", "DETECTED_NEGATIVE") 

50 ) 

51 

52 

53class DecorrelateALKernelTask(pipeBase.Task): 

54 """Decorrelate the effect of convolution by Alard-Lupton matching kernel in image difference 

55 

56 Notes 

57 ----- 

58 

59 Pipe-task that removes the neighboring-pixel covariance in an 

60 image difference that are added when the template image is 

61 convolved with the Alard-Lupton PSF matching kernel. 

62 

63 The image differencing pipeline task @link 

64 ip.diffim.psfMatch.PsfMatchTask PSFMatchTask@endlink and @link 

65 ip.diffim.psfMatch.PsfMatchConfigAL PSFMatchConfigAL@endlink uses 

66 the Alard and Lupton (1998) method for matching the PSFs of the 

67 template and science exposures prior to subtraction. The 

68 Alard-Lupton method identifies a matching kernel, which is then 

69 (typically) convolved with the template image to perform PSF 

70 matching. This convolution has the effect of adding covariance 

71 between neighboring pixels in the template image, which is then 

72 added to the image difference by subtraction. 

73 

74 The pixel covariance may be corrected by whitening the noise of 

75 the image difference. This task performs such a decorrelation by 

76 computing a decorrelation kernel (based upon the A&L matching 

77 kernel and variances in the template and science images) and 

78 convolving the image difference with it. This process is described 

79 in detail in [DMTN-021](http://dmtn-021.lsst.io). 

80 

81 This task has no standalone example, however it is applied as a 

82 subtask of pipe.tasks.imageDifference.ImageDifferenceTask. 

83 """ 

84 ConfigClass = DecorrelateALKernelConfig 

85 _DefaultName = "ip_diffim_decorrelateALKernel" 

86 

87 def __init__(self, *args, **kwargs): 

88 """Create the image decorrelation Task 

89 

90 Parameters 

91 ---------- 

92 args : 

93 arguments to be passed to ``lsst.pipe.base.task.Task.__init__`` 

94 kwargs : 

95 keyword arguments to be passed to ``lsst.pipe.base.task.Task.__init__`` 

96 """ 

97 pipeBase.Task.__init__(self, *args, **kwargs) 

98 

99 self.statsControl = afwMath.StatisticsControl() 

100 self.statsControl.setNumSigmaClip(3.) 

101 self.statsControl.setNumIter(3) 

102 self.statsControl.setAndMask(afwImage.Mask.getPlaneBitMask(self.config.ignoreMaskPlanes)) 

103 

104 def computeVarianceMean(self, exposure): 

105 statObj = afwMath.makeStatistics(exposure.getMaskedImage().getVariance(), 

106 exposure.getMaskedImage().getMask(), 

107 afwMath.MEANCLIP, self.statsControl) 

108 var = statObj.getValue(afwMath.MEANCLIP) 

109 return var 

110 

111 @pipeBase.timeMethod 

112 def run(self, scienceExposure, templateExposure, subtractedExposure, psfMatchingKernel, 

113 preConvKernel=None, xcen=None, ycen=None, svar=None, tvar=None, templateMatched=True): 

114 """Perform decorrelation of an image difference exposure. 

115 

116 Decorrelates the diffim due to the convolution of the templateExposure with the 

117 A&L PSF matching kernel. Currently can accept a spatially varying matching kernel but in 

118 this case it simply uses a static kernel from the center of the exposure. The decorrelation 

119 is described in [DMTN-021, Equation 1](http://dmtn-021.lsst.io/#equation-1), where 

120 `exposure` is I_1; templateExposure is I_2; `subtractedExposure` is D(k); 

121 `psfMatchingKernel` is kappa; and svar and tvar are their respective 

122 variances (see below). 

123 

124 Parameters 

125 ---------- 

126 scienceExposure : `lsst.afw.image.Exposure` 

127 The original science exposure (before `preConvKernel` applied). 

128 templateExposure : `lsst.afw.image.Exposure` 

129 The original template exposure warped into the science exposure dimensions. 

130 subtractedExposure : `lsst.afw.iamge.Exposure` 

131 the subtracted exposure produced by 

132 `ip_diffim.ImagePsfMatchTask.subtractExposures()`. The `subtractedExposure` must 

133 inherit its PSF from `exposure`, see notes below. 

134 psfMatchingKernel : `lsst.afw.detection.Psf` 

135 An (optionally spatially-varying) PSF matching kernel produced 

136 by `ip_diffim.ImagePsfMatchTask.subtractExposures()`. 

137 preConvKernel : `lsst.afw.math.Kernel`, optional 

138 if not None, then the `scienceExposure` was pre-convolved with this kernel. 

139 Allowed only if ``templateMatched==True``. 

140 xcen : `float`, optional 

141 X-pixel coordinate to use for computing constant matching kernel to use 

142 If `None` (default), then use the center of the image. 

143 ycen : `float`, optional 

144 Y-pixel coordinate to use for computing constant matching kernel to use 

145 If `None` (default), then use the center of the image. 

146 svar : `float`, optional 

147 Image variance for science image 

148 If `None` (default) then compute the variance over the entire input science image. 

149 tvar : `float`, optional 

150 Image variance for template image 

151 If `None` (default) then compute the variance over the entire input template image. 

152 templateMatched : `bool`, optional 

153 If True, the template exposure was matched (convolved) to the science exposure. 

154 See also notes below. 

155 

156 Returns 

157 ------- 

158 result : `lsst.pipe.base.Struct` 

159 - ``correctedExposure`` : the decorrelated diffim 

160 

161 Notes 

162 ----- 

163 The `subtractedExposure` is NOT updated. The returned `correctedExposure` has an updated but 

164 spatially fixed PSF. It is calculated as the center of image PSF corrected by the center of 

165 image matching kernel. 

166 

167 If ``templateMatched==True``, the templateExposure was matched (convolved) 

168 to the ``scienceExposure`` by ``psfMatchingKernel``. Otherwise the ``scienceExposure`` 

169 was matched (convolved) by ``psfMatchingKernel``. 

170 

171 This task discards the variance plane of ``subtractedExposure`` and re-computes 

172 it from the variance planes of ``scienceExposure`` and ``templateExposure``. 

173 The image plane of ``subtractedExposure`` must be at the photometric level 

174 set by the AL PSF matching in `ImagePsfMatchTask.subtractExposures`. 

175 The assumptions about the photometric level are controlled by the 

176 `templateMatched` option in this task. 

177 

178 Here we currently convert a spatially-varying matching kernel into a constant kernel, 

179 just by computing it at the center of the image (tickets DM-6243, DM-6244). 

180 

181 We are also using a constant accross-the-image measure of sigma (sqrt(variance)) to compute 

182 the decorrelation kernel. 

183 

184 TODO DM-23857 As part of the spatially varying correction implementation 

185 consider whether returning a Struct is still necessary. 

186 """ 

187 if preConvKernel is not None and not templateMatched: 

188 raise ValueError("Pre-convolution and the matching of the " 

189 "science exposure is not supported.") 

190 

191 spatialKernel = psfMatchingKernel 

192 kimg = afwImage.ImageD(spatialKernel.getDimensions()) 

193 bbox = subtractedExposure.getBBox() 

194 if xcen is None: 

195 xcen = (bbox.getBeginX() + bbox.getEndX()) / 2. 

196 if ycen is None: 

197 ycen = (bbox.getBeginY() + bbox.getEndY()) / 2. 

198 self.log.info("Using matching kernel computed at (%d, %d)", xcen, ycen) 

199 spatialKernel.computeImage(kimg, False, xcen, ycen) 

200 

201 if svar is None: 

202 svar = self.computeVarianceMean(scienceExposure) 

203 if tvar is None: 

204 tvar = self.computeVarianceMean(templateExposure) 

205 self.log.infof("Variance (science, template): ({:.5e}, {:.5e})", svar, tvar) 

206 

207 if templateMatched: 

208 # Regular subtraction, we convolved the template 

209 self.log.info("Decorrelation after template image convolution") 

210 expVar = svar 

211 matchedVar = tvar 

212 exposure = scienceExposure 

213 matchedExposure = templateExposure 

214 else: 

215 # We convolved the science image 

216 self.log.info("Decorrelation after science image convolution") 

217 expVar = tvar 

218 matchedVar = svar 

219 exposure = templateExposure 

220 matchedExposure = scienceExposure 

221 

222 # Should not happen unless entire image has been masked, which could happen 

223 # if this is a small subimage of the main exposure. In this case, just return a full NaN 

224 # exposure 

225 if np.isnan(expVar) or np.isnan(matchedVar): 

226 # Double check that one of the exposures is all NaNs 

227 if (np.all(np.isnan(exposure.image.array)) 

228 or np.all(np.isnan(matchedExposure.image.array))): 

229 self.log.warn('Template or science image is entirely NaNs: skipping decorrelation.') 

230 outExposure = subtractedExposure.clone() 

231 return pipeBase.Struct(correctedExposure=outExposure, ) 

232 

233 # The maximal correction value converges to sqrt(matchedVar/expVar). 

234 # Correction divergence warning if the correction exceeds 4 orders of magnitude. 

235 mOverExpVar = matchedVar/expVar 

236 if mOverExpVar > 1e8: 

237 self.log.warn("Diverging correction: matched image variance is " 

238 " much larger than the unconvolved one's" 

239 f", matchedVar/expVar:{mOverExpVar:.2e}") 

240 

241 oldVarMean = self.computeVarianceMean(subtractedExposure) 

242 self.log.info("Variance (uncorrected diffim): %f", oldVarMean) 

243 

244 if preConvKernel is not None: 

245 self.log.info('Using a pre-convolution kernel as part of decorrelation correction.') 

246 kimg2 = afwImage.ImageD(preConvKernel.getDimensions()) 

247 preConvKernel.computeImage(kimg2, False) 

248 pckArr = kimg2.array 

249 

250 kArr = kimg.array 

251 diffExpArr = subtractedExposure.image.array 

252 psfImg = subtractedExposure.getPsf().computeKernelImage(geom.Point2D(xcen, ycen)) 

253 psfDim = psfImg.getDimensions() 

254 psfArr = psfImg.array 

255 

256 # Determine the common shape 

257 kSum = np.sum(kArr) 

258 kSumSq = kSum*kSum 

259 self.log.debugf("Matching kernel sum: {:.2e}", kSum) 

260 preSum = 1. 

261 if preConvKernel is None: 

262 self.computeCommonShape(kArr.shape, psfArr.shape, diffExpArr.shape) 

263 corrft = self.computeCorrection(kArr, expVar, matchedVar) 

264 else: 

265 preSum = np.sum(pckArr) 

266 self.log.debugf("pre-convolution kernel sum: {:.2e}", preSum) 

267 self.computeCommonShape(pckArr.shape, kArr.shape, 

268 psfArr.shape, diffExpArr.shape) 

269 corrft = self.computeCorrection(kArr, expVar, matchedVar, preConvArr=pckArr) 

270 

271 diffExpArr = self.computeCorrectedImage(corrft, diffExpArr) 

272 corrPsfArr = self.computeCorrectedDiffimPsf(corrft, psfArr) 

273 

274 psfcI = afwImage.ImageD(psfDim) 

275 psfcI.array = corrPsfArr 

276 psfcK = afwMath.FixedKernel(psfcI) 

277 psfNew = measAlg.KernelPsf(psfcK) 

278 

279 correctedExposure = subtractedExposure.clone() 

280 correctedExposure.image.array[...] = diffExpArr # Allow for numpy type casting 

281 # The subtracted exposure variance plane is already correlated, we cannot propagate 

282 # it through another convolution; instead we need to use the uncorrelated originals 

283 # The whitening should scale it to expVar + matchedVar on average 

284 varImg = correctedExposure.variance.array 

285 # Allow for numpy type casting 

286 varImg[...] = preSum*preSum*exposure.variance.array + kSumSq*matchedExposure.variance.array 

287 if not templateMatched: 

288 # ImagePsfMatch.subtractExposures re-scales the difference in 

289 # the science image convolution mode 

290 varImg /= kSumSq 

291 correctedExposure.setPsf(psfNew) 

292 

293 newVarMean = self.computeVarianceMean(correctedExposure) 

294 self.log.infof("Variance (corrected diffim): {:.5e}", newVarMean) 

295 

296 # TODO DM-23857 As part of the spatially varying correction implementation 

297 # consider whether returning a Struct is still necessary. 

298 return pipeBase.Struct(correctedExposure=correctedExposure, ) 

299 

300 def computeCommonShape(self, *shapes): 

301 """Calculate the common shape for FFT operations. Set `self.freqSpaceShape` 

302 internally. 

303 

304 Parameters 

305 ---------- 

306 shapes : one or more `tuple` of `int` 

307 Shapes of the arrays. All must have the same dimensionality. 

308 At least one shape must be provided. 

309 

310 Returns 

311 ------- 

312 None. 

313 

314 Notes 

315 ----- 

316 For each dimension, gets the smallest even number greater than or equal to 

317 `N1+N2-1` where `N1` and `N2` are the two largest values. 

318 In case of only one shape given, rounds up to even each dimension value. 

319 """ 

320 S = np.array(shapes, dtype=int) 

321 if len(shapes) > 2: 

322 S.sort(axis=0) 

323 S = S[-2:] 

324 if len(shapes) > 1: 

325 commonShape = np.sum(S, axis=0) - 1 

326 else: 

327 commonShape = S[0] 

328 commonShape[commonShape % 2 != 0] += 1 

329 self.freqSpaceShape = tuple(commonShape) 

330 self.log.info(f"Common frequency space shape {self.freqSpaceShape}") 

331 

332 @staticmethod 

333 def padCenterOriginArray(A, newShape: tuple, useInverse=False): 

334 """Zero pad an image where the origin is at the center and replace the 

335 origin to the corner as required by the periodic input of FFT. Implement also 

336 the inverse operation, crop the padding and re-center data. 

337 

338 Parameters 

339 ---------- 

340 A : `numpy.ndarray` 

341 An array to copy from. 

342 newShape : `tuple` of `int` 

343 The dimensions of the resulting array. For padding, the resulting array 

344 must be larger than A in each dimension. For the inverse operation this 

345 must be the original, before padding size of the array. 

346 useInverse : bool, optional 

347 Selector of forward, add padding, operation (False) 

348 or its inverse, crop padding, operation (True). 

349 

350 Returns 

351 ------- 

352 R : `numpy.ndarray` 

353 The padded or unpadded array with shape of `newShape` and the same dtype as A. 

354 

355 Notes 

356 ----- 

357 For odd dimensions, the splitting is rounded to 

358 put the center pixel into the new corner origin (0,0). This is to be consistent 

359 e.g. for a dirac delta kernel that is originally located at the center pixel. 

360 """ 

361 

362 # The forward and inverse operations should round odd dimension halves at the opposite 

363 # sides to get the pixels back to their original positions. 

364 if not useInverse: 

365 # Forward operation: First and second halves with respect to the axes of A. 

366 firstHalves = [x//2 for x in A.shape] 

367 secondHalves = [x-y for x, y in zip(A.shape, firstHalves)] 

368 else: 

369 # Inverse operation: Opposite rounding 

370 secondHalves = [x//2 for x in newShape] 

371 firstHalves = [x-y for x, y in zip(newShape, secondHalves)] 

372 

373 R = np.zeros_like(A, shape=newShape) 

374 R[-firstHalves[0]:, -firstHalves[1]:] = A[:firstHalves[0], :firstHalves[1]] 

375 R[:secondHalves[0], -firstHalves[1]:] = A[-secondHalves[0]:, :firstHalves[1]] 

376 R[:secondHalves[0], :secondHalves[1]] = A[-secondHalves[0]:, -secondHalves[1]:] 

377 R[-firstHalves[0]:, :secondHalves[1]] = A[:firstHalves[0], -secondHalves[1]:] 

378 return R 

379 

380 def computeCorrection(self, kappa, svar, tvar, preConvArr=None): 

381 """Compute the Lupton decorrelation post-convolution kernel for decorrelating an 

382 image difference, based on the PSF-matching kernel. 

383 

384 Parameters 

385 ---------- 

386 kappa : `numpy.ndarray` 

387 A matching kernel 2-d numpy.array derived from Alard & Lupton PSF matching. 

388 svar : `float` 

389 Average variance of science image used for PSF matching. 

390 tvar : `float` 

391 Average variance of the template (matched) image used for PSF matching. 

392 preConvArr : `numpy.ndarray`, optional 

393 If not None, then pre-filtering was applied 

394 to science exposure, and this is the pre-convolution kernel. 

395 

396 Returns 

397 ------- 

398 corrft : `numpy.ndarray` of `float` 

399 The frequency space representation of the correction. The array is real (dtype float). 

400 Shape is `self.freqSpaceShape`. 

401 

402 Notes 

403 ----- 

404 The maximum correction factor converges to `sqrt(tvar/svar)` towards high frequencies. 

405 This should be a plausible value. 

406 """ 

407 kSum = np.sum(kappa) 

408 kappa = self.padCenterOriginArray(kappa, self.freqSpaceShape) 

409 kft = np.fft.fft2(kappa) 

410 kftAbsSq = np.real(np.conj(kft) * kft) 

411 # If there is no pre-convolution kernel, use placeholder scalars 

412 if preConvArr is None: 

413 preSum = 1. 

414 preAbsSq = 1. 

415 else: 

416 preSum = np.sum(preConvArr) 

417 preConvArr = self.padCenterOriginArray(preConvArr, self.freqSpaceShape) 

418 preK = np.fft.fft2(preConvArr) 

419 preAbsSq = np.real(np.conj(preK)*preK) 

420 

421 denom = svar * preAbsSq + tvar * kftAbsSq 

422 # Division by zero protection, though we don't expect to hit it 

423 # (rather we'll have numerical noise) 

424 tiny = np.finfo(kftAbsSq.dtype).tiny * 1000. 

425 flt = denom < tiny 

426 sumFlt = np.sum(flt) 

427 if sumFlt > 0: 

428 self.log.warnf("Avoid zero division. Skip decorrelation " 

429 "at {} divergent frequencies.", sumFlt) 

430 denom[flt] = 1. 

431 kft = np.sqrt((svar * preSum*preSum + tvar * kSum*kSum) / denom) 

432 # Don't do any correction at these frequencies 

433 # the difference image should be close to zero anyway, so can't be decorrelated 

434 if sumFlt > 0: 

435 kft[flt] = 1. 

436 return kft 

437 

438 def computeCorrectedDiffimPsf(self, corrft, psfOld): 

439 """Compute the (decorrelated) difference image's new PSF. 

440 

441 Parameters 

442 ---------- 

443 corrft : `numpy.ndarray` 

444 The frequency space representation of the correction calculated by 

445 `computeCorrection`. Shape must be `self.freqSpaceShape`. 

446 psfOld : `numpy.ndarray` 

447 The psf of the difference image to be corrected. 

448 

449 Returns 

450 ------- 

451 psfNew : `numpy.ndarray` 

452 The corrected psf, same shape as `psfOld`, sum normed to 1. 

453 

454 Notes 

455 ----- 

456 There is no algorithmic guarantee that the corrected psf can 

457 meaningfully fit to the same size as the original one. 

458 """ 

459 psfShape = psfOld.shape 

460 psfNew = self.padCenterOriginArray(psfOld, self.freqSpaceShape) 

461 psfNew = np.fft.fft2(psfNew) 

462 psfNew *= corrft 

463 psfNew = np.fft.ifft2(psfNew) 

464 psfNew = psfNew.real 

465 psfNew = self.padCenterOriginArray(psfNew, psfShape, useInverse=True) 

466 psfNew = psfNew/psfNew.sum() 

467 return psfNew 

468 

469 def computeCorrectedImage(self, corrft, imgOld): 

470 """Compute the decorrelated difference image. 

471 

472 Parameters 

473 ---------- 

474 corrft : `numpy.ndarray` 

475 The frequency space representation of the correction calculated by 

476 `computeCorrection`. Shape must be `self.freqSpaceShape`. 

477 imgOld : `numpy.ndarray` 

478 The difference image to be corrected. 

479 

480 Returns 

481 ------- 

482 imgNew : `numpy.ndarray` 

483 The corrected image, same size as the input. 

484 """ 

485 expShape = imgOld.shape 

486 imgNew = np.copy(imgOld) 

487 filtInf = np.isinf(imgNew) 

488 filtNan = np.isnan(imgNew) 

489 imgNew[filtInf] = np.nan 

490 imgNew[filtInf | filtNan] = np.nanmean(imgNew) 

491 imgNew = self.padCenterOriginArray(imgNew, self.freqSpaceShape) 

492 imgNew = np.fft.fft2(imgNew) 

493 imgNew *= corrft 

494 imgNew = np.fft.ifft2(imgNew) 

495 imgNew = imgNew.real 

496 imgNew = self.padCenterOriginArray(imgNew, expShape, useInverse=True) 

497 imgNew[filtNan] = np.nan 

498 imgNew[filtInf] = np.inf 

499 return imgNew 

500 

501 

502class DecorrelateALKernelMapper(DecorrelateALKernelTask, ImageMapper): 

503 """Task to be used as an ImageMapper for performing 

504 A&L decorrelation on subimages on a grid across a A&L difference image. 

505 

506 This task subclasses DecorrelateALKernelTask in order to implement 

507 all of that task's configuration parameters, as well as its `run` method. 

508 """ 

509 

510 ConfigClass = DecorrelateALKernelConfig 

511 _DefaultName = 'ip_diffim_decorrelateALKernelMapper' 

512 

513 def __init__(self, *args, **kwargs): 

514 DecorrelateALKernelTask.__init__(self, *args, **kwargs) 

515 

516 def run(self, subExposure, expandedSubExposure, fullBBox, 

517 template, science, alTaskResult=None, psfMatchingKernel=None, 

518 preConvKernel=None, **kwargs): 

519 """Perform decorrelation operation on `subExposure`, using 

520 `expandedSubExposure` to allow for invalid edge pixels arising from 

521 convolutions. 

522 

523 This method performs A&L decorrelation on `subExposure` using 

524 local measures for image variances and PSF. `subExposure` is a 

525 sub-exposure of the non-decorrelated A&L diffim. It also 

526 requires the corresponding sub-exposures of the template 

527 (`template`) and science (`science`) exposures. 

528 

529 Parameters 

530 ---------- 

531 subExposure : `lsst.afw.image.Exposure` 

532 the sub-exposure of the diffim 

533 expandedSubExposure : `lsst.afw.image.Exposure` 

534 the expanded sub-exposure upon which to operate 

535 fullBBox : `lsst.geom.Box2I` 

536 the bounding box of the original exposure 

537 template : `lsst.afw.image.Exposure` 

538 the corresponding sub-exposure of the template exposure 

539 science : `lsst.afw.image.Exposure` 

540 the corresponding sub-exposure of the science exposure 

541 alTaskResult : `lsst.pipe.base.Struct` 

542 the result of A&L image differencing on `science` and 

543 `template`, importantly containing the resulting 

544 `psfMatchingKernel`. Can be `None`, only if 

545 `psfMatchingKernel` is not `None`. 

546 psfMatchingKernel : Alternative parameter for passing the 

547 A&L `psfMatchingKernel` directly. 

548 preConvKernel : If not None, then pre-filtering was applied 

549 to science exposure, and this is the pre-convolution 

550 kernel. 

551 kwargs : 

552 additional keyword arguments propagated from 

553 `ImageMapReduceTask.run`. 

554 

555 Returns 

556 ------- 

557 A `pipeBase.Struct` containing: 

558 

559 - ``subExposure`` : the result of the `subExposure` processing. 

560 - ``decorrelationKernel`` : the decorrelation kernel, currently 

561 not used. 

562 

563 Notes 

564 ----- 

565 This `run` method accepts parameters identical to those of 

566 `ImageMapper.run`, since it is called from the 

567 `ImageMapperTask`. See that class for more information. 

568 """ 

569 templateExposure = template # input template 

570 scienceExposure = science # input science image 

571 if alTaskResult is None and psfMatchingKernel is None: 

572 raise RuntimeError('Both alTaskResult and psfMatchingKernel cannot be None') 

573 psfMatchingKernel = alTaskResult.psfMatchingKernel if alTaskResult is not None else psfMatchingKernel 

574 

575 # subExp and expandedSubExp are subimages of the (un-decorrelated) diffim! 

576 # So here we compute corresponding subimages of templateExposure and scienceExposure 

577 subExp2 = scienceExposure.Factory(scienceExposure, expandedSubExposure.getBBox()) 

578 subExp1 = templateExposure.Factory(templateExposure, expandedSubExposure.getBBox()) 

579 

580 # Prevent too much log INFO verbosity from DecorrelateALKernelTask.run 

581 logLevel = self.log.getLevel() 

582 self.log.setLevel(lsst.log.WARN) 

583 res = DecorrelateALKernelTask.run(self, subExp2, subExp1, expandedSubExposure, 

584 psfMatchingKernel, preConvKernel) 

585 self.log.setLevel(logLevel) # reset the log level 

586 

587 diffim = res.correctedExposure.Factory(res.correctedExposure, subExposure.getBBox()) 

588 out = pipeBase.Struct(subExposure=diffim, ) 

589 return out 

590 

591 

592class DecorrelateALKernelMapReduceConfig(ImageMapReduceConfig): 

593 """Configuration parameters for the ImageMapReduceTask to direct it to use 

594 DecorrelateALKernelMapper as its mapper for A&L decorrelation. 

595 """ 

596 mapper = pexConfig.ConfigurableField( 

597 doc='A&L decorrelation task to run on each sub-image', 

598 target=DecorrelateALKernelMapper 

599 ) 

600 

601 

602class DecorrelateALKernelSpatialConfig(pexConfig.Config): 

603 """Configuration parameters for the DecorrelateALKernelSpatialTask. 

604 """ 

605 decorrelateConfig = pexConfig.ConfigField( 

606 dtype=DecorrelateALKernelConfig, 

607 doc='DecorrelateALKernel config to use when running on complete exposure (non spatially-varying)', 

608 ) 

609 

610 decorrelateMapReduceConfig = pexConfig.ConfigField( 

611 dtype=DecorrelateALKernelMapReduceConfig, 

612 doc='DecorrelateALKernelMapReduce config to use when running on each sub-image (spatially-varying)', 

613 ) 

614 

615 ignoreMaskPlanes = pexConfig.ListField( 

616 dtype=str, 

617 doc="""Mask planes to ignore for sigma-clipped statistics""", 

618 default=("INTRP", "EDGE", "DETECTED", "SAT", "CR", "BAD", "NO_DATA", "DETECTED_NEGATIVE") 

619 ) 

620 

621 def setDefaults(self): 

622 self.decorrelateMapReduceConfig.gridStepX = self.decorrelateMapReduceConfig.gridStepY = 40 

623 self.decorrelateMapReduceConfig.cellSizeX = self.decorrelateMapReduceConfig.cellSizeY = 41 

624 self.decorrelateMapReduceConfig.borderSizeX = self.decorrelateMapReduceConfig.borderSizeY = 8 

625 self.decorrelateMapReduceConfig.reducer.reduceOperation = 'average' 

626 

627 

628class DecorrelateALKernelSpatialTask(pipeBase.Task): 

629 """Decorrelate the effect of convolution by Alard-Lupton matching kernel in image difference 

630 

631 Notes 

632 ----- 

633 

634 Pipe-task that removes the neighboring-pixel covariance in an 

635 image difference that are added when the template image is 

636 convolved with the Alard-Lupton PSF matching kernel. 

637 

638 This task is a simple wrapper around @ref DecorrelateALKernelTask, 

639 which takes a `spatiallyVarying` parameter in its `run` method. If 

640 it is `False`, then it simply calls the `run` method of @ref 

641 DecorrelateALKernelTask. If it is True, then it uses the @ref 

642 ImageMapReduceTask framework to break the exposures into 

643 subExposures on a grid, and performs the `run` method of @ref 

644 DecorrelateALKernelTask on each subExposure. This enables it to 

645 account for spatially-varying PSFs and noise in the exposures when 

646 performing the decorrelation. 

647 

648 This task has no standalone example, however it is applied as a 

649 subtask of pipe.tasks.imageDifference.ImageDifferenceTask. 

650 There is also an example of its use in `tests/testImageDecorrelation.py`. 

651 """ 

652 ConfigClass = DecorrelateALKernelSpatialConfig 

653 _DefaultName = "ip_diffim_decorrelateALKernelSpatial" 

654 

655 def __init__(self, *args, **kwargs): 

656 """Create the image decorrelation Task 

657 

658 Parameters 

659 ---------- 

660 args : 

661 arguments to be passed to 

662 `lsst.pipe.base.task.Task.__init__` 

663 kwargs : 

664 additional keyword arguments to be passed to 

665 `lsst.pipe.base.task.Task.__init__` 

666 """ 

667 pipeBase.Task.__init__(self, *args, **kwargs) 

668 

669 self.statsControl = afwMath.StatisticsControl() 

670 self.statsControl.setNumSigmaClip(3.) 

671 self.statsControl.setNumIter(3) 

672 self.statsControl.setAndMask(afwImage.Mask.getPlaneBitMask(self.config.ignoreMaskPlanes)) 

673 

674 def computeVarianceMean(self, exposure): 

675 """Compute the mean of the variance plane of `exposure`. 

676 """ 

677 statObj = afwMath.makeStatistics(exposure.getMaskedImage().getVariance(), 

678 exposure.getMaskedImage().getMask(), 

679 afwMath.MEANCLIP, self.statsControl) 

680 var = statObj.getValue(afwMath.MEANCLIP) 

681 return var 

682 

683 def run(self, scienceExposure, templateExposure, subtractedExposure, psfMatchingKernel, 

684 spatiallyVarying=True, preConvKernel=None, templateMatched=True): 

685 """Perform decorrelation of an image difference exposure. 

686 

687 Decorrelates the diffim due to the convolution of the 

688 templateExposure with the A&L psfMatchingKernel. If 

689 `spatiallyVarying` is True, it utilizes the spatially varying 

690 matching kernel via the `imageMapReduce` framework to perform 

691 spatially-varying decorrelation on a grid of subExposures. 

692 

693 Parameters 

694 ---------- 

695 scienceExposure : `lsst.afw.image.Exposure` 

696 the science Exposure used for PSF matching 

697 templateExposure : `lsst.afw.image.Exposure` 

698 the template Exposure used for PSF matching 

699 subtractedExposure : `lsst.afw.image.Exposure` 

700 the subtracted Exposure produced by `ip_diffim.ImagePsfMatchTask.subtractExposures()` 

701 psfMatchingKernel : 

702 an (optionally spatially-varying) PSF matching kernel produced 

703 by `ip_diffim.ImagePsfMatchTask.subtractExposures()` 

704 spatiallyVarying : `bool` 

705 if True, perform the spatially-varying operation 

706 preConvKernel : `lsst.meas.algorithms.Psf` 

707 if not none, the scienceExposure has been pre-filtered with this kernel. (Currently 

708 this option is experimental.) 

709 templateMatched : `bool`, optional 

710 If True, the template exposure was matched (convolved) to the science exposure. 

711 

712 Returns 

713 ------- 

714 results : `lsst.pipe.base.Struct` 

715 a structure containing: 

716 

717 - ``correctedExposure`` : the decorrelated diffim 

718 

719 """ 

720 self.log.info('Running A&L decorrelation: spatiallyVarying=%r' % spatiallyVarying) 

721 

722 svar = self.computeVarianceMean(scienceExposure) 

723 tvar = self.computeVarianceMean(templateExposure) 

724 if np.isnan(svar) or np.isnan(tvar): # Should not happen unless entire image has been masked. 

725 # Double check that one of the exposures is all NaNs 

726 if (np.all(np.isnan(scienceExposure.image.array)) 

727 or np.all(np.isnan(templateExposure.image.array))): 

728 self.log.warn('Template or science image is entirely NaNs: skipping decorrelation.') 

729 if np.isnan(svar): 

730 svar = 1e-9 

731 if np.isnan(tvar): 

732 tvar = 1e-9 

733 

734 var = self.computeVarianceMean(subtractedExposure) 

735 

736 if spatiallyVarying: 

737 self.log.info("Variance (science, template): (%f, %f)", svar, tvar) 

738 self.log.info("Variance (uncorrected diffim): %f", var) 

739 config = self.config.decorrelateMapReduceConfig 

740 task = ImageMapReduceTask(config=config) 

741 results = task.run(subtractedExposure, science=scienceExposure, 

742 template=templateExposure, psfMatchingKernel=psfMatchingKernel, 

743 preConvKernel=preConvKernel, forceEvenSized=True, 

744 templateMatched=templateMatched) 

745 results.correctedExposure = results.exposure 

746 

747 # Make sure masks of input image are propagated to diffim 

748 def gm(exp): 

749 return exp.getMaskedImage().getMask() 

750 gm(results.correctedExposure)[:, :] = gm(subtractedExposure) 

751 

752 var = self.computeVarianceMean(results.correctedExposure) 

753 self.log.info("Variance (corrected diffim): %f", var) 

754 

755 else: 

756 config = self.config.decorrelateConfig 

757 task = DecorrelateALKernelTask(config=config) 

758 results = task.run(scienceExposure, templateExposure, 

759 subtractedExposure, psfMatchingKernel, preConvKernel=preConvKernel, 

760 templateMatched=templateMatched) 

761 

762 return results