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

26import lsst.afw.image as afwImage 

27import lsst.afw.math as afwMath 

28import lsst.meas.algorithms as measAlg 

29import lsst.pipe.base as pipeBase 

30import lsst.pex.config as pexConfig 

31 

32from .imagePsfMatch import (ImagePsfMatchTask, ImagePsfMatchConfig, 

33 subtractAlgorithmRegistry) 

34 

35__all__ = ["ZogyTask", "ZogyConfig", 

36 "ZogyImagePsfMatchConfig", "ZogyImagePsfMatchTask"] 

37 

38 

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

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

41 

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

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

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

45proscribed methodology, computing all convolutions in Fourier space, 

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

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

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

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

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

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

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

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

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

55eliminate these artifacts. 

56 

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

58already accurately photometrically and astrometrically registered. 

59 

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

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

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

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

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

65unit test. 

66""" 

67 

68 

69class ZogyConfig(pexConfig.Config): 

70 """Configuration parameters for the ZogyTask 

71 """ 

72 

73 templateFluxScaling = pexConfig.Field( 

74 dtype=float, 

75 default=1., 

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

77 ) 

78 

79 scienceFluxScaling = pexConfig.Field( 

80 dtype=float, 

81 default=1., 

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

83 ) 

84 

85 scaleByCalibration = pexConfig.Field( 

86 dtype=bool, 

87 default=True, 

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

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

90 ) 

91 

92 correctBackground = pexConfig.Field( 

93 dtype=bool, 

94 default=False, 

95 doc="Subtract exposure background mean to have zero expectation value." 

96 ) 

97 

98 ignoreMaskPlanes = pexConfig.ListField( 

99 dtype=str, 

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

101 doc="Mask planes to ignore for statistics" 

102 ) 

103 maxPsfCentroidDist = pexConfig.Field( 

104 dtype=float, 

105 default=0.2, 

106 doc="Maximum centroid difference allowed between the two exposure PSFs (pixels)." 

107 ) 

108 

109 

110class ZogyTask(pipeBase.Task): 

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

112 additional details. 

113 

114 """ 

115 ConfigClass = ZogyConfig 

116 _DefaultName = "imageDifferenceZogy" 

117 

118 def _computeVarianceMean(self, exposure): 

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

120 """ 

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

122 exposure.getMaskedImage().getMask(), 

123 afwMath.MEANCLIP, self.statsControl) 

124 var = statObj.getValue(afwMath.MEANCLIP) 

125 return var 

126 

127 @staticmethod 

128 def padCenterOriginArray(A, newShape, useInverse=False, dtype=None): 

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

130 origin to the corner as required by the periodic input of FFT. 

131 

132 Implement also the inverse operation, crop the padding and re-center data. 

133 

134 Parameters 

135 ---------- 

136 A : `numpy.ndarray` 

137 An array to copy from. 

138 newShape : `tuple` of `int` 

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

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

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

142 useInverse : bool, optional 

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

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

145 dtype: `numpy.dtype`, optional 

146 Dtype of output array. Values must be implicitly castable to this type. 

147 Use to get expected result type, e.g. single float (nympy.float32). 

148 If not specified, dtype is inherited from ``A``. 

149 

150 Returns 

151 ------- 

152 R : `numpy.ndarray` 

153 The padded or unpadded array with shape of `newShape` and dtype of ``dtype``. 

154 

155 Notes 

156 ----- 

157 For odd dimensions, the splitting is rounded to 

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

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

160 

161 

162 Raises 

163 ------ 

164 ValueError : ``newShape`` dimensions must be greater than or equal to the 

165 dimensions of ``A`` for the forward operation and less than or equal to 

166 for the inverse operation. 

167 """ 

168 

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

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

171 if not useInverse: 

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

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

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

175 for d1, d2 in zip(newShape, A.shape): 

176 if d1 < d2: 

177 raise ValueError("Newshape dimensions must be greater or equal") 

178 else: 

179 # Inverse operation: Opposite rounding 

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

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

182 for d1, d2 in zip(newShape, A.shape): 

183 if d1 > d2: 

184 raise ValueError("Newshape dimensions must be smaller or equal") 

185 

186 if dtype is None: 

187 dtype = A.dtype 

188 

189 R = np.zeros(newShape, dtype=dtype) 

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

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

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

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

194 return R 

195 

196 def computeCommonShape(self, *shapes): 

197 """Calculate the common shape for FFT operations. 

198 

199 Set ``self.freqSpaceShape`` internally. 

200 

201 Parameters 

202 ---------- 

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

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

205 At least one shape must be provided. 

206 

207 Returns 

208 ------- 

209 None 

210 

211 Notes 

212 ----- 

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

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

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

216 """ 

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

218 if len(shapes) > 2: 

219 S.sort(axis=0) 

220 S = S[-2:] 

221 if len(shapes) > 1: 

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

223 else: 

224 commonShape = S[0] 

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

226 self.freqSpaceShape = tuple(commonShape) 

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

228 

229 def padAndFftImage(self, imgArr): 

230 """Prepare and forward FFT an image array. 

231 

232 Parameters 

233 ---------- 

234 imgArr : `numpy.ndarray` of `float` 

235 Original array. In-place modified as `numpy.nan` and `numpy.inf` are replaced by 

236 array mean. 

237 

238 Returns 

239 ------- 

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

241 - ``imFft`` : `numpy.ndarray` of `numpy.complex`. 

242 FFT of image. 

243 - ``filtInf``, ``filtNaN`` : `numpy.ndarray` of `bool` 

244 

245 Notes 

246 ----- 

247 Save location of non-finite values for restoration, and replace them 

248 with image mean values. Re-center and zero pad array by `padCenterOriginArray`. 

249 """ 

250 filtInf = np.isinf(imgArr) 

251 filtNaN = np.isnan(imgArr) 

252 imgArr[filtInf] = np.nan 

253 imgArr[filtInf | filtNaN] = np.nanmean(imgArr) 

254 self.log.debug("Replacing {} Inf and {} NaN values.".format( 

255 np.sum(filtInf), np.sum(filtNaN))) 

256 imgArr = self.padCenterOriginArray(imgArr, self.freqSpaceShape) 

257 imgArr = np.fft.fft2(imgArr) 

258 return pipeBase.Struct(imFft=imgArr, filtInf=filtInf, filtNaN=filtNaN) 

259 

260 def inverseFftAndCropImage(self, imgArr, origSize, filtInf=None, filtNaN=None, dtype=None): 

261 """Inverse FFT and crop padding from image array. 

262 

263 Parameters 

264 ---------- 

265 imgArr : `numpy.ndarray` of `numpy.complex` 

266 Fourier space array representing a real image. 

267 

268 origSize : `tuple` of `int` 

269 Original unpadded shape tuple of the image to be cropped to. 

270 

271 filtInf, filtNan : `numpy.ndarray` of indices, optional 

272 If specified, they are used as index arrays for ``result`` to set values to 

273 `numpy.inf` and `numpy.nan` respectively at these positions. 

274 

275 dtype : `numpy.dtype`, optional 

276 Dtype of result array to cast return values to implicitly. This is to 

277 spare one array copy operation at reducing double precision to single. 

278 If `None` result inherits dtype of `imgArr`. 

279 

280 Returns 

281 ------- 

282 result : `numpy.ndarray` of `dtype` 

283 """ 

284 imgNew = np.fft.ifft2(imgArr) 

285 imgNew = imgNew.real 

286 imgNew = self.padCenterOriginArray(imgNew, origSize, useInverse=True, dtype=dtype) 

287 if filtInf is not None: 

288 imgNew[filtInf] = np.inf 

289 if filtNaN is not None: 

290 imgNew[filtNaN] = np.nan 

291 return imgNew 

292 

293 @staticmethod 

294 def computePsfAtCenter(exposure): 

295 """Computes the PSF image at the bbox center point. 

296 

297 This may be at a fractional pixel position. 

298 

299 Parameters 

300 ---------- 

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

302 Exposure with psf. 

303 

304 Returns 

305 ------- 

306 psfImg : `lsst.afw.image.Image` 

307 Calculated psf image. 

308 """ 

309 bbox = exposure.getBBox() 

310 cen = bbox.getCenter() 

311 psf = exposure.getPsf() 

312 psfImg = psf.computeKernelImage(cen) # Centered and normed 

313 return psfImg 

314 

315 @staticmethod 

316 def subtractImageMean(image, mask, statsControl): 

317 """In-place subtraction of sigma-clipped mean of the image. 

318 

319 Parameters 

320 ---------- 

321 image : `lsst.afw.image.Image` 

322 Image to manipulate. Its sigma clipped mean is in-place subtracted. 

323 

324 mask : `lsst.afw.image.Mask` 

325 Mask to use for ignoring pixels. 

326 

327 statsControl : `lsst.afw.math.StatisticsControl` 

328 Config of sigma clipped mean statistics calculation. 

329 

330 Returns 

331 ------- 

332 None 

333 

334 Raises 

335 ------ 

336 ValueError : If image mean is nan. 

337 """ 

338 statObj = afwMath.makeStatistics(image, mask, 

339 afwMath.MEANCLIP, statsControl) 

340 mean = statObj.getValue(afwMath.MEANCLIP) 

341 if not np.isnan(mean): 

342 image -= mean 

343 else: 

344 raise ValueError("Image mean is NaN.") 

345 

346 def prepareFullExposure(self, exposure1, exposure2, correctBackground=False): 

347 """Performs calculations that apply to the full exposures once only in the psf matching. 

348 

349 Parameters 

350 ---------- 

351 

352 correctBackground : `bool`, optional 

353 If True, subtracts sigma-clipped mean of exposures. The algorithm 

354 assumes zero expectation value at background pixels. 

355 

356 Returns 

357 ------- 

358 None 

359 

360 Notes 

361 ----- 

362 Set a number of instance fields with pre-calculated values. ``psfShape``, 

363 ``imgShape`` fields follow the numpy ndarray shape convention i.e. height, 

364 width. 

365 

366 Raises 

367 ------ 

368 ValueError : If photometric calibrations are not available while 

369 ``config.scaleByCalibration`` equals True. 

370 """ 

371 

372 self.statsControl = afwMath.StatisticsControl() 

373 self.statsControl.setNumSigmaClip(3.) 

374 self.statsControl.setNumIter(3) 

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

376 self.config.ignoreMaskPlanes)) 

377 

378 exposure1 = exposure1.clone() 

379 exposure2 = exposure2.clone() 

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

381 if self.config.scaleByCalibration: 

382 calibObj1 = exposure1.getPhotoCalib() 

383 calibObj2 = exposure2.getPhotoCalib() 

384 if calibObj1 is None or calibObj2 is None: 

385 raise ValueError("Photometric calibrations are not available for both exposures.") 

386 mImg1 = calibObj1.calibrateImage(exposure1.maskedImage) 

387 mImg2 = calibObj2.calibrateImage(exposure2.maskedImage) 

388 self.F1 = 1. 

389 self.F2 = 1. 

390 else: 

391 self.F1 = self.config.templateFluxScaling # default is 1 

392 self.F2 = self.config.scienceFluxScaling # default is 1 

393 mImg1 = exposure1.maskedImage 

394 mImg2 = exposure2.maskedImage 

395 

396 # mImgs can be in-place modified 

397 if correctBackground: 

398 self.subtractImageMean(mImg1.image, mImg1.mask, self.statsControl) 

399 self.subtractImageMean(mImg2.image, mImg2.mask, self.statsControl) 

400 

401 psfBBox1 = exposure1.getPsf().computeBBox() 

402 psfBBox2 = exposure2.getPsf().computeBBox() 

403 # Shapes for numpy arrays 

404 self.psfShape1 = (psfBBox1.getHeight(), psfBBox1.getWidth()) 

405 self.psfShape2 = (psfBBox2.getHeight(), psfBBox2.getWidth()) 

406 self.imgShape = (mImg1.getHeight(), mImg1.getWidth()) 

407 # We need the calibrated, full size original 

408 # MaskedImages for the variance plane calculations 

409 exposure1.maskedImage = mImg1 

410 exposure2.maskedImage = mImg2 

411 # TODO DM-25174 : Here we need actually not psfShape but an 

412 # estimation of the size of Pd and Ps 

413 # worst case scenario a padding of imgShape (? TBC) 

414 self.computeCommonShape(self.imgShape, self.psfShape1, self.psfShape2) 

415 

416 self.fullExp1 = exposure1 

417 self.fullExp2 = exposure2 

418 

419 self.fftFullIm1 = self.padAndFftImage(mImg1.image.array) 

420 self.fftVarPl1 = self.padAndFftImage(mImg1.variance.array) 

421 self.fftFullIm2 = self.padAndFftImage(mImg2.image.array) 

422 self.fftVarPl2 = self.padAndFftImage(mImg2.variance.array) 

423 

424 def prepareSubExposure(self, bbox1=None, bbox2=None, psf1=None, psf2=None, sig1=None, sig2=None): 

425 """Perform per-sub exposure preparations. 

426 

427 Parameters 

428 ---------- 

429 sig1, sig2 : `float`, optional 

430 For debug purposes only, copnsider that the image 

431 may already be rescaled by the photometric calibration. 

432 bbox1, bbox2 : `lsst.geom.Box2I`, optional 

433 If specified, the region of the full exposure to use. 

434 

435 psf1, psf2 : `lsst.afw.detection.Psf`, optional 

436 If specified, use given psf as the sub exposure psf. For debug purposes. 

437 

438 sig1, sig2 : `float`, optional 

439 If specified, use value as the sub-exposures' background noise sigma value. 

440 

441 Returns 

442 ------- 

443 None 

444 

445 Notes 

446 ----- 

447 TODO DM-23855: Performing ZOGY on a grid is not yet implemented. 

448 Set (replace) a number of instance fields with pre-calculated values 

449 about the current sub exposure including the FFT of the psfs. 

450 

451 Raises 

452 ------ 

453 ValueError: If sub-exposure dimensions do not match. 

454 """ 

455 if bbox1 is None: 

456 subExposure1 = self.fullExp1.clone() 

457 else: 

458 subExposure1 = self.fullExp1.Factory(self.exposure1, bbox1) 

459 if bbox2 is None: 

460 subExposure2 = self.fullExp2.clone() 

461 else: 

462 subExposure2 = self.fullExp2.Factory(self.exposure2, bbox2) 

463 

464 if subExposure1.getDimensions() != subExposure2.getDimensions(): 

465 raise ValueError("Subexposure dimensions do not match.") 

466 

467 if psf1 is None: 

468 self.subExpPsf1 = self.computePsfAtCenter(subExposure1) 

469 else: 

470 self.subExpPsf1 = psf1 

471 if psf2 is None: 

472 self.subExpPsf2 = self.computePsfAtCenter(subExposure2) 

473 else: 

474 self.subExpPsf2 = psf2 

475 self.checkCentroids(self.subExpPsf1.array, self.subExpPsf2.array) 

476 # sig1 and sig2 should not be set externally, just for debug purpose 

477 if sig1 is None: 

478 sig1 = np.sqrt(self._computeVarianceMean(subExposure1)) 

479 self.subExpVar1 = sig1*sig1 

480 if sig2 is None: 

481 sig2 = np.sqrt(self._computeVarianceMean(subExposure2)) 

482 self.subExpVar2 = sig2*sig2 

483 

484 D = self.padCenterOriginArray(self.subExpPsf1.array, self.freqSpaceShape) 

485 self.psfFft1 = np.fft.fft2(D) 

486 D = self.padCenterOriginArray(self.subExpPsf2.array, self.freqSpaceShape) 

487 self.psfFft2 = np.fft.fft2(D) 

488 

489 self.subExposure1 = subExposure1 

490 self.subExposure2 = subExposure2 

491 

492 @staticmethod 

493 def pixelSpaceSquare(D): 

494 """Square the argument in pixel space. 

495 

496 Parameters 

497 ---------- 

498 D : 2D `numpy.ndarray` of `numpy.complex` 

499 Fourier transform of a real valued array. 

500 

501 Returns 

502 ------- 

503 R : `numpy.ndarray` of `numpy.complex` 

504 

505 Notes 

506 ----- 

507 ``D`` is to be inverse Fourier transformed, squared and then 

508 forward Fourier transformed again, i.e. an autoconvolution in Fourier space. 

509 This operation is not distributive over multiplication. 

510 ``pixelSpaceSquare(A*B) != pixelSpaceSquare(A)*pixelSpaceSquare(B)`` 

511 """ 

512 R = np.real(np.fft.ifft2(D)) 

513 R *= R 

514 R = np.fft.fft2(R) 

515 return R 

516 

517 @staticmethod 

518 def getCentroid(A): 

519 """Calculate the centroid coordinates of a 2D array. 

520 

521 Parameters 

522 ---------- 

523 A : 2D `numpy.ndarray` of `float` 

524 The input array. Must not be all exact zero. 

525 

526 Notes 

527 ----- 

528 Calculates the centroid as if the array represented a 2D geometrical shape with 

529 weights per cell, allowing for "negative" weights. If sum equals to exact (float) zero, 

530 calculates centroid of absolute value array. 

531 

532 The geometrical center is defined as (0,0), independently of the array shape. 

533 For an odd dimension, this is the center of the center pixel, 

534 for an even dimension, this is between the two center pixels. 

535 

536 Returns 

537 ------- 

538 ycen, xcen : `tuple` of `float` 

539 

540 """ 

541 s = np.sum(A) 

542 if s == 0.: 

543 A = np.fabs(A) 

544 s = np.sum(A) 

545 w = np.arange(A.shape[0], dtype=float) - (A.shape[0] - 1.)/2 

546 ycen = np.sum(w[:, np.newaxis]*A)/s 

547 w = np.arange(A.shape[1], dtype=float) - (A.shape[1] - 1.)/2 

548 xcen = np.sum(w[np.newaxis, :]*A)/s 

549 

550 return ycen, xcen 

551 

552 def checkCentroids(self, psfArr1, psfArr2): 

553 """Check whether two PSF array centroids' distance is within tolerance. 

554 

555 Parameters 

556 ---------- 

557 psfArr1, psfArr2 : `numpy.ndarray` of `float` 

558 Input PSF arrays to check. 

559 

560 Returns 

561 ------- 

562 None 

563 

564 Raises 

565 ------ 

566 ValueError: 

567 Centroid distance exceeds `config.maxPsfCentroidDist` pixels. 

568 """ 

569 yc1, xc1 = self.getCentroid(psfArr1) 

570 yc2, xc2 = self.getCentroid(psfArr2) 

571 dy = yc2 - yc1 

572 dx = xc2 - xc1 

573 if dy*dy + dx*dx > self.config.maxPsfCentroidDist*self.config.maxPsfCentroidDist: 

574 raise ValueError( 

575 f"PSF centroids are offset by more than {self.config.maxPsfCentroidDist:.2f} pixels.") 

576 

577 def calculateFourierDiffim(self, psf1, im1, varPlane1, F1, varMean1, 

578 psf2, im2, varPlane2, F2, varMean2, calculateScore=True): 

579 """Convolve and subtract two images in Fourier space. 

580 

581 Calculate the ZOGY proper difference image, score image and their PSFs. 

582 All input and output arrays are in Fourier space. 

583 

584 Parameters 

585 ---------- 

586 psf1, psf2, im1, im2, varPlane1, varPlane2 : `numpy.ndarray` of `numpy.complex`, 

587 shape ``self.freqSpaceShape`` 

588 Psf, image and variance plane arrays respectively. 

589 All arrays must be already in Fourier space. 

590 

591 varMean1, varMean2: `numpy.float` > 0. 

592 Average per-pixel noise variance in im1, im2 respectively. Used as weighing 

593 of input images. Must be greater than zero. 

594 

595 F1, F2 : `numpy.float` > 0. 

596 Photometric scaling of the images. See eqs. (5)--(9) 

597 

598 calculateScore : `bool`, optional 

599 If True (default), calculate and return the detection significance (score) image. 

600 Otherwise, these return fields are `None`. 

601 

602 Returns 

603 ------- 

604 result : `pipe.base.Struct` 

605 All arrays are in Fourier space and have shape ``self.freqSpaceShape``. 

606 - ``Fd`` : `float` 

607 Photometric level of ``D``. 

608 - ``D`` : `numpy.ndarray` of `numpy.complex` 

609 The difference image. 

610 - ``varplaneD`` : `numpy.ndarray` of `numpy.complex` 

611 Variance plane of ``D``. 

612 - ``Pd`` : `numpy.ndarray` of `numpy.complex` 

613 PSF of ``D``. 

614 - ``S`` : `numpy.ndarray` of `numpy.complex` or `None` 

615 Significance (score) image. 

616 - ``varplaneS`` : `numpy.ndarray` of `numpy.complex` or `None` 

617 Variance plane of ``S``. 

618 - ``Ps`` : `numpy.ndarray` of `numpy.complex` 

619 PSF of ``S``. 

620 

621 Notes 

622 ----- 

623 All array inputs and outputs are Fourier-space images with size of 

624 `self.freqSpaceShape` in this method. 

625 

626 ``varMean1``, ``varMean2`` quantities are part of the noise model and not to be confused 

627 with the variance of image frequency components or with ``varPlane1``, ``varPlane2`` that 

628 are the Fourier transform of the variance planes. 

629 """ 

630 var1F2Sq = varMean1*F2*F2 

631 var2F1Sq = varMean2*F1*F1 

632 # We need reals for comparison, also real operations are usually faster 

633 psfAbsSq1 = np.real(np.conj(psf1)*psf1) 

634 psfAbsSq2 = np.real(np.conj(psf2)*psf2) 

635 FdDenom = np.sqrt(var1F2Sq + var2F1Sq) # one number 

636 

637 # Secure positive limit to avoid floating point operations resulting in exact zero 

638 tiny = np.finfo(psf1.dtype).tiny * 100 

639 sDenom = var1F2Sq*psfAbsSq2 + var2F1Sq*psfAbsSq1 # array, eq. (12) 

640 # Frequencies where both psfs are too close to zero. 

641 # We expect this only in cases when psf1, psf2 are identical, 

642 # and either having very well sampled Gaussian tails 

643 # or having "edges" such that some sinc-like zero crossings are found at symmetry points 

644 # 

645 # if sDenom < tiny then it can be == 0. -> `denom` = 0. and 0/0 occur at `c1` , `c2` 

646 # if we keep SDenom = tiny, denom ~ O(sqrt(tiny)), Pd ~ O(sqrt(tiny)), S ~ O(sqrt(tiny)*tiny) == 0 

647 # Where S = 0 then Pd = 0 and D should still yield the same variance ~ O(1) 

648 # For safety, we set S = 0 explicitly, too, though it should be unnecessary. 

649 fltZero = sDenom < tiny 

650 nZero = np.sum(fltZero) 

651 self.log.debug(f"There are {nZero} frequencies where both FFTd PSFs are close to zero.") 

652 if nZero > 0: 

653 # We expect only a small fraction of such frequencies 

654 fltZero = np.nonzero(fltZero) # Tuple of index arrays 

655 sDenom[fltZero] = tiny # Avoid division problem but overwrite result anyway 

656 denom = np.sqrt(sDenom) # array, eq. (13) 

657 

658 c1 = F2*psf2/denom 

659 c2 = F1*psf1/denom 

660 if nZero > 0: 

661 c1[fltZero] = F2/FdDenom 

662 c2[fltZero] = F1/FdDenom 

663 D = c1*im1 - c2*im2 # Difference image eq. (13) 

664 varPlaneD = self.pixelSpaceSquare(c1)*varPlane1 + self.pixelSpaceSquare(c2)*varPlane2 # eq. (26) 

665 

666 Pd = FdDenom*psf1*psf2/denom # Psf of D eq. (14) 

667 if nZero > 0: 

668 Pd[fltZero] = 0 

669 

670 Fd = F1*F2/FdDenom # Flux scaling of D eq. (15) 

671 if calculateScore: 

672 c1 = F1*F2*F2*np.conj(psf1)*psfAbsSq2/sDenom 

673 c2 = F2*F1*F1*np.conj(psf2)*psfAbsSq1/sDenom 

674 if nZero > 0: 

675 c1[fltZero] = 0 

676 c2[fltZero] = 0 

677 S = c1*im1 - c2*im2 # eq. (12) 

678 varPlaneS = self.pixelSpaceSquare(c1)*varPlane1 + self.pixelSpaceSquare(c2)*varPlane2 

679 Ps = np.conj(Pd)*Pd # eq. (17) Source detection expects a PSF 

680 else: 

681 S = None 

682 Ps = None 

683 varPlaneS = None 

684 return pipeBase.Struct(D=D, Pd=Pd, varPlaneD=varPlaneD, Fd=Fd, 

685 S=S, Ps=Ps, varPlaneS=varPlaneS) 

686 

687 @staticmethod 

688 def calculateMaskPlane(mask1, mask2, effPsf1=None, effPsf2=None): 

689 """Calculate the mask plane of the difference image. 

690 

691 Parameters 

692 ---------- 

693 mask1, maks2 : `lsst.afw.image.Mask` 

694 Mask planes of the two exposures. 

695 

696 

697 Returns 

698 ------- 

699 diffmask : `lsst.afw.image.Mask` 

700 Mask plane for the subtraction result. 

701 

702 Notes 

703 ----- 

704 TODO DM-25174 : Specification of effPsf1, effPsf2 are not yet supported. 

705 """ 

706 

707 # mask1 x effPsf2 | mask2 x effPsf1 

708 if effPsf1 is not None or effPsf2 is not None: 

709 # TODO: DM-25174 effPsf1, effPsf2: the effective psf for cross-blurring. 

710 # We need a "size" approximation of the c1 and c2 coefficients to make effPsfs 

711 # Also convolution not yet supports mask-only operation 

712 raise NotImplementedError("Mask plane only 'convolution' operation is not yet supported") 

713 R = mask1.clone() 

714 R |= mask2 

715 return R 

716 

717 @staticmethod 

718 def makeKernelPsfFromArray(A): 

719 """Create a non spatially varying PSF from a `numpy.ndarray`. 

720 

721 Parameters 

722 ---------- 

723 A : `numpy.ndarray` 

724 2D array to use as the new psf image. The pixels are copied. 

725 

726 Returns 

727 ------- 

728 psfNew : `lsst.meas.algorithms.KernelPsf` 

729 The constructed PSF. 

730 """ 

731 psfImg = afwImage.ImageD(A.astype(np.float64, copy=True), deep=False) 

732 psfNew = measAlg.KernelPsf(afwMath.FixedKernel(psfImg)) 

733 return psfNew 

734 

735 def makeDiffimSubExposure(self, ftDiff): 

736 """Wrap array results into Exposure objects. 

737 

738 Parameters 

739 ---------- 

740 ftDiff : `lsst.pipe.base.Struct` 

741 Result struct by `calculateFourierDiffim`. 

742 

743 Returns 

744 ------- 

745 resultName : `lsst.pipe.base.Struct` 

746 - ``diffSubExp`` : `lsst.afw.image.Exposure` 

747 The difference (sub)exposure. The exposure is calibrated 

748 in its pixel values, and has a constant `PhotoCalib` object of 1. 

749 - ``scoreSubExp`` : `lsst.afw.image.Exposure` or `None` 

750 The score (sub)exposure if it was calculated. 

751 """ 

752 D = self.inverseFftAndCropImage( 

753 ftDiff.D, self.imgShape, np.logical_or(self.fftFullIm1.filtInf, self.fftFullIm2.filtInf), 

754 np.logical_or(self.fftFullIm1.filtNaN, self.fftFullIm2.filtNaN), 

755 dtype=self.subExposure1.image.dtype) 

756 varPlaneD = self.inverseFftAndCropImage( 

757 ftDiff.varPlaneD, self.imgShape, np.logical_or(self.fftVarPl1.filtInf, self.fftVarPl2.filtInf), 

758 np.logical_or(self.fftVarPl1.filtNaN, self.fftVarPl2.filtNaN), 

759 dtype=self.subExposure1.variance.dtype) 

760 Pd = self.inverseFftAndCropImage( 

761 ftDiff.Pd, self.psfShape1, dtype=self.subExpPsf1.dtype) 

762 sumPd = np.sum(Pd) 

763 # If this is smaller than 1. it is an indicator that it does not fit its original dimensions 

764 self.log.info(f"Pd sum before normalization: {sumPd:.3f}") 

765 Pd /= sumPd 

766 

767 diffSubExposure = self.subExposure1.clone() 

768 # Indices of the subexposure bbox in the full image array 

769 bbox = self.subExposure1.getBBox() 

770 xy0 = self.fullExp1.getXY0() 

771 imgD = afwImage.Image(D, deep=False, xy0=xy0, dtype=self.subExposure1.image.dtype) 

772 diffSubExposure.image = imgD[bbox] 

773 imgVarPlaneD = afwImage.Image(varPlaneD, deep=False, xy0=xy0, 

774 dtype=self.subExposure1.variance.dtype) 

775 diffSubExposure.variance = imgVarPlaneD[bbox] 

776 diffSubExposure.mask = self.calculateMaskPlane(self.subExposure1.mask, self.subExposure2.mask) 

777 

778 # Calibrate the image; subexposures must be on the same photometric scale 

779 diffSubExposure.maskedImage /= ftDiff.Fd 

780 # Now the subExposure calibration is 1. everywhere 

781 calibOne = afwImage.PhotoCalib(1.) 

782 diffSubExposure.setPhotoCalib(calibOne) 

783 # Set the PSF of this subExposure 

784 diffSubExposure.setPsf(self.makeKernelPsfFromArray(Pd)) 

785 

786 if ftDiff.S is not None: 

787 S = self.inverseFftAndCropImage( 

788 ftDiff.S, self.imgShape, np.logical_or(self.fftFullIm1.filtInf, self.fftFullIm2.filtInf), 

789 np.logical_or(self.fftFullIm1.filtNaN, self.fftFullIm2.filtNaN), 

790 dtype=self.subExposure1.image.dtype) 

791 varPlaneS = self.inverseFftAndCropImage( 

792 ftDiff.varPlaneS, self.imgShape, 

793 np.logical_or(self.fftVarPl1.filtInf, self.fftVarPl2.filtInf), 

794 np.logical_or(self.fftVarPl1.filtNaN, self.fftVarPl2.filtNaN), 

795 dtype=self.subExposure1.variance.dtype) 

796 imgS = afwImage.Image(S, deep=False, xy0=xy0, dtype=self.subExposure1.image.dtype) 

797 imgVarPlaneS = afwImage.Image(varPlaneS, deep=False, xy0=xy0, 

798 dtype=self.subExposure1.variance.dtype) 

799 imgS = imgS[bbox] 

800 imgVarPlaneS = imgVarPlaneS[bbox] 

801 

802 # Ensure that no 0/0 occur in S/sigma(S). 

803 tiny = np.finfo(varPlaneS.dtype).tiny * 10 

804 fltZero = imgVarPlaneS.array < tiny 

805 imgVarPlaneS.array[fltZero] = tiny 

806 imgS.array[fltZero] = 0 

807 

808 # PSF of S 

809 Ps = self.inverseFftAndCropImage(ftDiff.Ps, self.psfShape1, dtype=self.subExpPsf1.dtype) 

810 sumPs = np.sum(Ps) 

811 self.log.info(f"Ps sum before normalization: {sumPs:.3f}") 

812 Ps /= sumPs 

813 

814 # TODO DM-23855 : Additional score image corrections may be done here 

815 

816 scoreSubExposure = self.subExposure1.clone() 

817 scoreSubExposure.image = imgS 

818 scoreSubExposure.variance = imgVarPlaneS 

819 scoreSubExposure.mask = diffSubExposure.mask 

820 scoreSubExposure.setPhotoCalib(None) 

821 scoreSubExposure.setPsf(self.makeKernelPsfFromArray(Ps)) 

822 else: 

823 scoreSubExposure = None 

824 

825 return pipeBase.Struct(diffSubExp=diffSubExposure, scoreSubExp=scoreSubExposure) 

826 

827 def run(self, exposure1, exposure2, calculateScore=True): 

828 """Task entry point to perform the zogy subtraction 

829 of ``exposure1-exposure2``. 

830 

831 Parameters 

832 ---------- 

833 exposure1, exposure2 : `lsst.afw.image.Exposure` 

834 Two exposures warped and matched into matching pixel dimensions. 

835 calculateScore : `bool`, optional 

836 If True (default), calculate the score image and return in ``scoreExp``. 

837 

838 

839 Returns 

840 ------- 

841 resultName : `lsst.pipe.base.Struct` 

842 - ``diffExp`` : `lsst.afw.image.Exposure` 

843 The Zogy difference exposure (``exposure1-exposure2``). 

844 - ``scoreExp`` : `lsst.afw.image.Exposure` or `None` 

845 The Zogy significance or score (S) exposure if ``calculateScore==True``. 

846 - ``ftDiff`` : `lsst.pipe.base.Struct` 

847 Lower level return struct by `calculateFourierDiffim` with added 

848 fields from the task instance. For debug purposes. 

849 

850 Notes 

851 ----- 

852 

853 The score image (``S``) is defined in the ZOGY paper as the detection 

854 statistic value at each pixel. In the ZOGY image model, the input images 

855 have uniform variance noises and thus ``S`` has uniform per pixel 

856 variance (though it is not scaled to 1). In Section 3.3 of the paper, 

857 there are "corrections" defined to the score image to correct the 

858 significance values for some deviations from the image model. The first 

859 of these corrections is the calculation of the _variance plane_ of ``S`` 

860 allowing for different per pixel variance values by following the 

861 overall convolution operation on the pixels of the input images. ``S`` 

862 scaled (divided) by its corrected per pixel noise is referred as 

863 ``Scorr`` in the paper. 

864 

865 In the current implementation, ``scoreExp`` contains ``S`` in its image 

866 plane and the calculated (non-uniform) variance plane of ``S`` in its 

867 variance plane. ``scoreExp`` can be used directly for source detection 

868 as a likelihood image by respecting its variance plane or can be divided 

869 by the square root of the variance plane to scale detection significance 

870 values into units of sigma. 

871 

872 TODO DM-23855 : Implement further correction tags to the variance of 

873 ``scoreExp``. As of DM-25174 it is not determined how important these 

874 further correction tags are. 

875 

876 TODO DM-23855 : spatially varying solution on a grid is not yet implemented 

877 """ 

878 # We use the dimensions of the 1st image only in the code 

879 if exposure1.getDimensions() != exposure2.getDimensions(): 

880 raise ValueError("Exposure dimensions do not match ({} != {} )".format( 

881 exposure1.getDimensions(), exposure2.getDimensions())) 

882 

883 self.prepareFullExposure(exposure1, exposure2, correctBackground=self.config.correctBackground) 

884 

885 # TODO DM-23855: Add grid splitting support here for spatially varying PSF support 

886 # Passing exposure1,2 won't be ok here: they're not photometrically scaled. 

887 # Use the modified full maskedImages here 

888 self.prepareSubExposure() 

889 ftDiff = self.calculateFourierDiffim( 

890 self.psfFft1, self.fftFullIm1.imFft, self.fftVarPl1.imFft, self.F1, self.subExpVar1, 

891 self.psfFft2, self.fftFullIm2.imFft, self.fftVarPl2.imFft, self.F2, self.subExpVar2, 

892 calculateScore=calculateScore) 

893 diffExp = self.makeDiffimSubExposure(ftDiff) 

894 # Add debug info from the task instance 

895 ftDiff.freqSpaceShape = self.freqSpaceShape 

896 ftDiff.imgShape = self.imgShape 

897 ftDiff.psfShape1 = self.psfShape1 

898 ftDiff.psfShape2 = self.psfShape2 

899 return pipeBase.Struct(diffExp=diffExp.diffSubExp, 

900 scoreExp=diffExp.scoreSubExp, 

901 ftDiff=ftDiff) 

902 

903 

904class ZogyImagePsfMatchConfig(ImagePsfMatchConfig): 

905 """Config for the ZogyImagePsfMatchTask""" 

906 

907 zogyConfig = pexConfig.ConfigField( 

908 dtype=ZogyConfig, 

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

910 ) 

911 

912 

913class ZogyImagePsfMatchTask(ImagePsfMatchTask): 

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

915 

916 This class inherits from ImagePsfMatchTask to contain the _warper 

917 subtask and related methods. 

918 """ 

919 

920 ConfigClass = ZogyImagePsfMatchConfig 

921 

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

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

924 

925 def run(self, scienceExposure, templateExposure, doWarping=True, spatiallyVarying=False): 

926 """Register, PSF-match, and subtract two Exposures, ``scienceExposure - templateExposure`` 

927 using the ZOGY algorithm. 

928 

929 Parameters 

930 ---------- 

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

932 exposure to be warped to scienceExposure. 

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

934 reference Exposure. 

935 doWarping : `bool` 

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

937 - if True then warp templateExposure to match scienceExposure 

938 - if False then raise an Exception 

939 spatiallyVarying : `bool` 

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

941 

942 Notes 

943 ----- 

944 Do the following, in order: 

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

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

947 

948 This is the new entry point of the task as of DM-25115. 

949 

950 

951 Returns 

952 ------- 

953 results : `lsst.pipe.base.Struct` containing these fields: 

954 - subtractedExposure: `lsst.afw.image.Exposure` 

955 The subtraction result. 

956 - warpedExposure: `lsst.afw.image.Exposure` or `None` 

957 templateExposure after warping to match scienceExposure 

958 """ 

959 

960 if spatiallyVarying: 

961 raise NotImplementedError( 

962 "DM-25115 Spatially varying zogy subtraction is not implemented.") 

963 

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

965 if doWarping: 

966 self.log.info("Warping templateExposure to scienceExposure") 

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

968 scienceExposure.getWcs()) 

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

970 templateExposure = self._warper.warpExposure( 

971 scienceExposure.getWcs(), templateExposure, destBBox=scienceExposure.getBBox()) 

972 templateExposure.setPsf(psfWarped) 

973 else: 

974 raise RuntimeError("Input images are not registered. Consider setting doWarping=True.") 

975 

976 config = self.config.zogyConfig 

977 task = ZogyTask(config=config) 

978 results = task.run(scienceExposure, templateExposure) 

979 results.warpedExposure = templateExposure 

980 return results 

981 

982 def subtractExposures(self, templateExposure, scienceExposure, 

983 doWarping=True, spatiallyVarying=True, inImageSpace=False, 

984 doPreConvolve=False): 

985 raise NotImplementedError 

986 

987 def subtractMaskedImages(self, templateExposure, scienceExposure, 

988 doWarping=True, spatiallyVarying=True, inImageSpace=False, 

989 doPreConvolve=False): 

990 raise NotImplementedError 

991 

992 

993subtractAlgorithmRegistry.register('zogy', ZogyImagePsfMatchTask)