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

355 statements  

« prev     ^ index     » next       coverage.py v7.4.4, created at 2024-03-20 05:19 -0700

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 

23from lmfit.models import GaussianModel 

24import scipy.stats 

25import warnings 

26 

27import lsst.afw.math as afwMath 

28import lsst.pex.config as pexConfig 

29import lsst.pipe.base as pipeBase 

30from lsst.geom import (Box2I, Point2I, Extent2I) 

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

32 arrangeFlatsByExpFlux, sigmaClipCorrection, 

33 CovFastFourierTransform) 

34 

35import lsst.pipe.base.connectionTypes as cT 

36 

37from lsst.ip.isr import PhotonTransferCurveDataset 

38from lsst.ip.isr import IsrTask 

39 

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

41 

42 

43class PhotonTransferCurveExtractConnections(pipeBase.PipelineTaskConnections, 

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

45 

46 inputExp = cT.Input( 

47 name="ptcInputExposurePairs", 

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

49 "measure covariances from.", 

50 storageClass="Exposure", 

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

52 multiple=True, 

53 deferLoad=True, 

54 ) 

55 inputPhotodiodeData = cT.Input( 

56 name="photodiode", 

57 doc="Photodiode readings data.", 

58 storageClass="IsrCalib", 

59 dimensions=("instrument", "exposure"), 

60 multiple=True, 

61 deferLoad=True, 

62 ) 

63 taskMetadata = cT.Input( 

64 name="isr_metadata", 

65 doc="Input task metadata to extract statistics from.", 

66 storageClass="TaskMetadata", 

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

68 multiple=True, 

69 ) 

70 outputCovariances = cT.Output( 

71 name="ptcCovariances", 

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

73 storageClass="PhotonTransferCurveDataset", 

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

75 isCalibration=True, 

76 multiple=True, 

77 ) 

78 

79 def __init__(self, *, config=None): 

80 if not config.doExtractPhotodiodeData: 

81 del self.inputPhotodiodeData 

82 

83 

84class PhotonTransferCurveExtractConfig(pipeBase.PipelineTaskConfig, 

85 pipelineConnections=PhotonTransferCurveExtractConnections): 

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

87 """ 

88 matchExposuresType = pexConfig.ChoiceField( 

89 dtype=str, 

90 doc="Match input exposures by time, flux, or expId", 

91 default='TIME', 

92 allowed={ 

93 "TIME": "Match exposures by exposure time.", 

94 "FLUX": "Match exposures by target flux. Use header keyword" 

95 " in matchExposuresByFluxKeyword to find the flux.", 

96 "EXPID": "Match exposures by exposure ID." 

97 } 

98 ) 

99 matchExposuresByFluxKeyword = pexConfig.Field( 

100 dtype=str, 

101 doc="Header keyword for flux if matchExposuresType is FLUX.", 

102 default='CCOBFLUX', 

103 ) 

104 maximumRangeCovariancesAstier = pexConfig.Field( 

105 dtype=int, 

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

107 default=8, 

108 ) 

109 binSize = pexConfig.Field( 

110 dtype=int, 

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

112 default=1, 

113 ) 

114 minMeanSignal = pexConfig.DictField( 

115 keytype=str, 

116 itemtype=float, 

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

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

119 " {'ALL_AMPS': value}", 

120 default={'ALL_AMPS': 0.0}, 

121 deprecated="This config has been moved to cpSolvePtcTask, and will be removed after v26.", 

122 ) 

123 maxMeanSignal = pexConfig.DictField( 

124 keytype=str, 

125 itemtype=float, 

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

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

128 " {'ALL_AMPS': value}", 

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

130 deprecated="This config has been moved to cpSolvePtcTask, and will be removed after v26.", 

131 ) 

132 maskNameList = pexConfig.ListField( 

133 dtype=str, 

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

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

136 ) 

137 nSigmaClipPtc = pexConfig.Field( 

138 dtype=float, 

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

140 default=5.5, 

141 ) 

142 nIterSigmaClipPtc = pexConfig.Field( 

143 dtype=int, 

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

145 default=3, 

146 ) 

147 minNumberGoodPixelsForCovariance = pexConfig.Field( 

148 dtype=int, 

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

150 " direclty).", 

151 default=10000, 

152 ) 

153 thresholdDiffAfwVarVsCov00 = pexConfig.Field( 

154 dtype=float, 

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

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

157 "a warning will be issued.", 

158 default=1., 

159 ) 

160 detectorMeasurementRegion = pexConfig.ChoiceField( 

161 dtype=str, 

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

163 default='AMP', 

164 allowed={ 

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

166 "FULL": "Full image." 

167 } 

168 ) 

169 numEdgeSuspect = pexConfig.Field( 

170 dtype=int, 

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

172 default=0, 

173 ) 

174 edgeMaskLevel = pexConfig.ChoiceField( 

175 dtype=str, 

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

177 default="DETECTOR", 

178 allowed={ 

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

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

181 }, 

182 ) 

183 doGain = pexConfig.Field( 

184 dtype=bool, 

185 doc="Calculate a gain per input flat pair.", 

186 default=True, 

187 ) 

188 gainCorrectionType = pexConfig.ChoiceField( 

189 dtype=str, 

190 doc="Correction type for the gain.", 

191 default='FULL', 

192 allowed={ 

193 'NONE': 'No correction.', 

194 'SIMPLE': 'First order correction.', 

195 'FULL': 'Second order correction.' 

196 } 

197 ) 

198 ksHistNBins = pexConfig.Field( 

199 dtype=int, 

200 doc="Number of bins for the KS test histogram.", 

201 default=100, 

202 ) 

203 ksHistLimitMultiplier = pexConfig.Field( 

204 dtype=float, 

205 doc="Number of sigma (as predicted from the mean value) to compute KS test histogram.", 

206 default=8.0, 

207 ) 

208 ksHistMinDataValues = pexConfig.Field( 

209 dtype=int, 

210 doc="Minimum number of good data values to compute KS test histogram.", 

211 default=100, 

212 ) 

213 auxiliaryHeaderKeys = pexConfig.ListField( 

214 dtype=str, 

215 doc="Auxiliary header keys to store with the PTC dataset.", 

216 default=[], 

217 ) 

218 doExtractPhotodiodeData = pexConfig.Field( 

219 dtype=bool, 

220 doc="Extract photodiode data?", 

221 default=False, 

222 ) 

223 photodiodeIntegrationMethod = pexConfig.ChoiceField( 

224 dtype=str, 

225 doc="Integration method for photodiode monitoring data.", 

226 default="CHARGE_SUM", 

227 allowed={ 

228 "DIRECT_SUM": ("Use numpy's trapz integrator on all photodiode " 

229 "readout entries"), 

230 "TRIMMED_SUM": ("Use numpy's trapz integrator, clipping the " 

231 "leading and trailing entries, which are " 

232 "nominally at zero baseline level."), 

233 "CHARGE_SUM": ("Treat the current values as integrated charge " 

234 "over the sampling interval and simply sum " 

235 "the values, after subtracting a baseline level."), 

236 }, 

237 ) 

238 photodiodeCurrentScale = pexConfig.Field( 

239 dtype=float, 

240 doc="Scale factor to apply to photodiode current values for the " 

241 "``CHARGE_SUM`` integration method.", 

242 default=-1.0, 

243 ) 

244 

245 

246class PhotonTransferCurveExtractTask(pipeBase.PipelineTask): 

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

248 

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

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

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

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

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

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

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

256 the same exposure time (alternatively, all input images could have 

257 the same exposure time but their flux changed). 

258 

259 The variance is calculated via afwMath, and the covariance 

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

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

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

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

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

265 be issued. 

266 

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

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

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

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

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

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

273 

274 The number of partially-filled 

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

276 than the number of input exposures because the task combines 

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

278 that the number of input dimensions matches 

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

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

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

282 PTC-measurement pipeline, `PhotonTransferCurveSolveTask`, 

283 which will assemble the multiple `PhotonTransferCurveDataset` 

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

285 as a function of flux to one of three models 

286 (see `PhotonTransferCurveSolveTask` for details). 

287 

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

289 sensors", arXiv:1905.08677. 

290 """ 

291 

292 ConfigClass = PhotonTransferCurveExtractConfig 

293 _DefaultName = 'cpPtcExtract' 

294 

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

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

297 

298 Parameters 

299 ---------- 

300 butlerQC : `~lsst.daf.butler.QuantumContext` 

301 Butler to operate on. 

302 inputRefs : `~lsst.pipe.base.InputQuantizedConnection` 

303 Input data refs to load. 

304 ouptutRefs : `~lsst.pipe.base.OutputQuantizedConnection` 

305 Output data refs to persist. 

306 """ 

307 inputs = butlerQC.get(inputRefs) 

308 # Ids of input list of exposure references 

309 # (deferLoad=True in the input connections) 

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

311 

312 # Dictionary, keyed by expTime (or expFlux or expId), with tuples 

313 # containing flat exposures and their IDs. 

314 matchType = self.config.matchExposuresType 

315 if matchType == 'TIME': 

316 inputs['inputExp'] = arrangeFlatsByExpTime(inputs['inputExp'], inputs['inputDims'], log=self.log) 

317 elif matchType == 'FLUX': 

318 inputs['inputExp'] = arrangeFlatsByExpFlux( 

319 inputs['inputExp'], 

320 inputs['inputDims'], 

321 self.config.matchExposuresByFluxKeyword, 

322 log=self.log, 

323 ) 

324 else: 

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

326 

327 outputs = self.run(**inputs) 

328 outputs = self._guaranteeOutputs(inputs['inputDims'], outputs, outputRefs) 

329 butlerQC.put(outputs, outputRefs) 

330 

331 def _guaranteeOutputs(self, inputDims, outputs, outputRefs): 

332 """Ensure that all outputRefs have a matching output, and if they do 

333 not, fill the output with dummy PTC datasets. 

334 

335 Parameters 

336 ---------- 

337 inputDims : `dict` [`str`, `int`] 

338 Input exposure dimensions. 

339 outputs : `lsst.pipe.base.Struct` 

340 Outputs from the ``run`` method. Contains the entry: 

341 

342 ``outputCovariances`` 

343 Output PTC datasets (`list` [`lsst.ip.isr.IsrCalib`]) 

344 outputRefs : `~lsst.pipe.base.OutputQuantizedConnection` 

345 Container with all of the outputs expected to be generated. 

346 

347 Returns 

348 ------- 

349 outputs : `lsst.pipe.base.Struct` 

350 Dummy dataset padded version of the input ``outputs`` with 

351 the same entries. 

352 """ 

353 newCovariances = [] 

354 for ref in outputRefs.outputCovariances: 

355 outputExpId = ref.dataId['exposure'] 

356 if outputExpId in inputDims: 

357 entry = inputDims.index(outputExpId) 

358 newCovariances.append(outputs.outputCovariances[entry]) 

359 else: 

360 newPtc = PhotonTransferCurveDataset(['no amp'], 'DUMMY', covMatrixSide=1) 

361 newPtc.setAmpValuesPartialDataset('no amp') 

362 newCovariances.append(newPtc) 

363 return pipeBase.Struct(outputCovariances=newCovariances) 

364 

365 def run(self, inputExp, inputDims, taskMetadata, inputPhotodiodeData=None): 

366 

367 """Measure covariances from difference of flat pairs 

368 

369 Parameters 

370 ---------- 

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

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

373 Dictionary that groups references to flat-field exposures that 

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

375 sequentially by their exposure id. 

376 inputDims : `list` 

377 List of exposure IDs. 

378 taskMetadata : `list` [`lsst.pipe.base.TaskMetadata`] 

379 List of exposures metadata from ISR. 

380 inputPhotodiodeData : `dict` [`str`, `lsst.ip.isr.PhotodiodeCalib`] 

381 Photodiode readings data (optional). 

382 

383 Returns 

384 ------- 

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

386 The resulting Struct contains: 

387 

388 ``outputCovariances`` 

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

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

391 """ 

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

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

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

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

396 detNum = detector.getId() 

397 amps = detector.getAmplifiers() 

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

399 

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

401 # specified in the config. 

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

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

404 for ampName in ampNames: 

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

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

407 elif ampName in self.config.maxMeanSignal: 

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

409 

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

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

412 elif ampName in self.config.minMeanSignal: 

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

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

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

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

417 # Create a dummy ptcDataset. Dummy datasets will be 

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

419 # dimensions match. 

420 dummyPtcDataset = PhotonTransferCurveDataset( 

421 ampNames, 'DUMMY', 

422 covMatrixSide=self.config.maximumRangeCovariancesAstier) 

423 for ampName in ampNames: 

424 dummyPtcDataset.setAmpValuesPartialDataset(ampName) 

425 

426 # Extract the photodiode data if requested. 

427 if self.config.doExtractPhotodiodeData: 

428 # Compute the photodiode integrals once, at the start. 

429 monitorDiodeCharge = {} 

430 for handle in inputPhotodiodeData: 

431 expId = handle.dataId['exposure'] 

432 pdCalib = handle.get() 

433 pdCalib.integrationMethod = self.config.photodiodeIntegrationMethod 

434 pdCalib.currentScale = self.config.photodiodeCurrentScale 

435 monitorDiodeCharge[expId] = pdCalib.integrate() 

436 

437 # Output list with PTC datasets. 

438 partialPtcDatasetList = [] 

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

440 # references: initialize outputlist with dummy PTC datasets. 

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

442 partialPtcDatasetList.append(dummyPtcDataset) 

443 

444 if self.config.numEdgeSuspect > 0: 

445 isrTask = IsrTask() 

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

447 self.config.numEdgeSuspect, self.config.edgeMaskLevel) 

448 

449 # Depending on the value of config.matchExposuresType 

450 # 'expTime' can stand for exposure time, flux, or ID. 

451 for expTime in inputExp: 

452 exposures = inputExp[expTime] 

453 if not np.isfinite(expTime): 

454 self.log.warning("Illegal/missing %s found (%s). Dropping exposure %d", 

455 self.config.matchExposuresType, expTime, exposures[0][1]) 

456 continue 

457 elif len(exposures) == 1: 

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

459 self.config.matchExposuresType, expTime, exposures[0][1]) 

460 continue 

461 else: 

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

463 # element is a tuple (exposure, expId) 

464 expRef1, expId1 = exposures[0] 

465 expRef2, expId2 = exposures[1] 

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

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

468 

469 if len(exposures) > 2: 

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

471 self.config.matchExposuresType, expTime, 

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

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

474 if self.config.numEdgeSuspect > 0: 

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

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

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

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

479 

480 # Extract any metadata keys from the headers. 

481 auxDict = {} 

482 metadata = exp1.getMetadata() 

483 for key in self.config.auxiliaryHeaderKeys: 

484 if key not in metadata: 

485 self.log.warning( 

486 "Requested auxiliary keyword %s not found in exposure metadata for %d", 

487 key, 

488 expId1, 

489 ) 

490 value = np.nan 

491 else: 

492 value = metadata[key] 

493 

494 auxDict[key] = value 

495 

496 nAmpsNan = 0 

497 partialPtcDataset = PhotonTransferCurveDataset( 

498 ampNames, 'PARTIAL', 

499 covMatrixSide=self.config.maximumRangeCovariancesAstier) 

500 for ampNumber, amp in enumerate(detector): 

501 ampName = amp.getName() 

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

503 region = amp.getBBox() 

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

505 region = None 

506 

507 # Get masked image regions, masking planes, statistic control 

508 # objects, and clipped means. Calculate once to reuse in 

509 # `measureMeanVarCov` and `getGainFromFlatPair`. 

510 im1Area, im2Area, imStatsCtrl, mu1, mu2 = self.getImageAreasMasksStats(exp1, exp2, 

511 region=region) 

512 

513 # Get the read noise for each exposure 

514 readNoise1 = dict() 

515 readNoise2 = dict() 

516 meanReadNoise = dict() 

517 

518 expMetadata1 = expRef1.get(component="metadata") 

519 metadataIndex1 = inputDims.index(expId1) 

520 thisTaskMetadata1 = taskMetadata[metadataIndex1] 

521 

522 expMetadata2 = expRef2.get(component="metadata") 

523 metadataIndex2 = inputDims.index(expId2) 

524 thisTaskMetadata2 = taskMetadata[metadataIndex2] 

525 

526 readNoise1[ampName] = self.getReadNoise(expMetadata1, thisTaskMetadata1, ampName) 

527 readNoise2[ampName] = self.getReadNoise(expMetadata2, thisTaskMetadata2, ampName) 

528 

529 meanReadNoise[ampName] = np.nanmean([readNoise1[ampName], readNoise2[ampName]]) 

530 

531 # We demand that both mu1 and mu2 be finite and greater than 0. 

532 if not np.isfinite(mu1) or not np.isfinite(mu2) \ 

533 or ((np.nan_to_num(mu1) + np.nan_to_num(mu2)/2.) <= 0.0): 

534 self.log.warning( 

535 "Illegal mean value(s) detected for amp %s on exposure pair %d/%d", 

536 ampName, 

537 expId1, 

538 expId2, 

539 ) 

540 partialPtcDataset.setAmpValuesPartialDataset( 

541 ampName, 

542 inputExpIdPair=(expId1, expId2), 

543 rawExpTime=expTime, 

544 expIdMask=False, 

545 ) 

546 continue 

547 

548 # `measureMeanVarCov` is the function that measures 

549 # the variance and covariances from a region of 

550 # the difference image of two flats at the same 

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

552 # returned is of the form: 

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

554 # {maxLag, maxLag}^2]. 

555 muDiff, varDiff, covAstier, rowMeanVariance = self.measureMeanVarCov(im1Area, 

556 im2Area, 

557 imStatsCtrl, 

558 mu1, 

559 mu2) 

560 # Estimate the gain from the flat pair 

561 if self.config.doGain: 

562 gain = self.getGainFromFlatPair(im1Area, im2Area, imStatsCtrl, mu1, mu2, 

563 correctionType=self.config.gainCorrectionType, 

564 readNoise=meanReadNoise[ampName]) 

565 else: 

566 gain = np.nan 

567 

568 # Correction factor for bias introduced by sigma 

569 # clipping. 

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

571 # to be squared. varDiff is calculated via 

572 # afwMath.VARIANCECLIP. 

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

574 varDiff *= varFactor 

575 

576 expIdMask = True 

577 # Mask data point at this mean signal level if 

578 # the signal, variance, or covariance calculations 

579 # from `measureMeanVarCov` resulted in NaNs. 

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

581 or (rowMeanVariance is np.nan)): 

582 self.log.warning("NaN mean, var or rowmeanVariance, or None cov in amp %s " 

583 "in exposure pair %d, %d of " 

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

585 nAmpsNan += 1 

586 expIdMask = False 

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

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

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

590 

591 # Mask data point if it is outside of the 

592 # specified mean signal range. 

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

594 expIdMask = False 

595 

596 if covAstier is not None: 

597 # Turn the tuples with the measured information 

598 # into covariance arrays. 

599 # covrow: (i, j, var (cov[0,0]), cov, npix) 

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

601 ampName) for covRow in covAstier] 

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

603 

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

605 self.config.maximumRangeCovariancesAstier) 

606 

607 # The returned covArray should only have 1 entry; 

608 # raise if this is not the case. 

609 if covArray.shape[0] != 1: 

610 raise RuntimeError("Serious programming error in covArray shape.") 

611 

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

613 

614 # Correct covArray for sigma clipping: 

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

616 covArray *= varFactor**2 

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

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

619 # (the first 0 is because this is a 3D array for insertion into 

620 # the combined dataset). 

621 covArray[0, 0, 0] /= varFactor 

622 

623 if expIdMask: 

624 # Run the Gaussian histogram only if this is a legal 

625 # amplifier. 

626 histVar, histChi2Dof, kspValue = self.computeGaussianHistogramParameters( 

627 im1Area, 

628 im2Area, 

629 imStatsCtrl, 

630 mu1, 

631 mu2, 

632 ) 

633 else: 

634 histVar = np.nan 

635 histChi2Dof = np.nan 

636 kspValue = 0.0 

637 

638 if self.config.doExtractPhotodiodeData: 

639 nExps = 0 

640 photoCharge = 0.0 

641 for expId in [expId1, expId2]: 

642 if expId in monitorDiodeCharge: 

643 photoCharge += monitorDiodeCharge[expId] 

644 nExps += 1 

645 if nExps > 0: 

646 photoCharge /= nExps 

647 else: 

648 photoCharge = np.nan 

649 else: 

650 photoCharge = np.nan 

651 

652 partialPtcDataset.setAmpValuesPartialDataset( 

653 ampName, 

654 inputExpIdPair=(expId1, expId2), 

655 rawExpTime=expTime, 

656 rawMean=muDiff, 

657 rawVar=varDiff, 

658 photoCharge=photoCharge, 

659 expIdMask=expIdMask, 

660 covariance=covArray[0, :, :], 

661 covSqrtWeights=covSqrtWeights[0, :, :], 

662 gain=gain, 

663 noise=meanReadNoise[ampName], 

664 histVar=histVar, 

665 histChi2Dof=histChi2Dof, 

666 kspValue=kspValue, 

667 rowMeanVariance=rowMeanVariance, 

668 ) 

669 

670 partialPtcDataset.setAuxValuesPartialDataset(auxDict) 

671 

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

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

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

675 # is necessary to extract the required index. 

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

677 # `partialPtcDatasetList` is a list of 

678 # `PhotonTransferCurveDataset` objects. Some of them 

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

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

681 # datasets with the mean signal, variance, and 

682 # covariance measurements at a given exposure 

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

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

685 # and assemble the measurements in the datasets 

686 # in an addecuate manner for fitting a PTC 

687 # model. 

688 partialPtcDataset.updateMetadataFromExposures([exp1, exp2]) 

689 partialPtcDataset.updateMetadata(setDate=True, detector=detector) 

690 partialPtcDatasetList[datasetIndex] = partialPtcDataset 

691 

692 if nAmpsNan == len(ampNames): 

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

694 self.log.warning(msg) 

695 

696 return pipeBase.Struct( 

697 outputCovariances=partialPtcDatasetList, 

698 ) 

699 

700 def makeCovArray(self, inputTuple, maxRangeFromTuple): 

701 """Make covariances array from tuple. 

702 

703 Parameters 

704 ---------- 

705 inputTuple : `numpy.ndarray` 

706 Structured array with rows with at least 

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

708 mu : `float` 

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

710 and mu2 is the mean value of flat2. 

711 afwVar : `float` 

712 Variance of difference flat, calculated with afw. 

713 cov : `float` 

714 Covariance value at lag(i, j) 

715 var : `float` 

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

717 i : `int` 

718 Lag in dimension "x". 

719 j : `int` 

720 Lag in dimension "y". 

721 npix : `int` 

722 Number of pixels used for covariance calculation. 

723 maxRangeFromTuple : `int` 

724 Maximum range to select from tuple. 

725 

726 Returns 

727 ------- 

728 cov : `numpy.array` 

729 Covariance arrays, indexed by mean signal mu. 

730 vCov : `numpy.array` 

731 Variance of the [co]variance arrays, indexed by mean signal mu. 

732 muVals : `numpy.array` 

733 List of mean signal values. 

734 """ 

735 if maxRangeFromTuple is not None: 

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

737 cutTuple = inputTuple[cut] 

738 else: 

739 cutTuple = inputTuple 

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

741 # same mu 

742 muTemp = cutTuple['mu'] 

743 ind = np.argsort(muTemp) 

744 

745 cutTuple = cutTuple[ind] 

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

747 mu = cutTuple['mu'] 

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

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

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

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

752 ind[steps] = 1 

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

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

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

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

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

758 c = 0.5*cutTuple['cov'] 

759 n = cutTuple['npix'] 

760 v = 0.5*cutTuple['var'] 

761 # book and fill 

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

763 var = np.zeros_like(cov) 

764 cov[ind, i, j] = c 

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

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

767 

768 return cov, var, muVals 

769 

770 def measureMeanVarCov(self, im1Area, im2Area, imStatsCtrl, mu1, mu2): 

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

772 and covariance of their difference. The variance is calculated 

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

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

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

776 one (covariance). 

777 

778 Parameters 

779 ---------- 

780 im1Area : `lsst.afw.image.maskedImage.MaskedImageF` 

781 Masked image from exposure 1. 

782 im2Area : `lsst.afw.image.maskedImage.MaskedImageF` 

783 Masked image from exposure 2. 

784 imStatsCtrl : `lsst.afw.math.StatisticsControl` 

785 Statistics control object. 

786 mu1: `float` 

787 Clipped mean of im1Area (ADU). 

788 mu2: `float` 

789 Clipped mean of im2Area (ADU). 

790 

791 Returns 

792 ------- 

793 mu : `float` 

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

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

796 NaN's, the returned value is NaN. 

797 varDiff : `float` 

798 Half of the clipped variance of the difference of the 

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

800 NaN's, the returned value is NaN. 

801 covDiffAstier : `list` or `NaN` 

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

803 dx : `int` 

804 Lag in x 

805 dy : `int` 

806 Lag in y 

807 var : `float` 

808 Variance at (dx, dy). 

809 cov : `float` 

810 Covariance at (dx, dy). 

811 nPix : `int` 

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

813 rowMeanVariance : `float` 

814 Variance of the means of each row in the difference image. 

815 Taken from `github.com/lsst-camera-dh/eo_pipe`. 

816 

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

818 """ 

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

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

821 return np.nan, np.nan, None, np.nan 

822 mu = 0.5*(mu1 + mu2) 

823 

824 # Take difference of pairs 

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

826 temp = im2Area.clone() 

827 temp *= mu1 

828 diffIm = im1Area.clone() 

829 diffIm *= mu2 

830 diffIm -= temp 

831 diffIm /= mu 

832 

833 if self.config.binSize > 1: 

834 diffIm = afwMath.binImage(diffIm, self.config.binSize) 

835 

836 # Calculate the variance (in ADU^2) of the means of rows for diffIm. 

837 # Taken from eo_pipe 

838 region = diffIm.getBBox() 

839 rowMeans = [] 

840 for row in range(region.minY, region.maxY): 

841 regionRow = Box2I(Point2I(region.minX, row), 

842 Extent2I(region.width, 1)) 

843 rowMeans.append(afwMath.makeStatistics(diffIm[regionRow], 

844 afwMath.MEAN, 

845 imStatsCtrl).getValue()) 

846 rowMeanVariance = afwMath.makeStatistics( 

847 np.array(rowMeans), afwMath.VARIANCECLIP, 

848 imStatsCtrl).getValue() 

849 

850 # Variance calculation via afwMath 

851 varDiff = 0.5*(afwMath.makeStatistics(diffIm, afwMath.VARIANCECLIP, imStatsCtrl).getValue()) 

852 

853 # Covariances calculations 

854 # Get the pixels that were not clipped 

855 varClip = afwMath.makeStatistics(diffIm, afwMath.VARIANCECLIP, imStatsCtrl).getValue() 

856 meanClip = afwMath.makeStatistics(diffIm, afwMath.MEANCLIP, imStatsCtrl).getValue() 

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

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

859 

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

861 # that were ignored by the clipping algorithm 

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

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

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

865 # calculations below. 

866 w = unmasked*wDiff 

867 

868 if np.sum(w) < self.config.minNumberGoodPixelsForCovariance/(self.config.binSize**2): 

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

870 "(than threshold %s)", np.sum(w), 

871 self.config.minNumberGoodPixelsForCovariance/(self.config.binSize**2)) 

872 return np.nan, np.nan, None, np.nan 

873 

874 maxRangeCov = self.config.maximumRangeCovariancesAstier 

875 

876 # Calculate covariances via FFT. 

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

878 # Calculate the sizes of FFT dimensions. 

879 s = shapeDiff + maxRangeCov 

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

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

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

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

884 # np.sum(w) is the same as npix[0][0] returned in covDiffAstier 

885 try: 

886 covDiffAstier = c.reportCovFastFourierTransform(maxRangeCov) 

887 except ValueError: 

888 # This is raised if there are not enough pixels. 

889 self.log.warning("Not enough pixels covering the requested covariance range in x/y (%d)", 

890 self.config.maximumRangeCovariancesAstier) 

891 return np.nan, np.nan, None, np.nan 

892 

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

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

895 # factor of 0.5 difference with afwMath.VARIANCECLIP. 

896 thresholdPercentage = self.config.thresholdDiffAfwVarVsCov00 

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

898 if fractionalDiff >= thresholdPercentage: 

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

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

901 

902 return mu, varDiff, covDiffAstier, rowMeanVariance 

903 

904 def getImageAreasMasksStats(self, exposure1, exposure2, region=None): 

905 """Get image areas in a region as well as masks and statistic objects. 

906 

907 Parameters 

908 ---------- 

909 exposure1 : `lsst.afw.image.ExposureF` 

910 First exposure of flat field pair. 

911 exposure2 : `lsst.afw.image.ExposureF` 

912 Second exposure of flat field pair. 

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

914 Region of each exposure where to perform the calculations 

915 (e.g, an amplifier). 

916 

917 Returns 

918 ------- 

919 im1Area : `lsst.afw.image.MaskedImageF` 

920 Masked image from exposure 1. 

921 im2Area : `lsst.afw.image.MaskedImageF` 

922 Masked image from exposure 2. 

923 imStatsCtrl : `lsst.afw.math.StatisticsControl` 

924 Statistics control object. 

925 mu1 : `float` 

926 Clipped mean of im1Area (ADU). 

927 mu2 : `float` 

928 Clipped mean of im2Area (ADU). 

929 """ 

930 if region is not None: 

931 im1Area = exposure1.maskedImage[region] 

932 im2Area = exposure2.maskedImage[region] 

933 else: 

934 im1Area = exposure1.maskedImage 

935 im2Area = exposure2.maskedImage 

936 

937 # Get mask planes and construct statistics control object from one 

938 # of the exposures 

939 imMaskVal = exposure1.getMask().getPlaneBitMask(self.config.maskNameList) 

940 imStatsCtrl = afwMath.StatisticsControl(self.config.nSigmaClipPtc, 

941 self.config.nIterSigmaClipPtc, 

942 imMaskVal) 

943 imStatsCtrl.setNanSafe(True) 

944 imStatsCtrl.setAndMask(imMaskVal) 

945 

946 mu1 = afwMath.makeStatistics(im1Area, afwMath.MEANCLIP, imStatsCtrl).getValue() 

947 mu2 = afwMath.makeStatistics(im2Area, afwMath.MEANCLIP, imStatsCtrl).getValue() 

948 

949 return (im1Area, im2Area, imStatsCtrl, mu1, mu2) 

950 

951 def getGainFromFlatPair(self, im1Area, im2Area, imStatsCtrl, mu1, mu2, 

952 correctionType='NONE', readNoise=None): 

953 """Estimate the gain from a single pair of flats. 

954 

955 The basic premise is 1/g = <(I1 - I2)^2/(I1 + I2)> = 1/const, 

956 where I1 and I2 correspond to flats 1 and 2, respectively. 

957 Corrections for the variable QE and the read-noise are then 

958 made following the derivation in Robert Lupton's forthcoming 

959 book, which gets 

960 

961 1/g = <(I1 - I2)^2/(I1 + I2)> - 1/mu(sigma^2 - 1/2g^2). 

962 

963 This is a quadratic equation, whose solutions are given by: 

964 

965 g = mu +/- sqrt(2*sigma^2 - 2*const*mu + mu^2)/(2*const*mu*2 

966 - 2*sigma^2) 

967 

968 where 'mu' is the average signal level and 'sigma' is the 

969 amplifier's readnoise. The positive solution will be used. 

970 The way the correction is applied depends on the value 

971 supplied for correctionType. 

972 

973 correctionType is one of ['NONE', 'SIMPLE' or 'FULL'] 

974 'NONE' : uses the 1/g = <(I1 - I2)^2/(I1 + I2)> formula. 

975 'SIMPLE' : uses the gain from the 'NONE' method for the 

976 1/2g^2 term. 

977 'FULL' : solves the full equation for g, discarding the 

978 non-physical solution to the resulting quadratic. 

979 

980 Parameters 

981 ---------- 

982 im1Area : `lsst.afw.image.maskedImage.MaskedImageF` 

983 Masked image from exposure 1. 

984 im2Area : `lsst.afw.image.maskedImage.MaskedImageF` 

985 Masked image from exposure 2. 

986 imStatsCtrl : `lsst.afw.math.StatisticsControl` 

987 Statistics control object. 

988 mu1: `float` 

989 Clipped mean of im1Area (ADU). 

990 mu2: `float` 

991 Clipped mean of im2Area (ADU). 

992 correctionType : `str`, optional 

993 The correction applied, one of ['NONE', 'SIMPLE', 'FULL'] 

994 readNoise : `float`, optional 

995 Amplifier readout noise (ADU). 

996 

997 Returns 

998 ------- 

999 gain : `float` 

1000 Gain, in e/ADU. 

1001 

1002 Raises 

1003 ------ 

1004 RuntimeError 

1005 Raise if `correctionType` is not one of 'NONE', 

1006 'SIMPLE', or 'FULL'. 

1007 """ 

1008 if correctionType not in ['NONE', 'SIMPLE', 'FULL']: 

1009 raise RuntimeError("Unknown correction type: %s" % correctionType) 

1010 

1011 if correctionType != 'NONE' and not np.isfinite(readNoise): 

1012 self.log.warning("'correctionType' in 'getGainFromFlatPair' is %s, " 

1013 "but 'readNoise' is NaN. Setting 'correctionType' " 

1014 "to 'NONE', so a gain value will be estimated without " 

1015 "corrections." % correctionType) 

1016 correctionType = 'NONE' 

1017 

1018 mu = 0.5*(mu1 + mu2) 

1019 

1020 # ratioIm = (I1 - I2)^2 / (I1 + I2) 

1021 temp = im2Area.clone() 

1022 ratioIm = im1Area.clone() 

1023 ratioIm -= temp 

1024 ratioIm *= ratioIm 

1025 

1026 # Sum of pairs 

1027 sumIm = im1Area.clone() 

1028 sumIm += temp 

1029 

1030 ratioIm /= sumIm 

1031 

1032 const = afwMath.makeStatistics(ratioIm, afwMath.MEAN, imStatsCtrl).getValue() 

1033 gain = 1. / const 

1034 

1035 if correctionType == 'SIMPLE': 

1036 gain = 1/(const - (1/mu)*(readNoise**2 - (1/2*gain**2))) 

1037 elif correctionType == 'FULL': 

1038 root = np.sqrt(mu**2 - 2*mu*const + 2*readNoise**2) 

1039 denom = (2*const*mu - 2*readNoise**2) 

1040 positiveSolution = (root + mu)/denom 

1041 gain = positiveSolution 

1042 

1043 return gain 

1044 

1045 def getReadNoise(self, exposureMetadata, taskMetadata, ampName): 

1046 """Gets readout noise for an amp from ISR metadata. 

1047 

1048 If possible, this attempts to get the now-standard headers 

1049 added to the exposure itself. If not found there, the ISR 

1050 TaskMetadata is searched. If neither of these has the value, 

1051 warn and set the read noise to NaN. 

1052 

1053 Parameters 

1054 ---------- 

1055 exposureMetadata : `lsst.daf.base.PropertySet` 

1056 Metadata to check for read noise first. 

1057 taskMetadata : `lsst.pipe.base.TaskMetadata` 

1058 List of exposures metadata from ISR for this exposure. 

1059 ampName : `str` 

1060 Amplifier name. 

1061 

1062 Returns 

1063 ------- 

1064 readNoise : `float` 

1065 The read noise for this set of exposure/amplifier. 

1066 """ 

1067 # Try from the exposure first. 

1068 expectedKey = f"LSST ISR OVERSCAN RESIDUAL SERIAL STDEV {ampName}" 

1069 if expectedKey in exposureMetadata: 

1070 return exposureMetadata[expectedKey] 

1071 

1072 # If not, try getting it from the task metadata. 

1073 expectedKey = f"RESIDUAL STDEV {ampName}" 

1074 if "isr" in taskMetadata: 

1075 if expectedKey in taskMetadata["isr"]: 

1076 return taskMetadata["isr"][expectedKey] 

1077 

1078 self.log.warning("Median readout noise from ISR metadata for amp %s " 

1079 "could not be calculated." % ampName) 

1080 return np.nan 

1081 

1082 def computeGaussianHistogramParameters(self, im1Area, im2Area, imStatsCtrl, mu1, mu2): 

1083 """Compute KS test for a Gaussian model fit to a histogram of the 

1084 difference image. 

1085 

1086 Parameters 

1087 ---------- 

1088 im1Area : `lsst.afw.image.MaskedImageF` 

1089 Masked image from exposure 1. 

1090 im2Area : `lsst.afw.image.MaskedImageF` 

1091 Masked image from exposure 2. 

1092 imStatsCtrl : `lsst.afw.math.StatisticsControl` 

1093 Statistics control object. 

1094 mu1 : `float` 

1095 Clipped mean of im1Area (ADU). 

1096 mu2 : `float` 

1097 Clipped mean of im2Area (ADU). 

1098 

1099 Returns 

1100 ------- 

1101 varFit : `float` 

1102 Variance from the Gaussian fit. 

1103 chi2Dof : `float` 

1104 Chi-squared per degree of freedom of Gaussian fit. 

1105 kspValue : `float` 

1106 The KS test p-value for the Gaussian fit. 

1107 

1108 Notes 

1109 ----- 

1110 The algorithm here was originally developed by Aaron Roodman. 

1111 Tests on the full focal plane of LSSTCam during testing has shown 

1112 that a KS test p-value cut of 0.01 is a good discriminant for 

1113 well-behaved flat pairs (p>0.01) and poorly behaved non-Gaussian 

1114 flat pairs (p<0.01). 

1115 """ 

1116 diffExp = im1Area.clone() 

1117 diffExp -= im2Area 

1118 

1119 sel = (((diffExp.mask.array & imStatsCtrl.getAndMask()) == 0) 

1120 & np.isfinite(diffExp.mask.array)) 

1121 diffArr = diffExp.image.array[sel] 

1122 

1123 numOk = len(diffArr) 

1124 

1125 if numOk >= self.config.ksHistMinDataValues and np.isfinite(mu1) and np.isfinite(mu2): 

1126 # Create a histogram symmetric around zero, with a bin size 

1127 # determined from the expected variance given by the average of 

1128 # the input signal levels. 

1129 lim = self.config.ksHistLimitMultiplier * np.sqrt((mu1 + mu2)/2.) 

1130 yVals, binEdges = np.histogram(diffArr, bins=self.config.ksHistNBins, range=[-lim, lim]) 

1131 

1132 # Fit the histogram with a Gaussian model. 

1133 model = GaussianModel() 

1134 yVals = yVals.astype(np.float64) 

1135 xVals = ((binEdges[0: -1] + binEdges[1:])/2.).astype(np.float64) 

1136 errVals = np.sqrt(yVals) 

1137 errVals[(errVals == 0.0)] = 1.0 

1138 pars = model.guess(yVals, x=xVals) 

1139 with warnings.catch_warnings(): 

1140 warnings.simplefilter("ignore") 

1141 # The least-squares fitter sometimes spouts (spurious) warnings 

1142 # when the model is very bad. Swallow these warnings now and 

1143 # let the KS test check the model below. 

1144 out = model.fit( 

1145 yVals, 

1146 pars, 

1147 x=xVals, 

1148 weights=1./errVals, 

1149 calc_covar=True, 

1150 method="least_squares", 

1151 ) 

1152 

1153 # Calculate chi2. 

1154 chiArr = out.residual 

1155 nDof = len(yVals) - 3 

1156 chi2Dof = np.sum(chiArr**2.)/nDof 

1157 sigmaFit = out.params["sigma"].value 

1158 

1159 # Calculate KS test p-value for the fit. 

1160 ksResult = scipy.stats.ks_1samp( 

1161 diffArr, 

1162 scipy.stats.norm.cdf, 

1163 (out.params["center"].value, sigmaFit), 

1164 ) 

1165 

1166 kspValue = ksResult.pvalue 

1167 if kspValue < 1e-15: 

1168 kspValue = 0.0 

1169 

1170 varFit = sigmaFit**2. 

1171 

1172 else: 

1173 varFit = np.nan 

1174 chi2Dof = np.nan 

1175 kspValue = 0.0 

1176 

1177 return varFit, chi2Dof, kspValue