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.geom as geom 

26import lsst.afw.image as afwImage 

27import lsst.afw.geom as afwGeom 

28import lsst.meas.algorithms as measAlg 

29import lsst.afw.math as afwMath 

30import lsst.pex.config as pexConfig 

31import lsst.pipe.base as pipeBase 

32 

33from .imageMapReduce import (ImageMapReduceConfig, ImageMapper, 

34 ImageMapReduceTask) 

35from .imagePsfMatch import (ImagePsfMatchTask, ImagePsfMatchConfig, 

36 subtractAlgorithmRegistry) 

37 

38__all__ = ["ZogyTask", "ZogyConfig", 

39 "ZogyMapper", "ZogyMapReduceConfig", 

40 "ZogyImagePsfMatchConfig", "ZogyImagePsfMatchTask"] 

41 

42 

43"""Tasks for performing the "Proper image subtraction" algorithm of 

44Zackay, et al. (2016), hereafter simply referred to as 'ZOGY (2016)'. 

45 

46`ZogyTask` contains methods to perform the basic estimation of the 

47ZOGY diffim `D`, its updated PSF, and the variance-normalized 

48likelihood image `S_corr`. We have implemented ZOGY using the 

49proscribed methodology, computing all convolutions in Fourier space, 

50and also variants in which the convolutions are performed in real 

51(image) space. The former is faster and results in fewer artifacts 

52when the PSFs are noisy (i.e., measured, for example, via 

53`PsfEx`). The latter is presumed to be preferred as it can account for 

54masks correctly with fewer "ringing" artifacts from edge effects or 

55saturated stars, but noisy PSFs result in their own smaller 

56artifacts. Removal of these artifacts is a subject of continuing 

57research. Currently, we "pad" the PSFs when performing the 

58subtractions in real space, which reduces, but does not entirely 

59eliminate these artifacts. 

60 

61All methods in `ZogyTask` assume template and science images are 

62already accurately photometrically and astrometrically registered. 

63 

64`ZogyMapper` is a wrapper which runs `ZogyTask` in the 

65`ImageMapReduce` framework, computing of ZOGY diffim's on small, 

66overlapping sub-images, thereby enabling complete ZOGY diffim's which 

67account for spatially-varying noise and PSFs across the two input 

68exposures. An example of the use of this task is in the `testZogy.py` 

69unit test. 

70""" 

71 

72 

73class ZogyConfig(pexConfig.Config): 

74 """Configuration parameters for the ZogyTask 

75 """ 

76 inImageSpace = pexConfig.Field( 

77 dtype=bool, 

78 default=False, 

79 doc="Perform all convolutions in real (image) space rather than Fourier space. " 

80 "Currently if True, this results in artifacts when using real (noisy) PSFs." 

81 ) 

82 

83 padSize = pexConfig.Field( 

84 dtype=int, 

85 default=7, 

86 doc="Number of pixels to pad PSFs to avoid artifacts (when inImageSpace is True)" 

87 ) 

88 

89 templateFluxScaling = pexConfig.Field( 

90 dtype=float, 

91 default=1., 

92 doc="Template flux scaling factor (Fr in ZOGY paper)" 

93 ) 

94 

95 scienceFluxScaling = pexConfig.Field( 

96 dtype=float, 

97 default=1., 

98 doc="Science flux scaling factor (Fn in ZOGY paper)" 

99 ) 

100 

101 scaleByCalibration = pexConfig.Field( 

102 dtype=bool, 

103 default=True, 

104 doc="Compute the flux normalization scaling based on the image calibration." 

105 "This overrides 'templateFluxScaling' and 'scienceFluxScaling'." 

106 ) 

107 

108 doTrimKernels = pexConfig.Field( 

109 dtype=bool, 

110 default=False, 

111 doc="Trim kernels for image-space ZOGY. Speeds up convolutions and shrinks artifacts. " 

112 "Subject of future research." 

113 ) 

114 

115 doFilterPsfs = pexConfig.Field( 

116 dtype=bool, 

117 default=True, 

118 doc="Filter PSFs for image-space ZOGY. Aids in reducing artifacts. " 

119 "Subject of future research." 

120 ) 

121 

122 ignoreMaskPlanes = pexConfig.ListField( 

123 dtype=str, 

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

125 doc="Mask planes to ignore for statistics" 

126 ) 

127 

128 

129MIN_KERNEL = 1.0e-4 

130 

131 

132class ZogyTask(pipeBase.Task): 

133 """Task to perform ZOGY proper image subtraction. See module-level documentation for 

134 additional details. 

135 

136 In all methods, im1 is R (reference, or template) and im2 is N (new, or science). 

137 """ 

138 ConfigClass = ZogyConfig 

139 _DefaultName = "ip_diffim_Zogy" 

140 

141 def __init__(self, templateExposure=None, scienceExposure=None, sig1=None, sig2=None, 

142 psf1=None, psf2=None, *args, **kwargs): 

143 """Create the ZOGY task. 

144 

145 Parameters 

146 ---------- 

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

148 Template exposure ("Reference image" in ZOGY (2016)). 

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

150 Science exposure ("New image" in ZOGY (2016)). Must have already been 

151 registered and photmetrically matched to template. 

152 sig1 : `float` 

153 (Optional) sqrt(variance) of `templateExposure`. If `None`, it is 

154 computed from the sqrt(mean) of the `templateExposure` variance image. 

155 sig2 : `float` 

156 (Optional) sqrt(variance) of `scienceExposure`. If `None`, it is 

157 computed from the sqrt(mean) of the `scienceExposure` variance image. 

158 psf1 : 2D `numpy.array` 

159 (Optional) 2D array containing the PSF image for the template. If 

160 `None`, it is extracted from the PSF taken at the center of `templateExposure`. 

161 psf2 : 2D `numpy.array` 

162 (Optional) 2D array containing the PSF image for the science img. If 

163 `None`, it is extracted from the PSF taken at the center of `scienceExposure`. 

164 *args 

165 additional arguments to be passed to 

166 `lsst.pipe.base.Task` 

167 **kwargs 

168 additional keyword arguments to be passed to 

169 `lsst.pipe.base.Task` 

170 """ 

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

172 self.template = self.science = None 

173 self.setup(templateExposure=templateExposure, scienceExposure=scienceExposure, 

174 sig1=sig1, sig2=sig2, psf1=psf1, psf2=psf2, *args, **kwargs) 

175 

176 def setup(self, templateExposure=None, scienceExposure=None, sig1=None, sig2=None, 

177 psf1=None, psf2=None, correctBackground=False, *args, **kwargs): 

178 """Set up the ZOGY task. 

179 

180 Parameters 

181 ---------- 

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

183 Template exposure ("Reference image" in ZOGY (2016)). 

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

185 Science exposure ("New image" in ZOGY (2016)). Must have already been 

186 registered and photometrically matched to template. 

187 sig1 : `float` 

188 (Optional) sqrt(variance) of `templateExposure`. If `None`, it is 

189 computed from the sqrt(mean) of the `templateExposure` variance image. 

190 sig2 : `float` 

191 (Optional) sqrt(variance) of `scienceExposure`. If `None`, it is 

192 computed from the sqrt(mean) of the `scienceExposure` variance image. 

193 psf1 : 2D `numpy.array` 

194 (Optional) 2D array containing the PSF image for the template. If 

195 `None`, it is extracted from the PSF taken at the center of `templateExposure`. 

196 psf2 : 2D `numpy.array` 

197 (Optional) 2D array containing the PSF image for the science img. If 

198 `None`, it is extracted from the PSF taken at the center of `scienceExposure`. 

199 correctBackground : `bool` 

200 (Optional) subtract sigma-clipped mean of exposures. Zogy doesn't correct 

201 nonzero backgrounds (unlike AL) so subtract them here. 

202 *args 

203 additional arguments to be passed to 

204 `lsst.pipe.base.Task` 

205 **kwargs 

206 additional keyword arguments to be passed to 

207 `lsst.pipe.base.Task` 

208 """ 

209 if self.template is None and templateExposure is None: 

210 return 

211 if self.science is None and scienceExposure is None: 

212 return 

213 

214 self.template = templateExposure 

215 self.science = scienceExposure 

216 

217 self.statsControl = afwMath.StatisticsControl() 

218 self.statsControl.setNumSigmaClip(3.) 

219 self.statsControl.setNumIter(3) 

220 self.statsControl.setAndMask(afwImage.Mask.getPlaneBitMask( 

221 self.config.ignoreMaskPlanes)) 

222 

223 self.im1 = self.template.getMaskedImage().getImage().getArray() 

224 self.im2 = self.science.getMaskedImage().getImage().getArray() 

225 self.im1_var = self.template.getMaskedImage().getVariance().getArray() 

226 self.im2_var = self.science.getMaskedImage().getVariance().getArray() 

227 

228 def selectPsf(psf, exposure): 

229 if psf is not None: 

230 return psf 

231 else: 

232 bbox1 = self.template.getBBox() 

233 xcen = (bbox1.getBeginX() + bbox1.getEndX()) / 2. 

234 ycen = (bbox1.getBeginY() + bbox1.getEndY()) / 2. 

235 return exposure.getPsf().computeKernelImage(geom.Point2D(xcen, ycen)).getArray() 

236 

237 self.im1_psf = selectPsf(psf1, self.template) 

238 self.im2_psf = selectPsf(psf2, self.science) 

239 

240 # Make sure PSFs are the same size. Messy, but should work for all cases. 

241 psf1 = self.im1_psf 

242 psf2 = self.im2_psf 

243 pShape1 = psf1.shape 

244 pShape2 = psf2.shape 

245 if (pShape1[0] < pShape2[0]): 

246 psf1 = np.pad(psf1, ((0, pShape2[0] - pShape1[0]), (0, 0)), mode='constant', constant_values=0.) 

247 elif (pShape2[0] < pShape1[0]): 

248 psf2 = np.pad(psf2, ((0, pShape1[0] - pShape2[0]), (0, 0)), mode='constant', constant_values=0.) 

249 if (pShape1[1] < pShape2[1]): 

250 psf1 = np.pad(psf1, ((0, 0), (0, pShape2[1] - pShape1[1])), mode='constant', constant_values=0.) 

251 elif (pShape2[1] < pShape1[1]): 

252 psf2 = np.pad(psf2, ((0, 0), (0, pShape1[1] - pShape2[1])), mode='constant', constant_values=0.) 

253 

254 # PSFs' centers may be offset relative to each other; now fix that! 

255 maxLoc1 = np.unravel_index(np.argmax(psf1), psf1.shape) 

256 maxLoc2 = np.unravel_index(np.argmax(psf2), psf2.shape) 

257 # *Very* rarely happens but if they're off by >1 pixel, do it more than once. 

258 while (maxLoc1[0] != maxLoc2[0]) or (maxLoc1[1] != maxLoc2[1]): 

259 if maxLoc1[0] > maxLoc2[0]: 

260 psf2[1:, :] = psf2[:-1, :] 

261 elif maxLoc1[0] < maxLoc2[0]: 

262 psf1[1:, :] = psf1[:-1, :] 

263 if maxLoc1[1] > maxLoc2[1]: 

264 psf2[:, 1:] = psf2[:, :-1] 

265 elif maxLoc1[1] < maxLoc2[1]: 

266 psf1[:, 1:] = psf1[:, :-1] 

267 maxLoc1 = np.unravel_index(np.argmax(psf1), psf1.shape) 

268 maxLoc2 = np.unravel_index(np.argmax(psf2), psf2.shape) 

269 

270 # Make sure there are no div-by-zeros 

271 psf1[psf1 < MIN_KERNEL] = MIN_KERNEL 

272 psf2[psf2 < MIN_KERNEL] = MIN_KERNEL 

273 

274 self.im1_psf = psf1 

275 self.im2_psf = psf2 

276 

277 self.sig1 = np.sqrt(self._computeVarianceMean(self.template)) if sig1 is None else sig1 

278 self.sig2 = np.sqrt(self._computeVarianceMean(self.science)) if sig2 is None else sig2 

279 # if sig1 or sig2 are NaN, then the entire region being Zogy-ed is masked. 

280 # Don't worry about it - the result will be masked but avoid warning messages. 

281 if np.isnan(self.sig1) or self.sig1 == 0: 

282 self.sig1 = 1. 

283 if np.isnan(self.sig2) or self.sig2 == 0: 

284 self.sig2 = 1. 

285 

286 # Zogy doesn't correct nonzero backgrounds (unlike AL) so subtract them here. 

287 if correctBackground: 

288 def _subtractImageMean(exposure): 

289 """Compute the sigma-clipped mean of the image of `exposure`.""" 

290 mi = exposure.getMaskedImage() 

291 statObj = afwMath.makeStatistics(mi.getImage(), mi.getMask(), 

292 afwMath.MEANCLIP, self.statsControl) 

293 mean = statObj.getValue(afwMath.MEANCLIP) 

294 if not np.isnan(mean): 

295 mi -= mean 

296 

297 _subtractImageMean(self.template) 

298 _subtractImageMean(self.science) 

299 

300 # Define the normalization of each image from the config 

301 self.Fr = self.config.templateFluxScaling # default is 1 

302 self.Fn = self.config.scienceFluxScaling # default is 1 

303 # If 'scaleByCalibration' is True then these norms are overwritten 

304 if self.config.scaleByCalibration: 

305 calib_template = self.template.getPhotoCalib() 

306 calib_science = self.science.getPhotoCalib() 

307 if calib_template is None: 

308 self.log.warning("No calibration information available for template image.") 

309 if calib_science is None: 

310 self.log.warning("No calibration information available for science image.") 

311 if calib_template is None or calib_science is None: 

312 self.log.warning("Due to lack of calibration information, " 

313 "reverting to templateFluxScaling and scienceFluxScaling.") 

314 else: 

315 self.Fr = 1/calib_template.getCalibrationMean() 

316 self.Fn = 1/calib_science.getCalibrationMean() 

317 

318 self.log.info("Setting template image scaling to Fr=%f" % self.Fr) 

319 self.log.info("Setting science image scaling to Fn=%f" % self.Fn) 

320 

321 self.padSize = self.config.padSize # default is 7 

322 

323 def _computeVarianceMean(self, exposure): 

324 """Compute the sigma-clipped mean of the variance image of `exposure`. 

325 """ 

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

327 exposure.getMaskedImage().getMask(), 

328 afwMath.MEANCLIP, self.statsControl) 

329 var = statObj.getValue(afwMath.MEANCLIP) 

330 return var 

331 

332 @staticmethod 

333 def _padPsfToSize(psf, size): 

334 """Zero-pad `psf` to the dimensions given by `size`. 

335 

336 Parameters 

337 ---------- 

338 psf : 2D `numpy.array` 

339 Input psf to be padded 

340 size : `list` 

341 Two element list containing the dimensions to pad the `psf` to 

342 

343 Returns 

344 ------- 

345 psf : 2D `numpy.array` 

346 The padded copy of the input `psf`. 

347 """ 

348 newArr = np.zeros(size) 

349 # The center of the PSF sould be placed in the center-right. 

350 offset = [size[0]//2 - psf.shape[0]//2, size[1]//2 - psf.shape[1]//2] 

351 tmp = newArr[offset[0]:(psf.shape[0] + offset[0]), offset[1]:(psf.shape[1] + offset[1])] 

352 tmp[:, :] = psf 

353 return newArr 

354 

355 def computePrereqs(self, psf1=None, psf2=None, padSize=0): 

356 """Compute standard ZOGY quantities used by (nearly) all methods. 

357 

358 Many of the ZOGY calculations require similar quantities, including 

359 FFTs of the PSFs, and the "denominator" term (e.g. in eq. 13 of 

360 ZOGY manuscript (2016). This function consolidates many of those 

361 operations. 

362 

363 Parameters 

364 ---------- 

365 psf1 : 2D `numpy.array` 

366 (Optional) Input psf of template, override if already padded 

367 psf2 : 2D `numpy.array` 

368 (Optional) Input psf of science image, override if already padded 

369 padSize : `int`, optional 

370 Number of pixels to pad the image on each side with zeroes. 

371 

372 Returns 

373 ------- 

374 A `lsst.pipe.base.Struct` containing: 

375 - Pr : 2D `numpy.array`, the (possibly zero-padded) template PSF 

376 - Pn : 2D `numpy.array`, the (possibly zero-padded) science PSF 

377 - Pr_hat : 2D `numpy.array`, the FFT of `Pr` 

378 - Pn_hat : 2D `numpy.array`, the FFT of `Pn` 

379 - denom : 2D `numpy.array`, the denominator of equation (13) in ZOGY (2016) manuscript 

380 - Fd : `float`, the relative flux scaling factor between science and template 

381 """ 

382 psf1 = self.im1_psf if psf1 is None else psf1 

383 psf2 = self.im2_psf if psf2 is None else psf2 

384 padSize = self.padSize if padSize is None else padSize 

385 Pr, Pn = psf1, psf2 

386 if padSize > 0: 

387 Pr = ZogyTask._padPsfToSize(psf1, (psf1.shape[0] + padSize, psf1.shape[1] + padSize)) 

388 Pn = ZogyTask._padPsfToSize(psf2, (psf2.shape[0] + padSize, psf2.shape[1] + padSize)) 

389 # Make sure there are no div-by-zeros 

390 psf1[np.abs(psf1) <= MIN_KERNEL] = MIN_KERNEL 

391 psf2[np.abs(psf2) <= MIN_KERNEL] = MIN_KERNEL 

392 

393 sigR, sigN = self.sig1, self.sig2 

394 Pr_hat = np.fft.fft2(Pr) 

395 Pr_hat2 = np.conj(Pr_hat) * Pr_hat 

396 Pn_hat = np.fft.fft2(Pn) 

397 Pn_hat2 = np.conj(Pn_hat) * Pn_hat 

398 denom = np.sqrt((sigN**2 * self.Fr**2 * Pr_hat2) + (sigR**2 * self.Fn**2 * Pn_hat2)) 

399 Fd = self.Fr * self.Fn / np.sqrt(sigN**2 * self.Fr**2 + sigR**2 * self.Fn**2) 

400 

401 res = pipeBase.Struct( 

402 Pr=Pr, Pn=Pn, Pr_hat=Pr_hat, Pn_hat=Pn_hat, denom=denom, Fd=Fd 

403 ) 

404 return res 

405 

406 def computeDiffimFourierSpace(self, debug=False, returnMatchedTemplate=False, **kwargs): 

407 r"""Compute ZOGY diffim `D` as proscribed in ZOGY (2016) manuscript 

408 

409 Parameters 

410 ---------- 

411 debug : `bool`, optional 

412 If set to True, filter the kernels by setting the edges to zero. 

413 returnMatchedTemplate : `bool`, optional 

414 Calculate the template image. 

415 If not set, the returned template will be None. 

416 

417 Notes 

418 ----- 

419 In all functions, im1 is R (reference, or template) and im2 is N (new, or science) 

420 Compute the ZOGY eqn. (13): 

421 

422 .. math:: 

423 

424 \widehat{D} = \frac{Fr\widehat{Pr}\widehat{N} - 

425 F_n\widehat{Pn}\widehat{R}}{\sqrt{\sigma_n^2 Fr^2 

426 \|\widehat{Pr}\|^2 + \sigma_r^2 F_n^2 \|\widehat{Pn}\|^2}} 

427 

428 where :math:`D` is the optimal difference image, :math:`R` and :math:`N` are the 

429 reference and "new" image, respectively, :math:`Pr` and :math:`P_n` are their 

430 PSFs, :math:`Fr` and :math:`Fn` are their flux-based zero-points (which we 

431 will set to one here), :math:`\sigma_r^2` and :math:`\sigma_n^2` are their 

432 variance, and :math:`\widehat{D}` denotes the FT of :math:`D`. 

433 

434 Returns 

435 ------- 

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

437 Result struct with components: 

438 

439 - ``D`` : 2D `numpy.array`, the proper image difference 

440 - ``D_var`` : 2D `numpy.array`, the variance image for `D` 

441 """ 

442 # Do all in fourier space (needs image-sized PSFs) 

443 psf1 = ZogyTask._padPsfToSize(self.im1_psf, self.im1.shape) 

444 psf2 = ZogyTask._padPsfToSize(self.im2_psf, self.im2.shape) 

445 

446 preqs = self.computePrereqs(psf1, psf2, padSize=0) # already padded the PSFs 

447 

448 def _filterKernel(K, trim_amount): 

449 # Filter the wings of Kn, Kr, set to zero 

450 ps = trim_amount 

451 K[:ps, :] = K[-ps:, :] = 0 

452 K[:, :ps] = K[:, -ps:] = 0 

453 return K 

454 

455 Kr_hat = self.Fr * preqs.Pr_hat / preqs.denom 

456 Kn_hat = self.Fn * preqs.Pn_hat / preqs.denom 

457 if debug and self.config.doTrimKernels: # default False 

458 # Suggestion from Barak to trim Kr and Kn to remove artifacts 

459 # Here we just filter them (in image space) to keep them the same size 

460 ps = (Kn_hat.shape[1] - 80)//2 

461 Kn = _filterKernel(np.fft.ifft2(Kn_hat), ps) 

462 Kn_hat = np.fft.fft2(Kn) 

463 Kr = _filterKernel(np.fft.ifft2(Kr_hat), ps) 

464 Kr_hat = np.fft.fft2(Kr) 

465 

466 def processImages(im1, im2, doAdd=False): 

467 # Some masked regions are NaN or infinite!, and FFTs no likey. 

468 im1[np.isinf(im1)] = np.nan 

469 im1[np.isnan(im1)] = np.nanmean(im1) 

470 im2[np.isinf(im2)] = np.nan 

471 im2[np.isnan(im2)] = np.nanmean(im2) 

472 

473 R_hat = np.fft.fft2(im1) 

474 N_hat = np.fft.fft2(im2) 

475 

476 D_hat = Kr_hat * N_hat 

477 D_hat_R = Kn_hat * R_hat 

478 if not doAdd: 

479 D_hat -= D_hat_R 

480 else: 

481 D_hat += D_hat_R 

482 

483 D = np.fft.ifft2(D_hat) 

484 D = np.fft.ifftshift(D.real) / preqs.Fd 

485 

486 R = None 

487 if returnMatchedTemplate: 

488 R = np.fft.ifft2(D_hat_R) 

489 R = np.fft.ifftshift(R.real) / preqs.Fd 

490 

491 return D, R 

492 

493 # First do the image 

494 D, R = processImages(self.im1, self.im2, doAdd=False) 

495 # Do the exact same thing to the var images, except add them 

496 D_var, R_var = processImages(self.im1_var, self.im2_var, doAdd=True) 

497 

498 return pipeBase.Struct(D=D, D_var=D_var, R=R, R_var=R_var) 

499 

500 def _doConvolve(self, exposure, kernel, recenterKernel=False): 

501 """Convolve an Exposure with a decorrelation convolution kernel. 

502 

503 Parameters 

504 ---------- 

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

506 Input exposure to be convolved. 

507 kernel : `numpy.array` 

508 2D `numpy.array` to convolve the image with 

509 recenterKernel : `bool`, optional 

510 Force the kernel center to the pixel with the maximum value. 

511 

512 Returns 

513 ------- 

514 A new `lsst.afw.image.Exposure` with the convolved pixels and the (possibly 

515 re-centered) kernel. 

516 

517 Notes 

518 ----- 

519 - We optionally re-center the kernel if necessary and return the possibly 

520 re-centered kernel 

521 """ 

522 kernelImg = afwImage.ImageD(kernel.shape[1], kernel.shape[0]) 

523 kernelImg.getArray()[:, :] = kernel 

524 kern = afwMath.FixedKernel(kernelImg) 

525 if recenterKernel: 

526 maxloc = np.unravel_index(np.argmax(kernel), kernel.shape) 

527 kern.setCtr(geom.Point2I(maxloc)) 

528 outExp = exposure.clone() # Do this to keep WCS, PSF, masks, etc. 

529 convCntrl = afwMath.ConvolutionControl(doNormalize=False, doCopyEdge=False, 

530 maxInterpolationDistance=0) 

531 try: 

532 afwMath.convolve(outExp.getMaskedImage(), exposure.getMaskedImage(), kern, convCntrl) 

533 except AttributeError: 

534 # Allow exposure to actually be an image/maskedImage 

535 # (getMaskedImage will throw AttributeError in that case) 

536 afwMath.convolve(outExp, exposure, kern, convCntrl) 

537 

538 return outExp, kern 

539 

540 def computeDiffimImageSpace(self, padSize=None, debug=False, **kwargs): 

541 """Compute ZOGY diffim `D` using image-space convlutions 

542 

543 This method is still being debugged as it results in artifacts 

544 when the PSFs are noisy (see module-level docstring). Thus 

545 there are several options still enabled by the `debug` flag, 

546 which are disabled by defult. 

547 

548 Parameters 

549 ---------- 

550 padSize : `int` 

551 The amount to pad the PSFs by 

552 debug : `bool` 

553 Flag to enable debugging tests and options 

554 

555 Returns 

556 ------- 

557 D : `lsst.afw.Exposure` 

558 the proper image difference, including correct variance, 

559 masks, and PSF 

560 """ 

561 preqs = self.computePrereqs(padSize=padSize) 

562 

563 delta = 0. 

564 if debug: 

565 delta = 1. # Regularize the ratio, a possible option to remove artifacts 

566 Kr_hat = (preqs.Pr_hat + delta) / (preqs.denom + delta) 

567 Kn_hat = (preqs.Pn_hat + delta) / (preqs.denom + delta) 

568 Kr = np.fft.ifft2(Kr_hat).real 

569 Kr = np.roll(np.roll(Kr, -1, 0), -1, 1) 

570 Kn = np.fft.ifft2(Kn_hat).real 

571 Kn = np.roll(np.roll(Kn, -1, 0), -1, 1) 

572 

573 def _trimKernel(self, K, trim_amount): 

574 # Trim out the wings of Kn, Kr (see notebook #15) 

575 # only necessary if it's from a measured psf and PsfEx seems to always make PSFs of size 41x41 

576 ps = trim_amount 

577 K = K[ps:-ps, ps:-ps] 

578 return K 

579 

580 padSize = self.padSize if padSize is None else padSize 

581 # Enabling this block (debug=True) makes it slightly faster, but ~25% worse artifacts: 

582 if debug and self.config.doTrimKernels: # default False 

583 # Filtering also makes it slightly faster (zeros are ignored in convolution) 

584 # but a bit worse. Filter the wings of Kn, Kr (see notebook #15) 

585 Kn = _trimKernel(Kn, padSize) 

586 Kr = _trimKernel(Kr, padSize) 

587 

588 # Note these are reverse-labelled, this is CORRECT! 

589 exp1, _ = self._doConvolve(self.template, Kn) 

590 exp2, _ = self._doConvolve(self.science, Kr) 

591 D = exp2 

592 tmp = D.getMaskedImage() 

593 tmp -= exp1.getMaskedImage() 

594 tmp /= preqs.Fd 

595 return pipeBase.Struct(D=D, R=exp1) 

596 

597 def _setNewPsf(self, exposure, psfArr): 

598 """Utility method to set an exposure's PSF when provided as a 2-d numpy.array 

599 """ 

600 bbox = exposure.getBBox() 

601 center = ((bbox.getBeginX() + bbox.getEndX()) // 2., (bbox.getBeginY() + bbox.getEndY()) // 2.) 

602 center = geom.Point2D(center[0], center[1]) 

603 psfI = afwImage.ImageD(psfArr.shape[1], psfArr.shape[0]) 

604 psfI.getArray()[:, :] = psfArr 

605 psfK = afwMath.FixedKernel(psfI) 

606 psfNew = measAlg.KernelPsf(psfK, center) 

607 exposure.setPsf(psfNew) 

608 return exposure 

609 

610 def computeDiffim(self, inImageSpace=None, padSize=None, 

611 returnMatchedTemplate=False, **kwargs): 

612 """Wrapper method to compute ZOGY proper diffim 

613 

614 This method should be used as the public interface for 

615 computing the ZOGY diffim. 

616 

617 Parameters 

618 ---------- 

619 inImageSpace : `bool` 

620 Override config `inImageSpace` parameter 

621 padSize : `int` 

622 Override config `padSize` parameter 

623 returnMatchedTemplate : `bool` 

624 Include the PSF-matched template in the results Struct 

625 **kwargs 

626 additional keyword arguments to be passed to 

627 `computeDiffimFourierSpace` or `computeDiffimImageSpace`. 

628 

629 Returns 

630 ------- 

631 An lsst.pipe.base.Struct containing: 

632 - D : `lsst.afw.Exposure` 

633 the proper image difference, including correct variance, 

634 masks, and PSF 

635 - R : `lsst.afw.Exposure` 

636 If `returnMatchedTemplate` is True, the PSF-matched template 

637 exposure 

638 """ 

639 R = None 

640 inImageSpace = self.config.inImageSpace if inImageSpace is None else inImageSpace 

641 if inImageSpace: 

642 padSize = self.padSize if padSize is None else padSize 

643 res = self.computeDiffimImageSpace(padSize=padSize, **kwargs) 

644 D = res.D 

645 if returnMatchedTemplate: 

646 R = res.R 

647 else: 

648 res = self.computeDiffimFourierSpace(**kwargs) 

649 D = self.science.clone() 

650 D.getMaskedImage().getImage().getArray()[:, :] = res.D 

651 D.getMaskedImage().getVariance().getArray()[:, :] = res.D_var 

652 if returnMatchedTemplate: 

653 R = self.science.clone() 

654 R.getMaskedImage().getImage().getArray()[:, :] = res.R 

655 R.getMaskedImage().getVariance().getArray()[:, :] = res.R_var 

656 

657 psf = self.computeDiffimPsf() 

658 D = self._setNewPsf(D, psf) 

659 return pipeBase.Struct(D=D, R=R) 

660 

661 def computeDiffimPsf(self, padSize=0, keepFourier=False, psf1=None, psf2=None): 

662 """Compute the ZOGY diffim PSF (ZOGY manuscript eq. 14) 

663 

664 Parameters 

665 ---------- 

666 padSize : `int` 

667 Override config `padSize` parameter 

668 keepFourier : `bool` 

669 Return the FFT of the diffim PSF (do not inverse-FFT it) 

670 psf1 : 2D `numpy.array` 

671 (Optional) Input psf of template, override if already padded 

672 psf2 : 2D `numpy.array` 

673 (Optional) Input psf of science image, override if already padded 

674 

675 Returns 

676 ------- 

677 Pd : 2D `numpy.array` 

678 The diffim PSF (or FFT of PSF if `keepFourier=True`) 

679 """ 

680 preqs = self.computePrereqs(psf1=psf1, psf2=psf2, padSize=padSize) 

681 

682 Pd_hat_numerator = (self.Fr * self.Fn * preqs.Pr_hat * preqs.Pn_hat) 

683 Pd_hat = Pd_hat_numerator / (preqs.Fd * preqs.denom) 

684 

685 if keepFourier: 

686 return Pd_hat 

687 

688 Pd = np.fft.ifft2(Pd_hat) 

689 Pd = np.fft.ifftshift(Pd).real 

690 

691 return Pd 

692 

693 def _computeVarAstGradients(self, xVarAst=0., yVarAst=0., inImageSpace=False, 

694 R_hat=None, Kr_hat=None, Kr=None, 

695 N_hat=None, Kn_hat=None, Kn=None): 

696 """Compute the astrometric noise correction terms 

697 

698 Compute the correction for estimated astrometric noise as 

699 proscribed in ZOGY (2016), section 3.3. All convolutions 

700 performed either in real (image) or Fourier space. 

701 

702 Parameters 

703 ---------- 

704 xVarAst, yVarAst : `float` 

705 estimated astrometric noise (variance of astrometric registration errors) 

706 inImageSpace : `bool` 

707 Perform all convolutions in real (image) space rather than Fourier space 

708 R_hat : 2-D `numpy.array` 

709 (Optional) FFT of template image, only required if `inImageSpace=False` 

710 Kr_hat : 2-D `numpy.array` 

711 FFT of Kr kernel (eq. 28 of ZOGY (2016)), only required if `inImageSpace=False` 

712 Kr : 2-D `numpy.array` 

713 Kr kernel (eq. 28 of ZOGY (2016)), only required if `inImageSpace=True`. 

714 Kr is associated with the template (reference). 

715 N_hat : 2-D `numpy.array` 

716 FFT of science image, only required if `inImageSpace=False` 

717 Kn_hat : 2-D `numpy.array` 

718 FFT of Kn kernel (eq. 29 of ZOGY (2016)), only required if `inImageSpace=False` 

719 Kn : 2-D `numpy.array` 

720 Kn kernel (eq. 29 of ZOGY (2016)), only required if `inImageSpace=True`. 

721 Kn is associated with the science (new) image. 

722 

723 Returns 

724 ------- 

725 VastSR, VastSN : 2-D `numpy.array` 

726 Arrays containing the values in eqs. 30 and 32 of ZOGY (2016). 

727 """ 

728 VastSR = VastSN = 0. 

729 if xVarAst + yVarAst > 0: # Do the astrometric variance correction 

730 if inImageSpace: 

731 S_R, _ = self._doConvolve(self.template, Kr) 

732 S_R = S_R.getMaskedImage().getImage().getArray() 

733 else: 

734 S_R = np.fft.ifft2(R_hat * Kr_hat) 

735 gradRx, gradRy = np.gradient(S_R) 

736 VastSR = xVarAst * gradRx**2. + yVarAst * gradRy**2. 

737 

738 if inImageSpace: 

739 S_N, _ = self._doConvolve(self.science, Kn) 

740 S_N = S_N.getMaskedImage().getImage().getArray() 

741 else: 

742 S_N = np.fft.ifft2(N_hat * Kn_hat) 

743 gradNx, gradNy = np.gradient(S_N) 

744 VastSN = xVarAst * gradNx**2. + yVarAst * gradNy**2. 

745 

746 return VastSR, VastSN 

747 

748 def computeScorrFourierSpace(self, xVarAst=0., yVarAst=0., **kwargs): 

749 """Compute corrected likelihood image, optimal for source detection 

750 

751 Compute ZOGY S_corr image. This image can be thresholded for 

752 detection without optimal filtering, and the variance image is 

753 corrected to account for astrometric noise (errors in 

754 astrometric registration whether systematic or due to effects 

755 such as DCR). The calculations here are all performed in 

756 Fourier space, as proscribed in ZOGY (2016). 

757 

758 Parameters 

759 ---------- 

760 xVarAst, yVarAst : `float` 

761 estimated astrometric noise (variance of astrometric registration errors) 

762 

763 Returns 

764 ------- 

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

766 Result struct with components: 

767 

768 - ``S`` : `numpy.array`, the likelihood image S (eq. 12 of ZOGY (2016)) 

769 - ``S_var`` : the corrected variance image (denominator of eq. 25 of ZOGY (2016)) 

770 - ``Dpsf`` : the PSF of the diffim D, likely never to be used. 

771 """ 

772 # Some masked regions are NaN or infinite!, and FFTs no likey. 

773 def fix_nans(im): 

774 """Replace any NaNs or Infs with the mean of the image.""" 

775 isbad = ~np.isfinite(im) 

776 if np.any(isbad): 

777 im[isbad] = np.nan 

778 im[isbad] = np.nanmean(im) 

779 return im 

780 

781 self.im1 = fix_nans(self.im1) 

782 self.im2 = fix_nans(self.im2) 

783 self.im1_var = fix_nans(self.im1_var) 

784 self.im2_var = fix_nans(self.im2_var) 

785 

786 # Do all in fourier space (needs image-sized PSFs) 

787 psf1 = ZogyTask._padPsfToSize(self.im1_psf, self.im1.shape) 

788 psf2 = ZogyTask._padPsfToSize(self.im2_psf, self.im2.shape) 

789 

790 preqs = self.computePrereqs(psf1, psf2, padSize=0) # already padded the PSFs 

791 

792 # Compute D_hat here (don't need D then, for speed) 

793 R_hat = np.fft.fft2(self.im1) 

794 N_hat = np.fft.fft2(self.im2) 

795 D_hat = self.Fr * preqs.Pr_hat * N_hat - self.Fn * preqs.Pn_hat * R_hat 

796 D_hat /= preqs.denom 

797 

798 Pd_hat = self.computeDiffimPsf(padSize=0, keepFourier=True, psf1=psf1, psf2=psf2) 

799 Pd_bar = np.conj(Pd_hat) 

800 S = np.fft.ifft2(D_hat * Pd_bar) 

801 

802 # Adjust the variance planes of the two images to contribute to the final detection 

803 # (eq's 26-29). 

804 Pn_hat2 = np.conj(preqs.Pn_hat) * preqs.Pn_hat 

805 Kr_hat = self.Fr * self.Fn**2. * np.conj(preqs.Pr_hat) * Pn_hat2 / preqs.denom**2. 

806 Pr_hat2 = np.conj(preqs.Pr_hat) * preqs.Pr_hat 

807 Kn_hat = self.Fn * self.Fr**2. * np.conj(preqs.Pn_hat) * Pr_hat2 / preqs.denom**2. 

808 

809 Kr_hat2 = np.fft.fft2(np.fft.ifft2(Kr_hat)**2.) 

810 Kn_hat2 = np.fft.fft2(np.fft.ifft2(Kn_hat)**2.) 

811 var1c_hat = Kr_hat2 * np.fft.fft2(self.im1_var) 

812 var2c_hat = Kn_hat2 * np.fft.fft2(self.im2_var) 

813 

814 # Do the astrometric variance correction 

815 fGradR, fGradN = self._computeVarAstGradients(xVarAst, yVarAst, inImageSpace=False, 

816 R_hat=R_hat, Kr_hat=Kr_hat, 

817 N_hat=N_hat, Kn_hat=Kn_hat) 

818 

819 S_var = np.sqrt(np.fft.ifftshift(np.fft.ifft2(var1c_hat + var2c_hat)) + fGradR + fGradN) 

820 S_var *= preqs.Fd 

821 

822 S = np.fft.ifftshift(np.fft.ifft2(Kn_hat * N_hat - Kr_hat * R_hat)) 

823 S *= preqs.Fd 

824 

825 Pd = self.computeDiffimPsf(padSize=0) 

826 return pipeBase.Struct(S=S.real, S_var=S_var.real, Dpsf=Pd) 

827 

828 def computeScorrImageSpace(self, xVarAst=0., yVarAst=0., padSize=None, **kwargs): 

829 """Compute corrected likelihood image, optimal for source detection 

830 

831 Compute ZOGY S_corr image. This image can be thresholded for 

832 detection without optimal filtering, and the variance image is 

833 corrected to account for astrometric noise (errors in 

834 astrometric registration whether systematic or due to effects 

835 such as DCR). The calculations here are all performed in 

836 Real (image) space. 

837 

838 Parameters 

839 ---------- 

840 xVarAst, yVarAst : `float` 

841 estimated astrometric noise (variance of astrometric registration errors) 

842 

843 Returns 

844 ------- 

845 A `lsst.pipe.base.Struct` containing: 

846 - S : `lsst.afw.image.Exposure`, the likelihood exposure S (eq. 12 of ZOGY (2016)), 

847 including corrected variance, masks, and PSF 

848 - D : `lsst.afw.image.Exposure`, the proper image difference, including correct 

849 variance, masks, and PSF 

850 """ 

851 # Do convolutions in image space 

852 preqs = self.computePrereqs(padSize=0) 

853 

854 padSize = self.padSize if padSize is None else padSize 

855 D = self.computeDiffimImageSpace(padSize=padSize).D 

856 Pd = self.computeDiffimPsf() 

857 D = self._setNewPsf(D, Pd) 

858 Pd_bar = np.fliplr(np.flipud(Pd)) 

859 S, _ = self._doConvolve(D, Pd_bar) 

860 tmp = S.getMaskedImage() 

861 tmp *= preqs.Fd 

862 

863 # Adjust the variance planes of the two images to contribute to the final detection 

864 # (eq's 26-29). 

865 Pn_hat2 = np.conj(preqs.Pn_hat) * preqs.Pn_hat 

866 Kr_hat = self.Fr * self.Fn**2. * np.conj(preqs.Pr_hat) * Pn_hat2 / preqs.denom**2. 

867 Pr_hat2 = np.conj(preqs.Pr_hat) * preqs.Pr_hat 

868 Kn_hat = self.Fn * self.Fr**2. * np.conj(preqs.Pn_hat) * Pr_hat2 / preqs.denom**2. 

869 

870 Kr = np.fft.ifft2(Kr_hat).real 

871 Kr = np.roll(np.roll(Kr, -1, 0), -1, 1) 

872 Kn = np.fft.ifft2(Kn_hat).real 

873 Kn = np.roll(np.roll(Kn, -1, 0), -1, 1) 

874 var1c, _ = self._doConvolve(self.template.getMaskedImage().getVariance(), Kr**2.) 

875 var2c, _ = self._doConvolve(self.science.getMaskedImage().getVariance(), Kn**2.) 

876 

877 # Do the astrometric variance correction 

878 fGradR, fGradN = self._computeVarAstGradients(xVarAst, yVarAst, inImageSpace=True, 

879 Kr=Kr, Kn=Kn) 

880 

881 Smi = S.getMaskedImage() 

882 Smi *= preqs.Fd 

883 S_var = np.sqrt(var1c.getArray() + var2c.getArray() + fGradR + fGradN) 

884 S.getMaskedImage().getVariance().getArray()[:, :] = S_var 

885 S = self._setNewPsf(S, Pd) 

886 

887 # also return diffim since it was calculated and might be desired 

888 return pipeBase.Struct(S=S, D=D) 

889 

890 def computeScorr(self, xVarAst=0., yVarAst=0., inImageSpace=None, padSize=0, **kwargs): 

891 """Wrapper method to compute ZOGY corrected likelihood image, optimal for 

892 source detection 

893 

894 This method should be used as the public interface for 

895 computing the ZOGY S_corr. 

896 

897 Parameters 

898 ---------- 

899 xVarAst, yVarAst : `float` 

900 estimated astrometric noise (variance of astrometric registration errors) 

901 inImageSpace : `bool` 

902 Override config `inImageSpace` parameter 

903 padSize : `int` 

904 Override config `padSize` parameter 

905 

906 Returns 

907 ------- 

908 S : `lsst.afw.image.Exposure` 

909 The likelihood exposure S (eq. 12 of ZOGY (2016)), 

910 including corrected variance, masks, and PSF 

911 """ 

912 inImageSpace = self.config.inImageSpace if inImageSpace is None else inImageSpace 

913 if inImageSpace: 

914 res = self.computeScorrImageSpace(xVarAst=xVarAst, yVarAst=yVarAst, padSize=padSize) 

915 S = res.S 

916 else: 

917 res = self.computeScorrFourierSpace(xVarAst=xVarAst, yVarAst=yVarAst) 

918 

919 S = self.science.clone() 

920 S.getMaskedImage().getImage().getArray()[:, :] = res.S 

921 S.getMaskedImage().getVariance().getArray()[:, :] = res.S_var 

922 S = self._setNewPsf(S, res.Dpsf) 

923 

924 return pipeBase.Struct(S=S) 

925 

926 

927class ZogyMapper(ZogyTask, ImageMapper): 

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

929 ZOGY image subtraction on a grid of subimages. 

930 """ 

931 ConfigClass = ZogyConfig 

932 _DefaultName = 'ip_diffim_ZogyMapper' 

933 

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

935 ImageMapper.__init__(self, *args, **kwargs) 

936 

937 def run(self, subExposure, expandedSubExposure, fullBBox, template, 

938 **kwargs): 

939 """Perform ZOGY proper image subtraction on sub-images 

940 

941 This method performs ZOGY proper image subtraction on 

942 `subExposure` using local measures for image variances and 

943 PSF. `subExposure` is a sub-exposure of the science image. It 

944 also requires the corresponding sub-exposures of the template 

945 (`template`). The operations are actually performed on 

946 `expandedSubExposure` to allow for invalid edge pixels arising 

947 from convolutions, which are then removed. 

948 

949 Parameters 

950 ---------- 

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

952 the sub-exposure of the diffim 

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

954 the expanded sub-exposure upon which to operate 

955 fullBBox : `lsst.geom.Box2I` 

956 the bounding box of the original exposure 

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

958 the template exposure, from which a corresponding sub-exposure 

959 is extracted 

960 **kwargs 

961 additional keyword arguments propagated from 

962 `ImageMapReduceTask.run`. These include: 

963 

964 ``doScorr`` : `bool` 

965 Compute and return the corrected likelihood image S_corr 

966 rather than the proper image difference 

967 ``inImageSpace`` : `bool` 

968 Perform all convolutions in real (image) space rather than 

969 in Fourier space. This option currently leads to artifacts 

970 when using real (measured and noisy) PSFs, thus it is set 

971 to `False` by default. 

972 These kwargs may also include arguments to be propagated to 

973 `ZogyTask.computeDiffim` and `ZogyTask.computeScorr`. 

974 

975 Returns 

976 ------- 

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

978 Result struct with components: 

979 

980 ``subExposure``: Either the subExposure of the proper image difference ``D``, 

981 or (if `doScorr==True`) the corrected likelihood exposure ``S``. 

982 

983 Notes 

984 ----- 

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

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

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

988 """ 

989 bbox = subExposure.getBBox() 

990 center = ((bbox.getBeginX() + bbox.getEndX()) // 2., (bbox.getBeginY() + bbox.getEndY()) // 2.) 

991 center = geom.Point2D(center[0], center[1]) 

992 

993 imageSpace = kwargs.pop('inImageSpace', False) 

994 doScorr = kwargs.pop('doScorr', False) 

995 sigmas = kwargs.pop('sigmas', None) 

996 padSize = kwargs.pop('padSize', 7) 

997 

998 # Psf and image for science img (index 2) 

999 subExp2 = expandedSubExposure 

1000 

1001 # Psf and image for template img (index 1) 

1002 subExp1 = template.Factory(template, expandedSubExposure.getBBox()) 

1003 

1004 if sigmas is None: 

1005 sig1 = sig2 = None 

1006 else: 

1007 # for testing, can use the input sigma (e.g., global value for entire exposure) 

1008 sig1, sig2 = sigmas[0], sigmas[1] 

1009 

1010 def _makePsfSquare(psf): 

1011 # Sometimes CoaddPsf does this. Make it square. 

1012 if psf.shape[0] < psf.shape[1]: 

1013 psf = np.pad(psf, ((0, psf.shape[1] - psf.shape[0]), (0, 0)), mode='constant', 

1014 constant_values=0.) 

1015 elif psf.shape[0] > psf.shape[1]: 

1016 psf = np.pad(psf, ((0, 0), (0, psf.shape[0] - psf.shape[1])), mode='constant', 

1017 constant_values=0.) 

1018 return psf 

1019 

1020 psf2 = subExp2.getPsf().computeKernelImage(center).getArray() 

1021 psf2 = _makePsfSquare(psf2) 

1022 

1023 psf1 = template.getPsf().computeKernelImage(center).getArray() 

1024 psf1 = _makePsfSquare(psf1) 

1025 

1026 # from diffimTests.diffimTests ... 

1027 if subExp1.getDimensions()[0] < psf1.shape[0] or subExp1.getDimensions()[1] < psf1.shape[1]: 

1028 return pipeBase.Struct(subExposure=subExposure) 

1029 

1030 def _filterPsf(psf): 

1031 """Filter a noisy Psf to remove artifacts. Subject of future research.""" 

1032 # only necessary if it's from a measured psf and PsfEx seems to always make PSFs of size 41x41 

1033 if psf.shape[0] == 41: # its from a measured psf 

1034 psf = psf.copy() 

1035 psf[psf < 0] = 0 

1036 psf[0:10, :] = psf[:, 0:10] = psf[31:41, :] = psf[:, 31:41] = 0 

1037 psf /= psf.sum() 

1038 

1039 return psf 

1040 

1041 psf1b = psf2b = None 

1042 if self.config.doFilterPsfs: # default True 

1043 # Note this *really* helps for measured psfs. 

1044 psf1b = _filterPsf(psf1) 

1045 psf2b = _filterPsf(psf2) 

1046 

1047 config = ZogyConfig() 

1048 if imageSpace is True: 

1049 config.inImageSpace = imageSpace 

1050 config.padSize = padSize # Don't need padding if doing all in fourier space 

1051 task = ZogyTask(templateExposure=subExp1, scienceExposure=subExp2, 

1052 sig1=sig1, sig2=sig2, psf1=psf1b, psf2=psf2b, config=config) 

1053 

1054 if not doScorr: 

1055 res = task.computeDiffim(**kwargs) 

1056 D = res.D 

1057 else: 

1058 res = task.computeScorr(**kwargs) 

1059 D = res.S 

1060 

1061 outExp = D.Factory(D, subExposure.getBBox()) 

1062 out = pipeBase.Struct(subExposure=outExp) 

1063 return out 

1064 

1065 

1066class ZogyMapReduceConfig(ImageMapReduceConfig): 

1067 """Config to be passed to ImageMapReduceTask 

1068 

1069 This config targets the imageMapper to use the ZogyMapper. 

1070 """ 

1071 mapper = pexConfig.ConfigurableField( 

1072 doc='Zogy task to run on each sub-image', 

1073 target=ZogyMapper 

1074 ) 

1075 

1076 

1077class ZogyImagePsfMatchConfig(ImagePsfMatchConfig): 

1078 """Config for the ZogyImagePsfMatchTask""" 

1079 

1080 zogyConfig = pexConfig.ConfigField( 

1081 dtype=ZogyConfig, 

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

1083 ) 

1084 

1085 zogyMapReduceConfig = pexConfig.ConfigField( 

1086 dtype=ZogyMapReduceConfig, 

1087 doc='ZogyMapReduce config to use when running Zogy on each sub-image (spatially-varying)', 

1088 ) 

1089 

1090 def setDefaults(self): 

1091 self.zogyMapReduceConfig.gridStepX = self.zogyMapReduceConfig.gridStepY = 40 

1092 self.zogyMapReduceConfig.cellSizeX = self.zogyMapReduceConfig.cellSizeY = 41 

1093 self.zogyMapReduceConfig.borderSizeX = self.zogyMapReduceConfig.borderSizeY = 8 

1094 self.zogyMapReduceConfig.reducer.reduceOperation = 'average' 

1095 self.zogyConfig.inImageSpace = False 

1096 

1097 

1098class ZogyImagePsfMatchTask(ImagePsfMatchTask): 

1099 """Task to perform Zogy PSF matching and image subtraction. 

1100 

1101 This class inherits from ImagePsfMatchTask to contain the _warper 

1102 subtask and related methods. 

1103 """ 

1104 

1105 ConfigClass = ZogyImagePsfMatchConfig 

1106 

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

1108 ImagePsfMatchTask.__init__(self, *args, **kwargs) 

1109 

1110 def _computeImageMean(self, exposure): 

1111 """Compute the sigma-clipped mean of the pixels image of `exposure`. 

1112 """ 

1113 statsControl = afwMath.StatisticsControl() 

1114 statsControl.setNumSigmaClip(3.) 

1115 statsControl.setNumIter(3) 

1116 ignoreMaskPlanes = ("INTRP", "EDGE", "DETECTED", "SAT", "CR", "BAD", "NO_DATA", "DETECTED_NEGATIVE") 

1117 statsControl.setAndMask(afwImage.Mask.getPlaneBitMask(ignoreMaskPlanes)) 

1118 statObj = afwMath.makeStatistics(exposure.getMaskedImage().getImage(), 

1119 exposure.getMaskedImage().getMask(), 

1120 afwMath.MEANCLIP | afwMath.MEDIAN, statsControl) 

1121 mn = statObj.getValue(afwMath.MEANCLIP) 

1122 med = statObj.getValue(afwMath.MEDIAN) 

1123 return mn, med 

1124 

1125 def subtractExposures(self, templateExposure, scienceExposure, 

1126 doWarping=True, spatiallyVarying=True, inImageSpace=False, 

1127 doPreConvolve=False): 

1128 """Register, PSF-match, and subtract two Exposures using the ZOGY algorithm. 

1129 

1130 Do the following, in order: 

1131 - Warp templateExposure to match scienceExposure, if their WCSs do not already match 

1132 - Compute subtracted exposure ZOGY image subtraction algorithm on the two exposures 

1133 

1134 Parameters 

1135 ---------- 

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

1137 exposure to PSF-match to scienceExposure. The exposure's mean value is subtracted 

1138 in-place. 

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

1140 reference Exposure. The exposure's mean value is subtracted in-place. 

1141 doWarping : `bool` 

1142 what to do if templateExposure's and scienceExposure's WCSs do not match: 

1143 - if True then warp templateExposure to match scienceExposure 

1144 - if False then raise an Exception 

1145 spatiallyVarying : `bool` 

1146 If True, perform the operation over a grid of patches across the two exposures 

1147 inImageSpace : `bool` 

1148 If True, perform the Zogy convolutions in image space rather than in frequency space. 

1149 doPreConvolve : `bool` 

1150 ***Currently not implemented.*** If True assume we are to compute the match filter-convolved 

1151 exposure which can be thresholded for detection. In the case of Zogy this would mean 

1152 we compute the Scorr image. 

1153 

1154 Returns 

1155 ------- 

1156 A `lsst.pipe.base.Struct` containing these fields: 

1157 - subtractedExposure: subtracted Exposure 

1158 - warpedExposure: templateExposure after warping to match scienceExposure (if doWarping true) 

1159 """ 

1160 

1161 mn1 = self._computeImageMean(templateExposure) 

1162 mn2 = self._computeImageMean(scienceExposure) 

1163 self.log.info("Exposure means=%f, %f; median=%f, %f:" % (mn1[0], mn2[0], mn1[1], mn2[1])) 

1164 if not np.isnan(mn1[0]) and np.abs(mn1[0]) > 1: 

1165 mi = templateExposure.getMaskedImage() 

1166 mi -= mn1[0] 

1167 if not np.isnan(mn2[0]) and np.abs(mn2[0]) > 1: 

1168 mi = scienceExposure.getMaskedImage() 

1169 mi -= mn2[0] 

1170 

1171 self.log.info('Running Zogy algorithm: spatiallyVarying=%r' % spatiallyVarying) 

1172 

1173 if not self._validateWcs(templateExposure, scienceExposure): 

1174 if doWarping: 

1175 self.log.info("Astrometrically registering template to science image") 

1176 # Also warp the PSF 

1177 xyTransform = afwGeom.makeWcsPairTransform(templateExposure.getWcs(), 

1178 scienceExposure.getWcs()) 

1179 psfWarped = measAlg.WarpedPsf(templateExposure.getPsf(), xyTransform) 

1180 templateExposure = self._warper.warpExposure(scienceExposure.getWcs(), 

1181 templateExposure, 

1182 destBBox=scienceExposure.getBBox()) 

1183 

1184 templateExposure.setPsf(psfWarped) 

1185 else: 

1186 self.log.error("ERROR: Input images not registered") 

1187 raise RuntimeError("Input images not registered") 

1188 

1189 def gm(exp): 

1190 return exp.getMaskedImage().getMask() 

1191 

1192 def ga(exp): 

1193 return exp.getMaskedImage().getImage().getArray() 

1194 

1195 if self.config.zogyConfig.inImageSpace: 

1196 inImageSpace = True # Override 

1197 self.log.info('Running Zogy algorithm: inImageSpace=%r' % inImageSpace) 

1198 if spatiallyVarying: 

1199 config = self.config.zogyMapReduceConfig 

1200 task = ImageMapReduceTask(config=config) 

1201 results = task.run(scienceExposure, template=templateExposure, inImageSpace=inImageSpace, 

1202 doScorr=doPreConvolve, forceEvenSized=False) 

1203 results.D = results.exposure 

1204 # The CoaddPsf, when used for detection does not utilize its spatially-varying 

1205 # properties; it simply computes the PSF at its getAveragePosition(). 

1206 # TODO: we need to get it to return the matchedExposure (convolved template) 

1207 # too, for dipole fitting; but the imageMapReduce task might need to be engineered 

1208 # for this purpose. 

1209 else: 

1210 config = self.config.zogyConfig 

1211 task = ZogyTask(scienceExposure=scienceExposure, templateExposure=templateExposure, 

1212 config=config) 

1213 if not doPreConvolve: 

1214 results = task.computeDiffim(inImageSpace=inImageSpace) 

1215 results.matchedExposure = results.R 

1216 else: 

1217 results = task.computeScorr(inImageSpace=inImageSpace) 

1218 results.D = results.S 

1219 

1220 # Make sure masks of input images are propagated to diffim 

1221 mask = results.D.getMaskedImage().getMask() 

1222 mask |= scienceExposure.getMaskedImage().getMask() 

1223 mask |= templateExposure.getMaskedImage().getMask() 

1224 results.D.getMaskedImage().getMask()[:, :] = mask 

1225 badBitsNan = mask.addMaskPlane('UNMASKEDNAN') 

1226 resultsArr = results.D.getMaskedImage().getMask().getArray() 

1227 resultsArr[np.isnan(resultsArr)] |= badBitsNan 

1228 resultsArr[np.isnan(scienceExposure.getMaskedImage().getImage().getArray())] |= badBitsNan 

1229 resultsArr[np.isnan(templateExposure.getMaskedImage().getImage().getArray())] |= badBitsNan 

1230 

1231 results.subtractedExposure = results.D 

1232 results.warpedExposure = templateExposure 

1233 return results 

1234 

1235 def subtractMaskedImages(self, templateExposure, scienceExposure, 

1236 doWarping=True, spatiallyVarying=True, inImageSpace=False, 

1237 doPreConvolve=False): 

1238 raise NotImplementedError 

1239 

1240 

1241subtractAlgorithmRegistry.register('zogy', ZogyImagePsfMatchTask)