Coverage for python/lsst/cp/pipe/ptc/cpExtractPtcTask.py: 15%

Shortcuts 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

193 statements  

1# This file is part of cp_pipe. 

2# 

3# Developed for the LSST Data Management System. 

4# This product includes software developed by the LSST Project 

5# (https://www.lsst.org). 

6# See the COPYRIGHT file at the top-level directory of this distribution 

7# for details of code ownership. 

8# 

9# This program is free software: you can redistribute it and/or modify 

10# it under the terms of the GNU General Public License as published by 

11# the Free Software Foundation, either version 3 of the License, or 

12# (at your option) any later version. 

13# 

14# This program is distributed in the hope that it will be useful, 

15# but WITHOUT ANY WARRANTY; without even the implied warranty of 

16# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 

17# GNU General Public License for more details. 

18# 

19# You should have received a copy of the GNU General Public License 

20# along with this program. If not, see <https://www.gnu.org/licenses/>. 

21# 

22import numpy as np 

23 

24import lsst.afw.math as afwMath 

25import lsst.pex.config as pexConfig 

26import lsst.pipe.base as pipeBase 

27from lsst.cp.pipe.utils import (arrangeFlatsByExpTime, arrangeFlatsByExpId, 

28 sigmaClipCorrection, CovFastFourierTransform) 

29 

30import lsst.pipe.base.connectionTypes as cT 

31 

32from lsst.ip.isr import PhotonTransferCurveDataset 

33from lsst.ip.isr import IsrTask 

34 

35__all__ = ['PhotonTransferCurveExtractConfig', 'PhotonTransferCurveExtractTask'] 

36 

37 

38class PhotonTransferCurveExtractConnections(pipeBase.PipelineTaskConnections, 

39 dimensions=("instrument", "detector")): 

40 

41 inputExp = cT.Input( 

42 name="ptcInputExposurePairs", 

43 doc="Input post-ISR processed exposure pairs (flats) to" 

44 "measure covariances from.", 

45 storageClass="Exposure", 

46 dimensions=("instrument", "exposure", "detector"), 

47 multiple=True, 

48 deferLoad=True, 

49 ) 

50 

51 outputCovariances = cT.Output( 

52 name="ptcCovariances", 

53 doc="Extracted flat (co)variances.", 

54 storageClass="PhotonTransferCurveDataset", 

55 dimensions=("instrument", "exposure", "detector"), 

56 multiple=True, 

57 ) 

58 

59 

60class PhotonTransferCurveExtractConfig(pipeBase.PipelineTaskConfig, 

61 pipelineConnections=PhotonTransferCurveExtractConnections): 

62 """Configuration for the measurement of covariances from flats. 

63 """ 

64 

65 matchByExposureId = pexConfig.Field( 

66 dtype=bool, 

67 doc="Should exposures be matched by ID rather than exposure time?", 

68 default=False, 

69 ) 

70 maximumRangeCovariancesAstier = pexConfig.Field( 

71 dtype=int, 

72 doc="Maximum range of covariances as in Astier+19", 

73 default=8, 

74 ) 

75 binSize = pexConfig.Field( 

76 dtype=int, 

77 doc="Bin the image by this factor in both dimensions.", 

78 default=1, 

79 ) 

80 minMeanSignal = pexConfig.DictField( 

81 keytype=str, 

82 itemtype=float, 

83 doc="Minimum values (inclusive) of mean signal (in ADU) per amp to use." 

84 " The same cut is applied to all amps if this parameter [`dict`] is passed as " 

85 " {'ALL_AMPS': value}", 

86 default={'ALL_AMPS': 0.0}, 

87 ) 

88 maxMeanSignal = pexConfig.DictField( 

89 keytype=str, 

90 itemtype=float, 

91 doc="Maximum values (inclusive) of mean signal (in ADU) below which to consider, per amp." 

92 " The same cut is applied to all amps if this dictionary is of the form" 

93 " {'ALL_AMPS': value}", 

94 default={'ALL_AMPS': 1e6}, 

95 ) 

96 maskNameList = pexConfig.ListField( 

97 dtype=str, 

98 doc="Mask list to exclude from statistics calculations.", 

99 default=['SUSPECT', 'BAD', 'NO_DATA', 'SAT'], 

100 ) 

101 nSigmaClipPtc = pexConfig.Field( 

102 dtype=float, 

103 doc="Sigma cut for afwMath.StatisticsControl()", 

104 default=5.5, 

105 ) 

106 nIterSigmaClipPtc = pexConfig.Field( 

107 dtype=int, 

108 doc="Number of sigma-clipping iterations for afwMath.StatisticsControl()", 

109 default=3, 

110 ) 

111 minNumberGoodPixelsForCovariance = pexConfig.Field( 

112 dtype=int, 

113 doc="Minimum number of acceptable good pixels per amp to calculate the covariances (via FFT or" 

114 " direclty).", 

115 default=10000, 

116 ) 

117 thresholdDiffAfwVarVsCov00 = pexConfig.Field( 

118 dtype=float, 

119 doc="If the absolute fractional differece between afwMath.VARIANCECLIP and Cov00 " 

120 "for a region of a difference image is greater than this threshold (percentage), " 

121 "a warning will be issued.", 

122 default=1., 

123 ) 

124 detectorMeasurementRegion = pexConfig.ChoiceField( 

125 dtype=str, 

126 doc="Region of each exposure where to perform the calculations (amplifier or full image).", 

127 default='AMP', 

128 allowed={ 

129 "AMP": "Amplifier of the detector.", 

130 "FULL": "Full image." 

131 } 

132 ) 

133 numEdgeSuspect = pexConfig.Field( 

134 dtype=int, 

135 doc="Number of edge pixels to be flagged as untrustworthy.", 

136 default=0, 

137 ) 

138 edgeMaskLevel = pexConfig.ChoiceField( 

139 dtype=str, 

140 doc="Mask edge pixels in which coordinate frame: DETECTOR or AMP?", 

141 default="DETECTOR", 

142 allowed={ 

143 'DETECTOR': 'Mask only the edges of the full detector.', 

144 'AMP': 'Mask edges of each amplifier.', 

145 }, 

146 ) 

147 

148 

149class PhotonTransferCurveExtractTask(pipeBase.PipelineTask, 

150 pipeBase.CmdLineTask): 

151 """Task to measure covariances from flat fields. 

152 

153 This task receives as input a list of flat-field images 

154 (flats), and sorts these flats in pairs taken at the 

155 same time (the task will raise if there is one one flat 

156 at a given exposure time, and it will discard extra flats if 

157 there are more than two per exposure time). This task measures 

158 the mean, variance, and covariances from a region (e.g., 

159 an amplifier) of the difference image of the two flats with 

160 the same exposure time. 

161 

162 The variance is calculated via afwMath, and the covariance 

163 via the methods in Astier+19 (appendix A). In theory, 

164 var = covariance[0,0]. This should be validated, and in the 

165 future, we may decide to just keep one (covariance). 

166 At this moment, if the two values differ by more than the value 

167 of `thresholdDiffAfwVarVsCov00` (default: 1%), a warning will 

168 be issued. 

169 

170 The measured covariances at a given exposure time (along with 

171 other quantities such as the mean) are stored in a PTC dataset 

172 object (`~lsst.ip.isr.PhotonTransferCurveDataset`), which gets 

173 partially filled at this stage (the remainder of the attributes 

174 of the dataset will be filled after running the second task of 

175 the PTC-measurement pipeline, `~PhotonTransferCurveSolveTask`). 

176 

177 The number of partially-filled 

178 `~lsst.ip.isr.PhotonTransferCurveDataset` objects will be less 

179 than the number of input exposures because the task combines 

180 input flats in pairs. However, it is required at this moment 

181 that the number of input dimensions matches 

182 bijectively the number of output dimensions. Therefore, a number 

183 of "dummy" PTC datasets are inserted in the output list. This 

184 output list will then be used as input of the next task in the 

185 PTC-measurement pipeline, `PhotonTransferCurveSolveTask`, 

186 which will assemble the multiple `PhotonTransferCurveDataset` 

187 objects into a single one in order to fit the measured covariances 

188 as a function of flux to one of three models 

189 (see `PhotonTransferCurveSolveTask` for details). 

190 

191 Reference: Astier+19: "The Shape of the Photon Transfer Curve of CCD 

192 sensors", arXiv:1905.08677. 

193 """ 

194 

195 ConfigClass = PhotonTransferCurveExtractConfig 

196 _DefaultName = 'cpPtcExtract' 

197 

198 def runQuantum(self, butlerQC, inputRefs, outputRefs): 

199 """Ensure that the input and output dimensions are passed along. 

200 

201 Parameters 

202 ---------- 

203 butlerQC : `~lsst.daf.butler.butlerQuantumContext.ButlerQuantumContext` 

204 Butler to operate on. 

205 inputRefs : `~lsst.pipe.base.connections.InputQuantizedConnection` 

206 Input data refs to load. 

207 ouptutRefs : `~lsst.pipe.base.connections.OutputQuantizedConnection` 

208 Output data refs to persist. 

209 """ 

210 inputs = butlerQC.get(inputRefs) 

211 # Ids of input list of exposure references 

212 # (deferLoad=True in the input connections) 

213 

214 inputs['inputDims'] = [expRef.datasetRef.dataId['exposure'] for expRef in inputRefs.inputExp] 

215 

216 # Dictionary, keyed by expTime, with tuples containing flat 

217 # exposures and their IDs. 

218 if self.config.matchByExposureId: 

219 inputs['inputExp'] = arrangeFlatsByExpId(inputs['inputExp'], inputs['inputDims']) 

220 else: 

221 inputs['inputExp'] = arrangeFlatsByExpTime(inputs['inputExp'], inputs['inputDims']) 

222 

223 outputs = self.run(**inputs) 

224 butlerQC.put(outputs, outputRefs) 

225 

226 def run(self, inputExp, inputDims): 

227 """Measure covariances from difference of flat pairs 

228 

229 Parameters 

230 ---------- 

231 inputExp : `dict` [`float`, `list` 

232 [`~lsst.pipe.base.connections.DeferredDatasetRef`]] 

233 Dictionary that groups references to flat-field exposures that 

234 have the same exposure time (seconds), or that groups them 

235 sequentially by their exposure id. 

236 

237 inputDims : `list` 

238 List of exposure IDs. 

239 

240 Returns 

241 ------- 

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

243 The results struct containing: 

244 

245 ``outputCovariances`` 

246 A list containing the per-pair PTC measurements (`list` 

247 [`lsst.ip.isr.PhotonTransferCurveDataset`]) 

248 """ 

249 # inputExp.values() returns a view, which we turn into a list. We then 

250 # access the first exposure-ID tuple to get the detector. 

251 # The first "get()" retrieves the exposure from the exposure reference. 

252 detector = list(inputExp.values())[0][0][0].get(component='detector') 

253 detNum = detector.getId() 

254 amps = detector.getAmplifiers() 

255 ampNames = [amp.getName() for amp in amps] 

256 

257 # Each amp may have a different min and max ADU signal 

258 # specified in the config. 

259 maxMeanSignalDict = {ampName: 1e6 for ampName in ampNames} 

260 minMeanSignalDict = {ampName: 0.0 for ampName in ampNames} 

261 for ampName in ampNames: 

262 if 'ALL_AMPS' in self.config.maxMeanSignal: 

263 maxMeanSignalDict[ampName] = self.config.maxMeanSignal['ALL_AMPS'] 

264 elif ampName in self.config.maxMeanSignal: 

265 maxMeanSignalDict[ampName] = self.config.maxMeanSignal[ampName] 

266 

267 if 'ALL_AMPS' in self.config.minMeanSignal: 

268 minMeanSignalDict[ampName] = self.config.minMeanSignal['ALL_AMPS'] 

269 elif ampName in self.config.minMeanSignal: 

270 minMeanSignalDict[ampName] = self.config.minMeanSignal[ampName] 

271 # These are the column names for `tupleRows` below. 

272 tags = [('mu', '<f8'), ('afwVar', '<f8'), ('i', '<i8'), ('j', '<i8'), ('var', '<f8'), 

273 ('cov', '<f8'), ('npix', '<i8'), ('ext', '<i8'), ('expTime', '<f8'), ('ampName', '<U3')] 

274 # Create a dummy ptcDataset. Dummy datasets will be 

275 # used to ensure that the number of output and input 

276 # dimensions match. 

277 dummyPtcDataset = PhotonTransferCurveDataset(ampNames, 'DUMMY', 

278 self.config.maximumRangeCovariancesAstier) 

279 # Initialize amps of `dummyPtcDatset`. 

280 for ampName in ampNames: 

281 dummyPtcDataset.setAmpValues(ampName) 

282 # Output list with PTC datasets. 

283 partialPtcDatasetList = [] 

284 # The number of output references needs to match that of input 

285 # references: initialize outputlist with dummy PTC datasets. 

286 for i in range(len(inputDims)): 

287 partialPtcDatasetList.append(dummyPtcDataset) 

288 

289 if self.config.numEdgeSuspect > 0: 

290 isrTask = IsrTask() 

291 self.log.info("Masking %d pixels from the edges of all exposures as SUSPECT.", 

292 self.config.numEdgeSuspect) 

293 

294 for expTime in inputExp: 

295 exposures = inputExp[expTime] 

296 if len(exposures) == 1: 

297 self.log.warning("Only one exposure found at expTime %f. Dropping exposure %d.", 

298 expTime, exposures[0][1]) 

299 continue 

300 else: 

301 # Only use the first two exposures at expTime. Each 

302 # elements is a tuple (exposure, expId) 

303 expRef1, expId1 = exposures[0] 

304 expRef2, expId2 = exposures[1] 

305 # use get() to obtain `lsst.afw.image.Exposure` 

306 exp1, exp2 = expRef1.get(), expRef2.get() 

307 

308 if len(exposures) > 2: 

309 self.log.warning("Already found 2 exposures at expTime %f. Ignoring exposures: %s", 

310 expTime, ", ".join(str(i[1]) for i in exposures[2:])) 

311 # Mask pixels at the edge of the detector or of each amp 

312 if self.config.numEdgeSuspect > 0: 

313 isrTask.maskEdges(exp1, numEdgePixels=self.config.numEdgeSuspect, 

314 maskPlane="SUSPECT", level=self.config.edgeMaskLevel) 

315 isrTask.maskEdges(exp2, numEdgePixels=self.config.numEdgeSuspect, 

316 maskPlane="SUSPECT", level=self.config.edgeMaskLevel) 

317 

318 nAmpsNan = 0 

319 partialPtcDataset = PhotonTransferCurveDataset(ampNames, 'PARTIAL', 

320 self.config.maximumRangeCovariancesAstier) 

321 for ampNumber, amp in enumerate(detector): 

322 ampName = amp.getName() 

323 # covAstier: [(i, j, var (cov[0,0]), cov, npix) for 

324 # (i,j) in {maxLag, maxLag}^2] 

325 if self.config.detectorMeasurementRegion == 'AMP': 

326 region = amp.getBBox() 

327 elif self.config.detectorMeasurementRegion == 'FULL': 

328 region = None 

329 # `measureMeanVarCov` is the function that measures 

330 # the variance and covariances from a region of 

331 # the difference image of two flats at the same 

332 # exposure time. The variable `covAstier` that is 

333 # returned is of the form: 

334 # [(i, j, var (cov[0,0]), cov, npix) for (i,j) in 

335 # {maxLag, maxLag}^2]. 

336 muDiff, varDiff, covAstier = self.measureMeanVarCov(exp1, exp2, region=region) 

337 # Correction factor for bias introduced by sigma 

338 # clipping. 

339 # Function returns 1/sqrt(varFactor), so it needs 

340 # to be squared. varDiff is calculated via 

341 # afwMath.VARIANCECLIP. 

342 varFactor = sigmaClipCorrection(self.config.nSigmaClipPtc)**2 

343 varDiff *= varFactor 

344 

345 expIdMask = True 

346 # Mask data point at this mean signal level if 

347 # the signal, variance, or covariance calculations 

348 # from `measureMeanVarCov` resulted in NaNs. 

349 if np.isnan(muDiff) or np.isnan(varDiff) or (covAstier is None): 

350 self.log.warning("NaN mean or var, or None cov in amp %s in exposure pair %d, %d of " 

351 "detector %d.", ampName, expId1, expId2, detNum) 

352 nAmpsNan += 1 

353 expIdMask = False 

354 covArray = np.full((1, self.config.maximumRangeCovariancesAstier, 

355 self.config.maximumRangeCovariancesAstier), np.nan) 

356 covSqrtWeights = np.full_like(covArray, np.nan) 

357 

358 # Mask data point if it is outside of the 

359 # specified mean signal range. 

360 if (muDiff <= minMeanSignalDict[ampName]) or (muDiff >= maxMeanSignalDict[ampName]): 

361 expIdMask = False 

362 

363 if covAstier is not None: 

364 # Turn the tuples with the measured information 

365 # into covariance arrays. 

366 tupleRows = [(muDiff, varDiff) + covRow + (ampNumber, expTime, 

367 ampName) for covRow in covAstier] 

368 tempStructArray = np.array(tupleRows, dtype=tags) 

369 covArray, vcov, _ = self.makeCovArray(tempStructArray, 

370 self.config.maximumRangeCovariancesAstier) 

371 covSqrtWeights = np.nan_to_num(1./np.sqrt(vcov)) 

372 

373 # Correct covArray for sigma clipping: 

374 # 1) Apply varFactor twice for the whole covariance matrix 

375 covArray *= varFactor**2 

376 # 2) But, only once for the variance element of the 

377 # matrix, covArray[0,0] (so divide one factor out). 

378 covArray[0, 0] /= varFactor 

379 

380 partialPtcDataset.setAmpValues(ampName, rawExpTime=[expTime], rawMean=[muDiff], 

381 rawVar=[varDiff], inputExpIdPair=[(expId1, expId2)], 

382 expIdMask=[expIdMask], covArray=covArray, 

383 covSqrtWeights=covSqrtWeights) 

384 # Use location of exp1 to save PTC dataset from (exp1, exp2) pair. 

385 # Below, np.where(expId1 == np.array(inputDims)) returns a tuple 

386 # with a single-element array, so [0][0] 

387 # is necessary to extract the required index. 

388 datasetIndex = np.where(expId1 == np.array(inputDims))[0][0] 

389 # `partialPtcDatasetList` is a list of 

390 # `PhotonTransferCurveDataset` objects. Some of them 

391 # will be dummy datasets (to match length of input 

392 # and output references), and the rest will have 

393 # datasets with the mean signal, variance, and 

394 # covariance measurements at a given exposure 

395 # time. The next ppart of the PTC-measurement 

396 # pipeline, `solve`, will take this list as input, 

397 # and assemble the measurements in the datasets 

398 # in an addecuate manner for fitting a PTC 

399 # model. 

400 partialPtcDatasetList[datasetIndex] = partialPtcDataset 

401 

402 if nAmpsNan == len(ampNames): 

403 msg = f"NaN mean in all amps of exposure pair {expId1}, {expId2} of detector {detNum}." 

404 self.log.warning(msg) 

405 return pipeBase.Struct( 

406 outputCovariances=partialPtcDatasetList, 

407 ) 

408 

409 def makeCovArray(self, inputTuple, maxRangeFromTuple): 

410 """Make covariances array from tuple. 

411 

412 Parameters 

413 ---------- 

414 inputTuple : `numpy.ndarray` 

415 Structured array with rows with at least 

416 (mu, afwVar, cov, var, i, j, npix), where: 

417 

418 mu : `float` 

419 0.5*(m1 + m2), where mu1 is the mean value of flat1 

420 and mu2 is the mean value of flat2. 

421 afwVar : `float` 

422 Variance of difference flat, calculated with afw. 

423 cov : `float` 

424 Covariance value at lag(i, j) 

425 var : `float` 

426 Variance(covariance value at lag(0, 0)) 

427 i : `int` 

428 Lag in dimension "x". 

429 j : `int` 

430 Lag in dimension "y". 

431 npix : `int` 

432 Number of pixels used for covariance calculation. 

433 

434 maxRangeFromTuple : `int` 

435 Maximum range to select from tuple. 

436 

437 Returns 

438 ------- 

439 cov : `numpy.array` 

440 Covariance arrays, indexed by mean signal mu. 

441 

442 vCov : `numpy.array` 

443 Variance arrays, indexed by mean signal mu. 

444 

445 muVals : `numpy.array` 

446 List of mean signal values. 

447 """ 

448 if maxRangeFromTuple is not None: 

449 cut = (inputTuple['i'] < maxRangeFromTuple) & (inputTuple['j'] < maxRangeFromTuple) 

450 cutTuple = inputTuple[cut] 

451 else: 

452 cutTuple = inputTuple 

453 # increasing mu order, so that we can group measurements with the 

454 # same mu 

455 muTemp = cutTuple['mu'] 

456 ind = np.argsort(muTemp) 

457 

458 cutTuple = cutTuple[ind] 

459 # should group measurements on the same image pairs(same average) 

460 mu = cutTuple['mu'] 

461 xx = np.hstack(([mu[0]], mu)) 

462 delta = xx[1:] - xx[:-1] 

463 steps, = np.where(delta > 0) 

464 ind = np.zeros_like(mu, dtype=int) 

465 ind[steps] = 1 

466 ind = np.cumsum(ind) # this acts as an image pair index. 

467 # now fill the 3-d cov array(and variance) 

468 muVals = np.array(np.unique(mu)) 

469 i = cutTuple['i'].astype(int) 

470 j = cutTuple['j'].astype(int) 

471 c = 0.5*cutTuple['cov'] 

472 n = cutTuple['npix'] 

473 v = 0.5*cutTuple['var'] 

474 # book and fill 

475 cov = np.ndarray((len(muVals), np.max(i)+1, np.max(j)+1)) 

476 var = np.zeros_like(cov) 

477 cov[ind, i, j] = c 

478 var[ind, i, j] = v**2/n 

479 var[:, 0, 0] *= 2 # var(v) = 2*v**2/N 

480 

481 return cov, var, muVals 

482 

483 def measureMeanVarCov(self, exposure1, exposure2, region=None): 

484 """Calculate the mean of each of two exposures and the variance 

485 and covariance of their difference. The variance is calculated 

486 via afwMath, and the covariance via the methods in Astier+19 

487 (appendix A). In theory, var = covariance[0,0]. This should 

488 be validated, and in the future, we may decide to just keep 

489 one (covariance). 

490 

491 Parameters 

492 ---------- 

493 exposure1 : `lsst.afw.image.exposure.ExposureF` 

494 First exposure of flat field pair. 

495 exposure2 : `lsst.afw.image.exposure.ExposureF` 

496 Second exposure of flat field pair. 

497 region : `lsst.geom.Box2I`, optional 

498 Region of each exposure where to perform the calculations 

499 (e.g, an amplifier). 

500 

501 Returns 

502 ------- 

503 mu : `float` or `NaN` 

504 0.5*(mu1 + mu2), where mu1, and mu2 are the clipped means 

505 of the regions in both exposures. If either mu1 or m2 are 

506 NaN's, the returned value is NaN. 

507 varDiff : `float` or `NaN` 

508 Half of the clipped variance of the difference of the 

509 regions inthe two input exposures. If either mu1 or m2 are 

510 NaN's, the returned value is NaN. 

511 covDiffAstier : `list` or `NaN` 

512 List with tuples of the form (dx, dy, var, cov, npix), where: 

513 dx : `int` 

514 Lag in x 

515 dy : `int` 

516 Lag in y 

517 var : `float` 

518 Variance at (dx, dy). 

519 cov : `float` 

520 Covariance at (dx, dy). 

521 nPix : `int` 

522 Number of pixel pairs used to evaluate var and cov. 

523 

524 If either mu1 or m2 are NaN's, the returned value is NaN. 

525 """ 

526 if region is not None: 

527 im1Area = exposure1.maskedImage[region] 

528 im2Area = exposure2.maskedImage[region] 

529 else: 

530 im1Area = exposure1.maskedImage 

531 im2Area = exposure2.maskedImage 

532 

533 if self.config.binSize > 1: 

534 im1Area = afwMath.binImage(im1Area, self.config.binSize) 

535 im2Area = afwMath.binImage(im2Area, self.config.binSize) 

536 

537 im1MaskVal = exposure1.getMask().getPlaneBitMask(self.config.maskNameList) 

538 im1StatsCtrl = afwMath.StatisticsControl(self.config.nSigmaClipPtc, 

539 self.config.nIterSigmaClipPtc, 

540 im1MaskVal) 

541 im1StatsCtrl.setNanSafe(True) 

542 im1StatsCtrl.setAndMask(im1MaskVal) 

543 

544 im2MaskVal = exposure2.getMask().getPlaneBitMask(self.config.maskNameList) 

545 im2StatsCtrl = afwMath.StatisticsControl(self.config.nSigmaClipPtc, 

546 self.config.nIterSigmaClipPtc, 

547 im2MaskVal) 

548 im2StatsCtrl.setNanSafe(True) 

549 im2StatsCtrl.setAndMask(im2MaskVal) 

550 

551 # Clipped mean of images; then average of mean. 

552 mu1 = afwMath.makeStatistics(im1Area, afwMath.MEANCLIP, im1StatsCtrl).getValue() 

553 mu2 = afwMath.makeStatistics(im2Area, afwMath.MEANCLIP, im2StatsCtrl).getValue() 

554 if np.isnan(mu1) or np.isnan(mu2): 

555 self.log.warning("Mean of amp in image 1 or 2 is NaN: %f, %f.", mu1, mu2) 

556 return np.nan, np.nan, None 

557 mu = 0.5*(mu1 + mu2) 

558 

559 # Take difference of pairs 

560 # symmetric formula: diff = (mu2*im1-mu1*im2)/(0.5*(mu1+mu2)) 

561 temp = im2Area.clone() 

562 temp *= mu1 

563 diffIm = im1Area.clone() 

564 diffIm *= mu2 

565 diffIm -= temp 

566 diffIm /= mu 

567 

568 diffImMaskVal = diffIm.getMask().getPlaneBitMask(self.config.maskNameList) 

569 diffImStatsCtrl = afwMath.StatisticsControl(self.config.nSigmaClipPtc, 

570 self.config.nIterSigmaClipPtc, 

571 diffImMaskVal) 

572 diffImStatsCtrl.setNanSafe(True) 

573 diffImStatsCtrl.setAndMask(diffImMaskVal) 

574 

575 # Variance calculation via afwMath 

576 varDiff = 0.5*(afwMath.makeStatistics(diffIm, afwMath.VARIANCECLIP, diffImStatsCtrl).getValue()) 

577 

578 # Covariances calculations 

579 # Get the pixels that were not clipped 

580 varClip = afwMath.makeStatistics(diffIm, afwMath.VARIANCECLIP, diffImStatsCtrl).getValue() 

581 meanClip = afwMath.makeStatistics(diffIm, afwMath.MEANCLIP, diffImStatsCtrl).getValue() 

582 cut = meanClip + self.config.nSigmaClipPtc*np.sqrt(varClip) 

583 unmasked = np.where(np.fabs(diffIm.image.array) <= cut, 1, 0) 

584 

585 # Get the pixels in the mask planes of the difference image 

586 # that were ignored by the clipping algorithm 

587 wDiff = np.where(diffIm.getMask().getArray() == 0, 1, 0) 

588 # Combine the two sets of pixels ('1': use; '0': don't use) 

589 # into a final weight matrix to be used in the covariance 

590 # calculations below. 

591 w = unmasked*wDiff 

592 

593 if np.sum(w) < self.config.minNumberGoodPixelsForCovariance: 

594 self.log.warning("Number of good points for covariance calculation (%s) is less " 

595 "(than threshold %s)", np.sum(w), self.config.minNumberGoodPixelsForCovariance) 

596 return np.nan, np.nan, None 

597 

598 maxRangeCov = self.config.maximumRangeCovariancesAstier 

599 

600 # Calculate covariances via FFT. 

601 shapeDiff = np.array(diffIm.image.array.shape) 

602 # Calculate the sizes of FFT dimensions. 

603 s = shapeDiff + maxRangeCov 

604 tempSize = np.array(np.log(s)/np.log(2.)).astype(int) 

605 fftSize = np.array(2**(tempSize+1)).astype(int) 

606 fftShape = (fftSize[0], fftSize[1]) 

607 

608 c = CovFastFourierTransform(diffIm.image.array, w, fftShape, maxRangeCov) 

609 covDiffAstier = c.reportCovFastFourierTransform(maxRangeCov) 

610 

611 # Compare Cov[0,0] and afwMath.VARIANCECLIP covDiffAstier[0] 

612 # is the Cov[0,0] element, [3] is the variance, and there's a 

613 # factor of 0.5 difference with afwMath.VARIANCECLIP. 

614 thresholdPercentage = self.config.thresholdDiffAfwVarVsCov00 

615 fractionalDiff = 100*np.fabs(1 - varDiff/(covDiffAstier[0][3]*0.5)) 

616 if fractionalDiff >= thresholdPercentage: 

617 self.log.warning("Absolute fractional difference between afwMatch.VARIANCECLIP and Cov[0,0] " 

618 "is more than %f%%: %f", thresholdPercentage, fractionalDiff) 

619 

620 return mu, varDiff, covDiffAstier