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# 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# 

22 

23__all__ = ['MeasurePhotonTransferCurveTask', 

24 'MeasurePhotonTransferCurveTaskConfig', 

25 'PhotonTransferCurveDataset'] 

26 

27import numpy as np 

28import matplotlib.pyplot as plt 

29from sqlite3 import OperationalError 

30from collections import Counter 

31 

32import lsst.afw.math as afwMath 

33import lsst.pex.config as pexConfig 

34import lsst.pipe.base as pipeBase 

35from .utils import (NonexistentDatasetTaskDataIdContainer, PairedVisitListTaskRunner, 

36 checkExpLengthEqual, fitLeastSq, fitBootstrap, funcPolynomial, funcAstier) 

37from scipy.optimize import least_squares 

38 

39import datetime 

40 

41from .astierCovPtcUtils import (fftSize, CovFft, computeCovDirect, fitData) 

42from .linearity import LinearitySolveTask 

43 

44 

45class MeasurePhotonTransferCurveTaskConfig(pexConfig.Config): 

46 """Config class for photon transfer curve measurement task""" 

47 ccdKey = pexConfig.Field( 

48 dtype=str, 

49 doc="The key by which to pull a detector from a dataId, e.g. 'ccd' or 'detector'.", 

50 default='ccd', 

51 ) 

52 ptcFitType = pexConfig.ChoiceField( 

53 dtype=str, 

54 doc="Fit PTC to approximation in Astier+19 (Equation 16) or to a polynomial.", 

55 default="POLYNOMIAL", 

56 allowed={ 

57 "POLYNOMIAL": "n-degree polynomial (use 'polynomialFitDegree' to set 'n').", 

58 "EXPAPPROXIMATION": "Approximation in Astier+19 (Eq. 16).", 

59 "FULLCOVARIANCE": "Full covariances model in Astier+19 (Eq. 20)" 

60 } 

61 ) 

62 sigmaClipFullFitCovariancesAstier = pexConfig.Field( 

63 dtype=float, 

64 doc="sigma clip for full model fit for FULLCOVARIANCE ptcFitType ", 

65 default=5.0, 

66 ) 

67 maxIterFullFitCovariancesAstier = pexConfig.Field( 

68 dtype=int, 

69 doc="Maximum number of iterations in full model fit for FULLCOVARIANCE ptcFitType", 

70 default=3, 

71 ) 

72 maximumRangeCovariancesAstier = pexConfig.Field( 

73 dtype=int, 

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

75 default=8, 

76 ) 

77 covAstierRealSpace = pexConfig.Field( 

78 dtype=bool, 

79 doc="Calculate covariances in real space or via FFT? (see appendix A of Astier+19).", 

80 default=False, 

81 ) 

82 polynomialFitDegree = pexConfig.Field( 

83 dtype=int, 

84 doc="Degree of polynomial to fit the PTC, when 'ptcFitType'=POLYNOMIAL.", 

85 default=3, 

86 ) 

87 linearity = pexConfig.ConfigurableField( 

88 target=LinearitySolveTask, 

89 doc="Task to solve the linearity." 

90 ) 

91 

92 doCreateLinearizer = pexConfig.Field( 

93 dtype=bool, 

94 doc="Calculate non-linearity and persist linearizer?", 

95 default=False, 

96 ) 

97 

98 binSize = pexConfig.Field( 

99 dtype=int, 

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

101 default=1, 

102 ) 

103 minMeanSignal = pexConfig.Field( 

104 dtype=float, 

105 doc="Minimum value (inclusive) of mean signal (in DN) above which to consider.", 

106 default=0, 

107 ) 

108 maxMeanSignal = pexConfig.Field( 

109 dtype=float, 

110 doc="Maximum value (inclusive) of mean signal (in DN) below which to consider.", 

111 default=9e6, 

112 ) 

113 initialNonLinearityExclusionThresholdPositive = pexConfig.RangeField( 

114 dtype=float, 

115 doc="Initially exclude data points with a variance that are more than a factor of this from being" 

116 " linear in the positive direction, from the PTC fit. Note that these points will also be" 

117 " excluded from the non-linearity fit. This is done before the iterative outlier rejection," 

118 " to allow an accurate determination of the sigmas for said iterative fit.", 

119 default=0.12, 

120 min=0.0, 

121 max=1.0, 

122 ) 

123 initialNonLinearityExclusionThresholdNegative = pexConfig.RangeField( 

124 dtype=float, 

125 doc="Initially exclude data points with a variance that are more than a factor of this from being" 

126 " linear in the negative direction, from the PTC fit. Note that these points will also be" 

127 " excluded from the non-linearity fit. This is done before the iterative outlier rejection," 

128 " to allow an accurate determination of the sigmas for said iterative fit.", 

129 default=0.25, 

130 min=0.0, 

131 max=1.0, 

132 ) 

133 sigmaCutPtcOutliers = pexConfig.Field( 

134 dtype=float, 

135 doc="Sigma cut for outlier rejection in PTC.", 

136 default=5.0, 

137 ) 

138 maskNameList = pexConfig.ListField( 

139 dtype=str, 

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

141 default=['SUSPECT', 'BAD', 'NO_DATA'], 

142 ) 

143 nSigmaClipPtc = pexConfig.Field( 

144 dtype=float, 

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

146 default=5.5, 

147 ) 

148 nIterSigmaClipPtc = pexConfig.Field( 

149 dtype=int, 

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

151 default=1, 

152 ) 

153 maxIterationsPtcOutliers = pexConfig.Field( 

154 dtype=int, 

155 doc="Maximum number of iterations for outlier rejection in PTC.", 

156 default=2, 

157 ) 

158 doFitBootstrap = pexConfig.Field( 

159 dtype=bool, 

160 doc="Use bootstrap for the PTC fit parameters and errors?.", 

161 default=False, 

162 ) 

163 instrumentName = pexConfig.Field( 

164 dtype=str, 

165 doc="Instrument name.", 

166 default='', 

167 ) 

168 

169 

170class PhotonTransferCurveDataset: 

171 """A simple class to hold the output data from the PTC task. 

172 

173 The dataset is made up of a dictionary for each item, keyed by the 

174 amplifiers' names, which much be supplied at construction time. 

175 

176 New items cannot be added to the class to save accidentally saving to the 

177 wrong property, and the class can be frozen if desired. 

178 

179 inputVisitPairs records the visits used to produce the data. 

180 When fitPtc() or fitCovariancesAstier() is run, a mask is built up, which is by definition 

181 always the same length as inputVisitPairs, rawExpTimes, rawMeans 

182 and rawVars, and is a list of bools, which are incrementally set to False 

183 as points are discarded from the fits. 

184 

185 PTC fit parameters for polynomials are stored in a list in ascending order 

186 of polynomial term, i.e. par[0]*x^0 + par[1]*x + par[2]*x^2 etc 

187 with the length of the list corresponding to the order of the polynomial 

188 plus one. 

189 

190 Parameters 

191 ---------- 

192 ampNames : `list` 

193 List with the names of the amplifiers of the detector at hand. 

194 

195 ptcFitType : `str` 

196 Type of model fitted to the PTC: "POLYNOMIAL", "EXPAPPROXIMATION", or "FULLCOVARIANCE". 

197 

198 Returns 

199 ------- 

200 `lsst.cp.pipe.ptc.PhotonTransferCurveDataset` 

201 Output dataset from MeasurePhotonTransferCurveTask. 

202 """ 

203 

204 def __init__(self, ampNames, ptcFitType): 

205 # add items to __dict__ directly because __setattr__ is overridden 

206 

207 # instance variables 

208 self.__dict__["ptcFitType"] = ptcFitType 

209 self.__dict__["ampNames"] = ampNames 

210 self.__dict__["badAmps"] = [] 

211 

212 # raw data variables 

213 # visitMask is the mask produced after outlier rejection. The mask produced by "FULLCOVARIANCE" 

214 # may differ from the one produced in the other two PTC fit types. 

215 self.__dict__["inputVisitPairs"] = {ampName: [] for ampName in ampNames} 

216 self.__dict__["visitMask"] = {ampName: [] for ampName in ampNames} 

217 self.__dict__["rawExpTimes"] = {ampName: [] for ampName in ampNames} 

218 self.__dict__["rawMeans"] = {ampName: [] for ampName in ampNames} 

219 self.__dict__["rawVars"] = {ampName: [] for ampName in ampNames} 

220 

221 # Gain and noise 

222 self.__dict__["gain"] = {ampName: -1. for ampName in ampNames} 

223 self.__dict__["gainErr"] = {ampName: -1. for ampName in ampNames} 

224 self.__dict__["noise"] = {ampName: -1. for ampName in ampNames} 

225 self.__dict__["noiseErr"] = {ampName: -1. for ampName in ampNames} 

226 

227 # if ptcFitTye in ["POLYNOMIAL", "EXPAPPROXIMATION"] 

228 # fit information 

229 self.__dict__["ptcFitPars"] = {ampName: [] for ampName in ampNames} 

230 self.__dict__["ptcFitParsError"] = {ampName: [] for ampName in ampNames} 

231 self.__dict__["ptcFitReducedChiSquared"] = {ampName: [] for ampName in ampNames} 

232 

233 # if ptcFitTye in ["FULLCOVARIANCE"] 

234 # "covariancesTuple" is a numpy recarray with entries of the form 

235 # ['mu', 'i', 'j', 'var', 'cov', 'npix', 'ext', 'expTime', 'ampName'] 

236 # "covariancesFits" has CovFit objects that fit the measured covariances to Eq. 20 of Astier+19. 

237 # In "covariancesFitsWithNoB", "b"=0 in the model described by Eq. 20 of Astier+19. 

238 self.__dict__["covariancesTuple"] = {ampName: [] for ampName in ampNames} 

239 self.__dict__["covariancesFitsWithNoB"] = {ampName: [] for ampName in ampNames} 

240 self.__dict__["covariancesFits"] = {ampName: [] for ampName in ampNames} 

241 self.__dict__["aMatrix"] = {ampName: [] for ampName in ampNames} 

242 self.__dict__["bMatrix"] = {ampName: [] for ampName in ampNames} 

243 

244 # "final" means that the "raw" vectors above had "visitMask" applied. 

245 self.__dict__["finalVars"] = {ampName: [] for ampName in ampNames} 

246 self.__dict__["finalModelVars"] = {ampName: [] for ampName in ampNames} 

247 self.__dict__["finalMeans"] = {ampName: [] for ampName in ampNames} 

248 

249 def __setattr__(self, attribute, value): 

250 """Protect class attributes""" 

251 if attribute not in self.__dict__: 

252 raise AttributeError(f"{attribute} is not already a member of PhotonTransferCurveDataset, which" 

253 " does not support setting of new attributes.") 

254 else: 

255 self.__dict__[attribute] = value 

256 

257 def getVisitsUsed(self, ampName): 

258 """Get the visits used, i.e. not discarded, for a given amp. 

259 

260 If no mask has been created yet, all visits are returned. 

261 """ 

262 if len(self.visitMask[ampName]) == 0: 

263 return self.inputVisitPairs[ampName] 

264 

265 # if the mask exists it had better be the same length as the visitPairs 

266 assert len(self.visitMask[ampName]) == len(self.inputVisitPairs[ampName]) 

267 

268 pairs = self.inputVisitPairs[ampName] 

269 mask = self.visitMask[ampName] 

270 # cast to bool required because numpy 

271 return [(v1, v2) for ((v1, v2), m) in zip(pairs, mask) if bool(m) is True] 

272 

273 def getGoodAmps(self): 

274 return [amp for amp in self.ampNames if amp not in self.badAmps] 

275 

276 

277class MeasurePhotonTransferCurveTask(pipeBase.CmdLineTask): 

278 """A class to calculate, fit, and plot a PTC from a set of flat pairs. 

279 

280 The Photon Transfer Curve (var(signal) vs mean(signal)) is a standard tool 

281 used in astronomical detectors characterization (e.g., Janesick 2001, 

282 Janesick 2007). If ptcFitType is "EXPAPPROXIMATION" or "POLYNOMIAL", this task calculates the 

283 PTC from a series of pairs of flat-field images; each pair taken at identical exposure 

284 times. The difference image of each pair is formed to eliminate fixed pattern noise, 

285 and then the variance of the difference image and the mean of the average image 

286 are used to produce the PTC. An n-degree polynomial or the approximation in Equation 

287 16 of Astier+19 ("The Shape of the Photon Transfer Curve of CCD sensors", 

288 arXiv:1905.08677) can be fitted to the PTC curve. These models include 

289 parameters such as the gain (e/DN) and readout noise. 

290 

291 Linearizers to correct for signal-chain non-linearity are also calculated. 

292 The `Linearizer` class, in general, can support per-amp linearizers, but in this 

293 task this is not supported. 

294 

295 If ptcFitType is "FULLCOVARIANCE", the covariances of the difference images are calculated via the 

296 DFT methods described in Astier+19 and the variances for the PTC are given by the cov[0,0] elements 

297 at each signal level. The full model in Equation 20 of Astier+19 is fit to the PTC to get the gain 

298 and the noise. 

299 

300 Parameters 

301 ---------- 

302 

303 *args: `list` 

304 Positional arguments passed to the Task constructor. None used at this 

305 time. 

306 **kwargs: `dict` 

307 Keyword arguments passed on to the Task constructor. None used at this 

308 time. 

309 

310 """ 

311 

312 RunnerClass = PairedVisitListTaskRunner 

313 ConfigClass = MeasurePhotonTransferCurveTaskConfig 

314 _DefaultName = "measurePhotonTransferCurve" 

315 

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

317 pipeBase.CmdLineTask.__init__(self, *args, **kwargs) 

318 self.makeSubtask("linearity") 

319 plt.interactive(False) # stop windows popping up when plotting. When headless, use 'agg' backend too 

320 self.config.validate() 

321 self.config.freeze() 

322 

323 @classmethod 

324 def _makeArgumentParser(cls): 

325 """Augment argument parser for the MeasurePhotonTransferCurveTask.""" 

326 parser = pipeBase.ArgumentParser(name=cls._DefaultName) 

327 parser.add_argument("--visit-pairs", dest="visitPairs", nargs="*", 

328 help="Visit pairs to use. Each pair must be of the form INT,INT e.g. 123,456") 

329 parser.add_id_argument("--id", datasetType="photonTransferCurveDataset", 

330 ContainerClass=NonexistentDatasetTaskDataIdContainer, 

331 help="The ccds to use, e.g. --id ccd=0..100") 

332 return parser 

333 

334 @pipeBase.timeMethod 

335 def runDataRef(self, dataRef, visitPairs): 

336 """Run the Photon Transfer Curve (PTC) measurement task. 

337 

338 For a dataRef (which is each detector here), 

339 and given a list of visit pairs (postISR) at different exposure times, 

340 measure the PTC. 

341 

342 Parameters 

343 ---------- 

344 dataRef : list of lsst.daf.persistence.ButlerDataRef 

345 dataRef for the detector for the visits to be fit. 

346 

347 visitPairs : `iterable` of `tuple` of `int` 

348 Pairs of visit numbers to be processed together 

349 """ 

350 

351 # setup necessary objects 

352 detNum = dataRef.dataId[self.config.ccdKey] 

353 camera = dataRef.get('camera') 

354 detector = camera[dataRef.dataId[self.config.ccdKey]] 

355 # expand some missing fields that we need for lsstCam. This is a work-around 

356 # for Gen2 problems that I (RHL) don't feel like solving. The calibs pipelines 

357 # (which inherit from CalibTask) use addMissingKeys() to do basically the same thing 

358 # 

359 # Basically, the butler's trying to look up the fields in `raw_visit` which won't work 

360 for name in dataRef.getButler().getKeys('bias'): 

361 if name not in dataRef.dataId: 

362 try: 

363 dataRef.dataId[name] = \ 

364 dataRef.getButler().queryMetadata('raw', [name], detector=detNum)[0] 

365 except OperationalError: 

366 pass 

367 

368 amps = detector.getAmplifiers() 

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

370 datasetPtc = PhotonTransferCurveDataset(ampNames, self.config.ptcFitType) 

371 self.log.info('Measuring PTC using %s visits for detector %s' % (visitPairs, detector.getId())) 

372 

373 tupleRecords = [] 

374 allTags = [] 

375 for (v1, v2) in visitPairs: 

376 # Get postISR exposures. 

377 try: 

378 dataRef.dataId['expId'] = v1 

379 exp1 = dataRef.get("postISRCCD", immediate=True) 

380 dataRef.dataId['expId'] = v2 

381 exp2 = dataRef.get("postISRCCD", immediate=True) 

382 except RuntimeError: 

383 self.log.warn(f"postISR exposure for either expId {v1} or expId {v2} could not be retreived. " 

384 "Ignoring flat pair.") 

385 continue 

386 del dataRef.dataId['expId'] 

387 

388 checkExpLengthEqual(exp1, exp2, v1, v2, raiseWithMessage=True) 

389 expTime = exp1.getInfo().getVisitInfo().getExposureTime() 

390 tupleRows = [] 

391 nAmpsNan = 0 

392 for ampNumber, amp in enumerate(detector): 

393 ampName = amp.getName() 

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

395 doRealSpace = self.config.covAstierRealSpace 

396 muDiff, varDiff, covAstier = self.measureMeanVarCov(exp1, exp2, region=amp.getBBox(), 

397 covAstierRealSpace=doRealSpace) 

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

399 msg = (f"NaN mean or var, or None cov in amp {ampNumber} in visit pair {v1}, {v2} " 

400 "of detector {detNum}.") 

401 self.log.warn(msg) 

402 nAmpsNan += 1 

403 continue 

404 tags = ['mu', 'i', 'j', 'var', 'cov', 'npix', 'ext', 'expTime', 'ampName'] 

405 if (muDiff <= self.config.minMeanSignal) or (muDiff >= self.config.maxMeanSignal): 

406 continue 

407 datasetPtc.rawExpTimes[ampName].append(expTime) 

408 datasetPtc.rawMeans[ampName].append(muDiff) 

409 datasetPtc.rawVars[ampName].append(varDiff) 

410 datasetPtc.inputVisitPairs[ampName].append((v1, v2)) 

411 

412 tupleRows += [(muDiff, ) + covRow + (ampNumber, expTime, ampName) for covRow in covAstier] 

413 if nAmpsNan == len(ampNames): 

414 msg = f"NaN mean in all amps of visit pair {v1}, {v2} of detector {detNum}." 

415 self.log.warn(msg) 

416 continue 

417 allTags += tags 

418 tupleRecords += tupleRows 

419 covariancesWithTags = np.core.records.fromrecords(tupleRecords, names=allTags) 

420 

421 if self.config.ptcFitType in ["FULLCOVARIANCE", ]: 

422 # Calculate covariances and fit them, including the PTC, to Astier+19 full model (Eq. 20) 

423 datasetPtc = self.fitCovariancesAstier(datasetPtc, covariancesWithTags) 

424 elif self.config.ptcFitType in ["EXPAPPROXIMATION", "POLYNOMIAL"]: 

425 # Fit the PTC to a polynomial or to Astier+19 exponential approximation (Eq. 16) 

426 # Fill up PhotonTransferCurveDataset object. 

427 datasetPtc = self.fitPtc(datasetPtc, self.config.ptcFitType) 

428 

429 # Fit a poynomial to calculate non-linearity and persist linearizer. 

430 if self.config.doCreateLinearizer: 

431 # Fit (non)linearity of signal vs time curve. 

432 # Fill up PhotonTransferCurveDataset object. 

433 # Fill up array for LUT linearizer (tableArray). 

434 # Produce coefficients for Polynomial ans Squared linearizers. 

435 # Build linearizer objects. 

436 dimensions = {'camera': camera.getName(), 'detector': detector.getId()} 

437 linearityResults = self.linearity.run(datasetPtc, camera, dimensions) 

438 linearizer = linearityResults.outputLinearizer 

439 

440 butler = dataRef.getButler() 

441 self.log.info("Writing linearizer:") 

442 

443 detName = detector.getName() 

444 now = datetime.datetime.utcnow() 

445 calibDate = now.strftime("%Y-%m-%d") 

446 

447 butler.put(linearizer, datasetType='Linearizer', dataId={'detector': detNum, 

448 'detectorName': detName, 'calibDate': calibDate}) 

449 

450 self.log.info(f"Writing PTC data to {dataRef.getUri(write=True)}") 

451 dataRef.put(datasetPtc, datasetType="photonTransferCurveDataset") 

452 

453 return pipeBase.Struct(exitStatus=0) 

454 

455 def fitCovariancesAstier(self, dataset, covariancesWithTagsArray): 

456 """Fit measured flat covariances to full model in Astier+19. 

457 

458 Parameters 

459 ---------- 

460 dataset : `lsst.cp.pipe.ptc.PhotonTransferCurveDataset` 

461 The dataset containing information such as the means, variances and exposure times. 

462 

463 covariancesWithTagsArray : `numpy.recarray` 

464 Tuple with at least (mu, cov, var, i, j, npix), where: 

465 mu : 0.5*(m1 + m2), where: 

466 mu1: mean value of flat1 

467 mu2: mean value of flat2 

468 cov: covariance value at lag(i, j) 

469 var: variance(covariance value at lag(0, 0)) 

470 i: lag dimension 

471 j: lag dimension 

472 npix: number of pixels used for covariance calculation. 

473 

474 Returns 

475 ------- 

476 dataset: `lsst.cp.pipe.ptc.PhotonTransferCurveDataset` 

477 This is the same dataset as the input paramter, however, it has been modified 

478 to include information such as the fit vectors and the fit parameters. See 

479 the class `PhotonTransferCurveDatase`. 

480 """ 

481 

482 covFits, covFitsNoB = fitData(covariancesWithTagsArray, maxMu=self.config.maxMeanSignal, 

483 r=self.config.maximumRangeCovariancesAstier, 

484 nSigmaFullFit=self.config.sigmaClipFullFitCovariancesAstier, 

485 maxIterFullFit=self.config.maxIterFullFitCovariancesAstier) 

486 

487 dataset.covariancesTuple = covariancesWithTagsArray 

488 dataset.covariancesFits = covFits 

489 dataset.covariancesFitsWithNoB = covFitsNoB 

490 dataset = self.getOutputPtcDataCovAstier(dataset, covFits) 

491 

492 return dataset 

493 

494 def getOutputPtcDataCovAstier(self, dataset, covFits): 

495 """Get output data for PhotonTransferCurveCovAstierDataset from CovFit objects. 

496 

497 Parameters 

498 ---------- 

499 dataset : `lsst.cp.pipe.ptc.PhotonTransferCurveDataset` 

500 The dataset containing information such as the means, variances and exposure times. 

501 

502 covFits: `dict` 

503 Dictionary of CovFit objects, with amp names as keys. 

504 

505 Returns 

506 ------- 

507 dataset : `lsst.cp.pipe.ptc.PhotonTransferCurveDataset` 

508 This is the same dataset as the input paramter, however, it has been modified 

509 to include extra information such as the mask 1D array, gains, reoudout noise, measured signal, 

510 measured variance, modeled variance, a, and b coefficient matrices (see Astier+19) per amplifier. 

511 See the class `PhotonTransferCurveDatase`. 

512 """ 

513 

514 for i, amp in enumerate(covFits): 

515 fit = covFits[amp] 

516 (meanVecFinal, varVecFinal, varVecModel, 

517 wc, varMask) = fit.getFitData(0, 0, divideByMu=False, returnMasked=True) 

518 gain = fit.getGain() 

519 dataset.visitMask[amp] = varMask 

520 dataset.gain[amp] = gain 

521 dataset.gainErr[amp] = fit.getGainErr() 

522 dataset.noise[amp] = np.sqrt(np.fabs(fit.getRon())) 

523 dataset.noiseErr[amp] = fit.getRonErr() 

524 dataset.finalVars[amp].append(varVecFinal/(gain**2)) 

525 dataset.finalModelVars[amp].append(varVecModel/(gain**2)) 

526 dataset.finalMeans[amp].append(meanVecFinal/gain) 

527 dataset.aMatrix[amp].append(fit.getA()) 

528 dataset.bMatrix[amp].append(fit.getB()) 

529 

530 return dataset 

531 

532 def measureMeanVarCov(self, exposure1, exposure2, region=None, covAstierRealSpace=False): 

533 """Calculate the mean of each of two exposures and the variance and covariance of their difference. 

534 

535 The variance is calculated via afwMath, and the covariance via the methods in Astier+19 (appendix A). 

536 In theory, var = covariance[0,0]. This should be validated, and in the future, we may decide to just 

537 keep one (covariance). 

538 

539 Parameters 

540 ---------- 

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

542 First exposure of flat field pair. 

543 

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

545 Second exposure of flat field pair. 

546 

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

548 Region of each exposure where to perform the calculations (e.g, an amplifier). 

549 

550 covAstierRealSpace : `bool`, optional 

551 Should the covariannces in Astier+19 be calculated in real space or via FFT? 

552 See Appendix A of Astier+19. 

553 

554 Returns 

555 ------- 

556 mu : `float` or `NaN` 

557 0.5*(mu1 + mu2), where mu1, and mu2 are the clipped means of the regions in 

558 both exposures. If either mu1 or m2 are NaN's, the returned value is NaN. 

559 

560 varDiff : `float` or `NaN` 

561 Half of the clipped variance of the difference of the regions inthe two input 

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

563 

564 covDiffAstier : `list` or `NaN` 

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

566 dx : `int` 

567 Lag in x 

568 dy : `int` 

569 Lag in y 

570 var : `float` 

571 Variance at (dx, dy). 

572 cov : `float` 

573 Covariance at (dx, dy). 

574 nPix : `int` 

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

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

577 """ 

578 

579 if region is not None: 

580 im1Area = exposure1.maskedImage[region] 

581 im2Area = exposure2.maskedImage[region] 

582 else: 

583 im1Area = exposure1.maskedImage 

584 im2Area = exposure2.maskedImage 

585 

586 if self.config.binSize > 1: 

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

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

589 

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

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

592 self.config.nIterSigmaClipPtc, 

593 im1MaskVal) 

594 im1StatsCtrl.setNanSafe(True) 

595 im1StatsCtrl.setAndMask(im1MaskVal) 

596 

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

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

599 self.config.nIterSigmaClipPtc, 

600 im2MaskVal) 

601 im2StatsCtrl.setNanSafe(True) 

602 im2StatsCtrl.setAndMask(im2MaskVal) 

603 

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

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

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

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

608 return np.nan, np.nan, None 

609 mu = 0.5*(mu1 + mu2) 

610 

611 # Take difference of pairs 

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

613 temp = im2Area.clone() 

614 temp *= mu1 

615 diffIm = im1Area.clone() 

616 diffIm *= mu2 

617 diffIm -= temp 

618 diffIm /= mu 

619 

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

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

622 self.config.nIterSigmaClipPtc, 

623 diffImMaskVal) 

624 diffImStatsCtrl.setNanSafe(True) 

625 diffImStatsCtrl.setAndMask(diffImMaskVal) 

626 

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

628 

629 # Get the mask and identify good pixels as '1', and the rest as '0'. 

630 w1 = np.where(im1Area.getMask().getArray() == 0, 1, 0) 

631 w2 = np.where(im2Area.getMask().getArray() == 0, 1, 0) 

632 

633 w12 = w1*w2 

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

635 w = w12*wDiff 

636 

637 maxRangeCov = self.config.maximumRangeCovariancesAstier 

638 if covAstierRealSpace: 

639 covDiffAstier = computeCovDirect(diffIm.getImage().getArray(), w, maxRangeCov) 

640 else: 

641 shapeDiff = diffIm.getImage().getArray().shape 

642 fftShape = (fftSize(shapeDiff[0] + maxRangeCov), fftSize(shapeDiff[1]+maxRangeCov)) 

643 c = CovFft(diffIm.getImage().getArray(), w, fftShape, maxRangeCov) 

644 covDiffAstier = c.reportCovFft(maxRangeCov) 

645 

646 return mu, varDiff, covDiffAstier 

647 

648 def computeCovDirect(self, diffImage, weightImage, maxRange): 

649 """Compute covariances of diffImage in real space. 

650 

651 For lags larger than ~25, it is slower than the FFT way. 

652 Taken from https://github.com/PierreAstier/bfptc/ 

653 

654 Parameters 

655 ---------- 

656 diffImage : `numpy.array` 

657 Image to compute the covariance of. 

658 

659 weightImage : `numpy.array` 

660 Weight image of diffImage (1's and 0's for good and bad pixels, respectively). 

661 

662 maxRange : `int` 

663 Last index of the covariance to be computed. 

664 

665 Returns 

666 ------- 

667 outList : `list` 

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

669 dx : `int` 

670 Lag in x 

671 dy : `int` 

672 Lag in y 

673 var : `float` 

674 Variance at (dx, dy). 

675 cov : `float` 

676 Covariance at (dx, dy). 

677 nPix : `int` 

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

679 """ 

680 outList = [] 

681 var = 0 

682 # (dy,dx) = (0,0) has to be first 

683 for dy in range(maxRange + 1): 

684 for dx in range(0, maxRange + 1): 

685 if (dx*dy > 0): 

686 cov1, nPix1 = self.covDirectValue(diffImage, weightImage, dx, dy) 

687 cov2, nPix2 = self.covDirectValue(diffImage, weightImage, dx, -dy) 

688 cov = 0.5*(cov1 + cov2) 

689 nPix = nPix1 + nPix2 

690 else: 

691 cov, nPix = self.covDirectValue(diffImage, weightImage, dx, dy) 

692 if (dx == 0 and dy == 0): 

693 var = cov 

694 outList.append((dx, dy, var, cov, nPix)) 

695 

696 return outList 

697 

698 def covDirectValue(self, diffImage, weightImage, dx, dy): 

699 """Compute covariances of diffImage in real space at lag (dx, dy). 

700 

701 Taken from https://github.com/PierreAstier/bfptc/ (c.f., appendix of Astier+19). 

702 

703 Parameters 

704 ---------- 

705 diffImage : `numpy.array` 

706 Image to compute the covariance of. 

707 

708 weightImage : `numpy.array` 

709 Weight image of diffImage (1's and 0's for good and bad pixels, respectively). 

710 

711 dx : `int` 

712 Lag in x. 

713 

714 dy : `int` 

715 Lag in y. 

716 

717 Returns 

718 ------- 

719 cov : `float` 

720 Covariance at (dx, dy) 

721 

722 nPix : `int` 

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

724 """ 

725 (nCols, nRows) = diffImage.shape 

726 # switching both signs does not change anything: 

727 # it just swaps im1 and im2 below 

728 if (dx < 0): 

729 (dx, dy) = (-dx, -dy) 

730 # now, we have dx >0. We have to distinguish two cases 

731 # depending on the sign of dy 

732 if dy >= 0: 

733 im1 = diffImage[dy:, dx:] 

734 w1 = weightImage[dy:, dx:] 

735 im2 = diffImage[:nCols - dy, :nRows - dx] 

736 w2 = weightImage[:nCols - dy, :nRows - dx] 

737 else: 

738 im1 = diffImage[:nCols + dy, dx:] 

739 w1 = weightImage[:nCols + dy, dx:] 

740 im2 = diffImage[-dy:, :nRows - dx] 

741 w2 = weightImage[-dy:, :nRows - dx] 

742 # use the same mask for all 3 calculations 

743 wAll = w1*w2 

744 # do not use mean() because weightImage=0 pixels would then count 

745 nPix = wAll.sum() 

746 im1TimesW = im1*wAll 

747 s1 = im1TimesW.sum()/nPix 

748 s2 = (im2*wAll).sum()/nPix 

749 p = (im1TimesW*im2).sum()/nPix 

750 cov = p - s1*s2 

751 

752 return cov, nPix 

753 

754 @staticmethod 

755 def _initialParsForPolynomial(order): 

756 assert(order >= 2) 

757 pars = np.zeros(order, dtype=np.float) 

758 pars[0] = 10 

759 pars[1] = 1 

760 pars[2:] = 0.0001 

761 return pars 

762 

763 @staticmethod 

764 def _boundsForPolynomial(initialPars): 

765 lowers = [np.NINF for p in initialPars] 

766 uppers = [np.inf for p in initialPars] 

767 lowers[1] = 0 # no negative gains 

768 return (lowers, uppers) 

769 

770 @staticmethod 

771 def _boundsForAstier(initialPars): 

772 lowers = [np.NINF for p in initialPars] 

773 uppers = [np.inf for p in initialPars] 

774 return (lowers, uppers) 

775 

776 @staticmethod 

777 def _getInitialGoodPoints(means, variances, maxDeviationPositive, maxDeviationNegative): 

778 """Return a boolean array to mask bad points. 

779 

780 A linear function has a constant ratio, so find the median 

781 value of the ratios, and exclude the points that deviate 

782 from that by more than a factor of maxDeviationPositive/negative. 

783 Asymmetric deviations are supported as we expect the PTC to turn 

784 down as the flux increases, but sometimes it anomalously turns 

785 upwards just before turning over, which ruins the fits, so it 

786 is wise to be stricter about restricting positive outliers than 

787 negative ones. 

788 

789 Too high and points that are so bad that fit will fail will be included 

790 Too low and the non-linear points will be excluded, biasing the NL fit.""" 

791 ratios = [b/a for (a, b) in zip(means, variances)] 

792 medianRatio = np.median(ratios) 

793 ratioDeviations = [(r/medianRatio)-1 for r in ratios] 

794 

795 # so that it doesn't matter if the deviation is expressed as positive or negative 

796 maxDeviationPositive = abs(maxDeviationPositive) 

797 maxDeviationNegative = -1. * abs(maxDeviationNegative) 

798 

799 goodPoints = np.array([True if (r < maxDeviationPositive and r > maxDeviationNegative) 

800 else False for r in ratioDeviations]) 

801 return goodPoints 

802 

803 def _makeZeroSafe(self, array, warn=True, substituteValue=1e-9): 

804 """""" 

805 nBad = Counter(array)[0] 

806 if nBad == 0: 

807 return array 

808 

809 if warn: 

810 msg = f"Found {nBad} zeros in array at elements {[x for x in np.where(array==0)[0]]}" 

811 self.log.warn(msg) 

812 

813 array[array == 0] = substituteValue 

814 return array 

815 

816 def fitPtc(self, dataset, ptcFitType): 

817 """Fit the photon transfer curve to a polynimial or to Astier+19 approximation. 

818 

819 Fit the photon transfer curve with either a polynomial of the order 

820 specified in the task config, or using the Astier approximation. 

821 

822 Sigma clipping is performed iteratively for the fit, as well as an 

823 initial clipping of data points that are more than 

824 config.initialNonLinearityExclusionThreshold away from lying on a 

825 straight line. This other step is necessary because the photon transfer 

826 curve turns over catastrophically at very high flux (because saturation 

827 drops the variance to ~0) and these far outliers cause the initial fit 

828 to fail, meaning the sigma cannot be calculated to perform the 

829 sigma-clipping. 

830 

831 Parameters 

832 ---------- 

833 dataset : `lsst.cp.pipe.ptc.PhotonTransferCurveDataset` 

834 The dataset containing the means, variances and exposure times 

835 

836 ptcFitType : `str` 

837 Fit a 'POLYNOMIAL' (degree: 'polynomialFitDegree') or 

838 'EXPAPPROXIMATION' (Eq. 16 of Astier+19) to the PTC 

839 

840 Returns 

841 ------- 

842 dataset: `lsst.cp.pipe.ptc.PhotonTransferCurveDataset` 

843 This is the same dataset as the input paramter, however, it has been modified 

844 to include information such as the fit vectors and the fit parameters. See 

845 the class `PhotonTransferCurveDatase`. 

846 """ 

847 

848 def errFunc(p, x, y): 

849 return ptcFunc(p, x) - y 

850 

851 sigmaCutPtcOutliers = self.config.sigmaCutPtcOutliers 

852 maxIterationsPtcOutliers = self.config.maxIterationsPtcOutliers 

853 

854 for i, ampName in enumerate(dataset.ampNames): 

855 timeVecOriginal = np.array(dataset.rawExpTimes[ampName]) 

856 meanVecOriginal = np.array(dataset.rawMeans[ampName]) 

857 varVecOriginal = np.array(dataset.rawVars[ampName]) 

858 varVecOriginal = self._makeZeroSafe(varVecOriginal) 

859 

860 mask = ((meanVecOriginal >= self.config.minMeanSignal) & 

861 (meanVecOriginal <= self.config.maxMeanSignal)) 

862 

863 goodPoints = self._getInitialGoodPoints(meanVecOriginal, varVecOriginal, 

864 self.config.initialNonLinearityExclusionThresholdPositive, 

865 self.config.initialNonLinearityExclusionThresholdNegative) 

866 mask = mask & goodPoints 

867 

868 if ptcFitType == 'EXPAPPROXIMATION': 

869 ptcFunc = funcAstier 

870 parsIniPtc = [-1e-9, 1.0, 10.] # a00, gain, noise 

871 bounds = self._boundsForAstier(parsIniPtc) 

872 if ptcFitType == 'POLYNOMIAL': 

873 ptcFunc = funcPolynomial 

874 parsIniPtc = self._initialParsForPolynomial(self.config.polynomialFitDegree + 1) 

875 bounds = self._boundsForPolynomial(parsIniPtc) 

876 

877 # Before bootstrap fit, do an iterative fit to get rid of outliers 

878 count = 1 

879 while count <= maxIterationsPtcOutliers: 

880 # Note that application of the mask actually shrinks the array 

881 # to size rather than setting elements to zero (as we want) so 

882 # always update mask itself and re-apply to the original data 

883 meanTempVec = meanVecOriginal[mask] 

884 varTempVec = varVecOriginal[mask] 

885 res = least_squares(errFunc, parsIniPtc, bounds=bounds, args=(meanTempVec, varTempVec)) 

886 pars = res.x 

887 

888 # change this to the original from the temp because the masks are ANDed 

889 # meaning once a point is masked it's always masked, and the masks must 

890 # always be the same length for broadcasting 

891 sigResids = (varVecOriginal - ptcFunc(pars, meanVecOriginal))/np.sqrt(varVecOriginal) 

892 newMask = np.array([True if np.abs(r) < sigmaCutPtcOutliers else False for r in sigResids]) 

893 mask = mask & newMask 

894 

895 nDroppedTotal = Counter(mask)[False] 

896 self.log.debug(f"Iteration {count}: discarded {nDroppedTotal} points in total for {ampName}") 

897 count += 1 

898 # objects should never shrink 

899 assert (len(mask) == len(timeVecOriginal) == len(meanVecOriginal) == len(varVecOriginal)) 

900 

901 dataset.visitMask[ampName] = mask # store the final mask 

902 parsIniPtc = pars 

903 meanVecFinal = meanVecOriginal[mask] 

904 varVecFinal = varVecOriginal[mask] 

905 

906 if Counter(mask)[False] > 0: 

907 self.log.info((f"Number of points discarded in PTC of amplifier {ampName}:" + 

908 f" {Counter(mask)[False]} out of {len(meanVecOriginal)}")) 

909 

910 if (len(meanVecFinal) < len(parsIniPtc)): 

911 msg = (f"\nSERIOUS: Not enough data points ({len(meanVecFinal)}) compared to the number of" 

912 f"parameters of the PTC model({len(parsIniPtc)}). Setting {ampName} to BAD.") 

913 self.log.warn(msg) 

914 # The first and second parameters of initial fit are discarded (bias and gain) 

915 # for the final NL coefficients 

916 dataset.badAmps.append(ampName) 

917 dataset.gain[ampName] = np.nan 

918 dataset.gainErr[ampName] = np.nan 

919 dataset.noise[ampName] = np.nan 

920 dataset.noiseErr[ampName] = np.nan 

921 dataset.ptcFitPars[ampName] = np.nan 

922 dataset.ptcFitParsError[ampName] = np.nan 

923 dataset.ptcFitReducedChiSquared[ampName] = np.nan 

924 continue 

925 

926 # Fit the PTC 

927 if self.config.doFitBootstrap: 

928 parsFit, parsFitErr, reducedChiSqPtc = fitBootstrap(parsIniPtc, meanVecFinal, 

929 varVecFinal, ptcFunc, 

930 weightsY=1./np.sqrt(varVecFinal)) 

931 else: 

932 parsFit, parsFitErr, reducedChiSqPtc = fitLeastSq(parsIniPtc, meanVecFinal, 

933 varVecFinal, ptcFunc, 

934 weightsY=1./np.sqrt(varVecFinal)) 

935 dataset.ptcFitPars[ampName] = parsFit 

936 dataset.ptcFitParsError[ampName] = parsFitErr 

937 dataset.ptcFitReducedChiSquared[ampName] = reducedChiSqPtc 

938 

939 if ptcFitType == 'EXPAPPROXIMATION': 

940 ptcGain = parsFit[1] 

941 ptcGainErr = parsFitErr[1] 

942 ptcNoise = np.sqrt(np.fabs(parsFit[2])) 

943 ptcNoiseErr = 0.5*(parsFitErr[2]/np.fabs(parsFit[2]))*np.sqrt(np.fabs(parsFit[2])) 

944 if ptcFitType == 'POLYNOMIAL': 

945 ptcGain = 1./parsFit[1] 

946 ptcGainErr = np.fabs(1./parsFit[1])*(parsFitErr[1]/parsFit[1]) 

947 ptcNoise = np.sqrt(np.fabs(parsFit[0]))*ptcGain 

948 ptcNoiseErr = (0.5*(parsFitErr[0]/np.fabs(parsFit[0]))*(np.sqrt(np.fabs(parsFit[0]))))*ptcGain 

949 dataset.gain[ampName] = ptcGain 

950 dataset.gainErr[ampName] = ptcGainErr 

951 dataset.noise[ampName] = ptcNoise 

952 dataset.noiseErr[ampName] = ptcNoiseErr 

953 if not len(dataset.ptcFitType) == 0: 

954 dataset.ptcFitType = ptcFitType 

955 

956 return dataset