Coverage for python/lsst/cp/pipe/ptc/cpSolvePtcTask.py: 11%

493 statements  

« prev     ^ index     » next       coverage.py v7.4.1, created at 2024-01-30 12:35 +0000

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 collections import Counter 

24 

25import lsst.pex.config as pexConfig 

26import lsst.pipe.base as pipeBase 

27from lsst.cp.pipe.utils import (fitLeastSq, fitBootstrap, funcPolynomial, 

28 funcAstier, symmetrize, Pol2D) 

29 

30from scipy.signal import fftconvolve 

31from scipy.optimize import least_squares 

32from itertools import groupby 

33from operator import itemgetter 

34 

35import lsst.pipe.base.connectionTypes as cT 

36 

37from lsst.ip.isr import PhotonTransferCurveDataset 

38 

39import copy 

40 

41 

42__all__ = ['PhotonTransferCurveSolveConfig', 'PhotonTransferCurveSolveTask'] 

43 

44 

45class PhotonTransferCurveSolveConnections(pipeBase.PipelineTaskConnections, 

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

47 inputCovariances = cT.Input( 

48 name="ptcCovariances", 

49 doc="Tuple with measured covariances from flats.", 

50 storageClass="PhotonTransferCurveDataset", 

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

52 isCalibration=True, 

53 multiple=True, 

54 ) 

55 camera = cT.PrerequisiteInput( 

56 name="camera", 

57 doc="Camera the input data comes from.", 

58 storageClass="Camera", 

59 dimensions=("instrument",), 

60 isCalibration=True, 

61 ) 

62 outputPtcDataset = cT.Output( 

63 name="ptcDatsetProposal", 

64 doc="Output proposed ptc dataset.", 

65 storageClass="PhotonTransferCurveDataset", 

66 dimensions=("instrument", "detector"), 

67 multiple=False, 

68 isCalibration=True, 

69 ) 

70 

71 

72class PhotonTransferCurveSolveConfig(pipeBase.PipelineTaskConfig, 

73 pipelineConnections=PhotonTransferCurveSolveConnections): 

74 """Configuration for fitting measured covariances. 

75 """ 

76 

77 ptcFitType = pexConfig.ChoiceField( 

78 dtype=str, 

79 doc="Fit PTC to Eq. 16, Eq. 20 in Astier+19, or to a polynomial.", 

80 default="POLYNOMIAL", 

81 allowed={ 

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

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

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

85 } 

86 ) 

87 minMeanSignal = pexConfig.DictField( 

88 keytype=str, 

89 itemtype=float, 

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

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

92 " {'ALL_AMPS': value}", 

93 default={'ALL_AMPS': 0.0}, 

94 ) 

95 maxMeanSignal = pexConfig.DictField( 

96 keytype=str, 

97 itemtype=float, 

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

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

100 " {'ALL_AMPS': value}", 

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

102 ) 

103 maximumRangeCovariancesAstier = pexConfig.Field( 

104 dtype=int, 

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

106 default=8, 

107 ) 

108 maximumRangeCovariancesAstierFullCovFit = pexConfig.Field( 

109 dtype=int, 

110 doc="Maximum range up to where to fit covariances as in Astier+19, " 

111 "for the FULLCOVARIANCE model." 

112 "This is different from maximumRangeCovariancesAstier." 

113 "It should be less or equal than maximumRangeCovariancesAstier." 

114 "The number of parameters for this model is " 

115 "3*maximumRangeCovariancesAstierFullCovFit^2 + 1, so increase with care " 

116 "so that the fit is not too slow.", 

117 default=8, 

118 ) 

119 doSubtractLongRangeCovariances = pexConfig.Field( 

120 dtype=bool, 

121 doc="Subtract long-range covariances before FULLCOVARIANCE fit, " 

122 "beyond startLongRangeCovariances?", 

123 default=False, 

124 ) 

125 startLongRangeCovariances = pexConfig.Field( 

126 dtype=int, 

127 doc="If doSubtractLongRangeCovariances is True, subtract covariances " 

128 "beyond this range. It should be less than maximumRangeCovariancesAstier. ", 

129 default=4, 

130 ) 

131 polyDegLongRangeCovariances = pexConfig.Field( 

132 dtype=int, 

133 doc="If doSubtractLongRangeCovariances is True, polynomial " 

134 "degree to fit data beyond startLongRangeCovariances.", 

135 default=1, 

136 ) 

137 sigmaClipFullFitCovariancesAstier = pexConfig.Field( 

138 dtype=float, 

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

140 default=5.0, 

141 ) 

142 maxIterFullFitCovariancesAstier = pexConfig.Field( 

143 dtype=int, 

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

145 default=3, 

146 ) 

147 polynomialFitDegree = pexConfig.Field( 

148 dtype=int, 

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

150 default=3, 

151 ) 

152 doLegacyTurnoffSelection = pexConfig.Field( 

153 dtype=bool, 

154 doc="Use 'legacy' computation for PTC turnoff selection. If set " 

155 "to False, then the KS test p-value selection will be used instead.", 

156 default=False, 

157 ) 

158 sigmaCutPtcOutliers = pexConfig.Field( 

159 dtype=float, 

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

161 default=5.0, 

162 ) 

163 maxIterationsPtcOutliers = pexConfig.RangeField( 

164 dtype=int, 

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

166 default=2, 

167 min=0 

168 ) 

169 maxSignalInitialPtcOutlierFit = pexConfig.Field( 

170 dtype=float, 

171 doc="Maximum signal considered for intial outlier fit. This should be below " 

172 "the PTC turnoff to ensure accurate outlier rejection. If " 

173 "scaleMaxSignalInitialPtcOutlierFit=True then the units are electrons; " 

174 "otherwise ADU.", 

175 default=50_000., 

176 ) 

177 scaleMaxSignalInitialPtcOutlierFit = pexConfig.Field( 

178 dtype=bool, 

179 doc="Scale maxSignalInitialPtcOutlierFit by approximate gain? If yes then " 

180 "maxSignalInitialPtcOutlierFit is assumed to have units of electrons, " 

181 "otherwise ADU.", 

182 default=True, 

183 ) 

184 minVarPivotSearch = pexConfig.Field( 

185 dtype=float, 

186 doc="The code looks for a pivot signal point after which the variance starts decreasing at high-flux" 

187 " to exclude then from the PTC model fit. However, sometimes at low fluxes, the variance" 

188 " decreases slightly. Set this variable for the variance value, in ADU^2, after which the pivot " 

189 " should be sought. Only used if doLegacyTurnoffSelection is True.", 

190 default=10000, 

191 ) 

192 consecutivePointsVarDecreases = pexConfig.RangeField( 

193 dtype=int, 

194 doc="Required number of consecutive points/fluxes in the PTC where the variance " 

195 "decreases in order to find a first estimate of the PTC turn-off. " 

196 "Only used if doLegacyTurnoffSelection is True.", 

197 default=2, 

198 min=2 

199 ) 

200 ksTestMinPvalue = pexConfig.Field( 

201 dtype=float, 

202 doc="Minimum value of the Gaussian histogram KS test p-value to be used in PTC fit. " 

203 "Only used if doLegacyTurnoffSelection is False.", 

204 default=0.01, 

205 ) 

206 doFitBootstrap = pexConfig.Field( 

207 dtype=bool, 

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

209 default=False, 

210 ) 

211 binSize = pexConfig.Field( 

212 dtype=int, 

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

214 default=1, 

215 ) 

216 

217 def validate(self): 

218 super().validate() 

219 fitMatrixSide = self.maximumRangeCovariancesAstierFullCovFit 

220 measureMatrixSide = self.maximumRangeCovariancesAstier 

221 if self.ptcFitType == "FULLCOVARIANCE": 

222 if fitMatrixSide > measureMatrixSide: 

223 raise RuntimeError("Covariance fit size %s is larger than" 

224 "measurement size %s.", 

225 fitMatrixSide, measureMatrixSide) 

226 if self.doSubtractLongRangeCovariances: 

227 startLag = self.startLongRangeCovariances 

228 if measureMatrixSide < startLag: 

229 raise RuntimeError("Covariance measure size %s is smaller than long" 

230 "-range covariance starting point %s.", 

231 measureMatrixSide, startLag) 

232 

233 

234class PhotonTransferCurveSolveTask(pipeBase.PipelineTask): 

235 """Task to fit the PTC from flat covariances. 

236 

237 The first task of the PTC measurement pipeline, 

238 ``PhotonTransferCurveMeasureTask`` (and assumed to have been run 

239 before this task), produced a list of 

240 `~lsst.ip.isr.PhotonTransferCurveDataset` objects. Each dataset 

241 contains the mean signal and covariances of the 

242 difference image of the flat-field images taken at 

243 the same exposure time. The list also contains dummy 

244 datasets (with no measurements), whose purpose is to have 

245 the input and output dimensions of ``PhotonTransferCurveMeasureTask`` 

246 match. 

247 

248 This task, ``PhotonTransferCurveSolveTask``, assembles the list 

249 of individual PTC datasets produced 

250 by ``PhotonTransferCurveMeasureTask`` into one single final PTC 

251 dataset, discarding the dummy datset as appropiate. 

252 The task fits the measured (co)variances to one of three models: 

253 a polynomial model of a given order, or the models described 

254 in equations 16 and 20 of Astier+19. These options are referred 

255 to as ``POLYNOMIAL``, ``EXPAPPROXIMATION``, and ``FULLCOVARIANCE`` 

256 in the configuration options of the task, respectively). 

257 Parameters of interest such as the gain and noise are derived 

258 from the fits. The ``FULLCOVARIANCE`` model is fitted to the 

259 full covariance data (as oppossed to the other two models, which 

260 are fit to the variance vs mean measurements only). 

261 

262 Astier+19: "The Shape of the Photon Transfer Curve 

263 of CCD sensors", arXiv:1905.08677 

264 """ 

265 

266 ConfigClass = PhotonTransferCurveSolveConfig 

267 _DefaultName = 'cpPhotonTransferCurveSolve' 

268 

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

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

271 

272 Parameters 

273 ---------- 

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

275 Butler to operate on. 

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

277 Input data refs to load. 

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

279 Output data refs to persist. 

280 """ 

281 inputs = butlerQC.get(inputRefs) 

282 detId = inputRefs.inputCovariances[0].dataId['detector'] 

283 outputs = self.run(inputCovariances=inputs['inputCovariances'], camera=inputs['camera'], detId=detId) 

284 butlerQC.put(outputs, outputRefs) 

285 

286 def run(self, inputCovariances, camera=None, detId=0): 

287 """Fit measured covariances to different models. 

288 

289 Parameters 

290 ---------- 

291 inputCovariances : `list` [`lsst.ip.isr.PhotonTransferCurveDataset`] 

292 List of lsst.ip.isr.PhotonTransferCurveDataset datasets. 

293 camera : `lsst.afw.cameraGeom.Camera`, optional 

294 Input camera. 

295 detId : `int` 

296 Detector ID to locate the detector in the camera and 

297 populate the `lsst.ip.isr.PhotonTransferCurveDataset` 

298 metadata. 

299 Returns 

300 ------- 

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

302 The resultins structure contains: 

303 

304 ``outputPtcDatset`` 

305 Final PTC dataset, containing information such as the 

306 means, variances, and exposure times 

307 (`lsst.ip.isr.PhotonTransferCurveDataset`). 

308 """ 

309 # Find the ampNames from a non-dummy ptc. 

310 ampNames = [] 

311 for partialPtcDataset in inputCovariances: 

312 if partialPtcDataset.ptcFitType != 'DUMMY': 

313 ampNames = partialPtcDataset.ampNames 

314 break 

315 

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

317 # specified in the config. 

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

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

320 for ampName in ampNames: 

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

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

323 elif ampName in self.config.maxMeanSignal: 

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

325 

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

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

328 elif ampName in self.config.minMeanSignal: 

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

330 

331 # Assemble individual PTC datasets into a single PTC dataset. 

332 datasetPtc = PhotonTransferCurveDataset( 

333 ampNames=ampNames, 

334 ptcFitType=self.config.ptcFitType, 

335 covMatrixSide=self.config.maximumRangeCovariancesAstier, 

336 covMatrixSideFullCovFit=self.config.maximumRangeCovariancesAstierFullCovFit) 

337 for partialPtcDataset in inputCovariances: 

338 # Ignore dummy datasets 

339 if partialPtcDataset.ptcFitType == 'DUMMY': 

340 continue 

341 for ampName in ampNames: 

342 # The partial dataset consists of lists of values for each 

343 # quantity. In the case of the input exposure pairs, this is a 

344 # list of tuples. In all cases we only want the first 

345 # (and only) element of the list. 

346 datasetPtc.inputExpIdPairs[ampName].append(partialPtcDataset.inputExpIdPairs[ampName][0]) 

347 datasetPtc.rawExpTimes[ampName] = np.append(datasetPtc.rawExpTimes[ampName], 

348 partialPtcDataset.rawExpTimes[ampName][0]) 

349 datasetPtc.rawMeans[ampName] = np.append(datasetPtc.rawMeans[ampName], 

350 partialPtcDataset.rawMeans[ampName][0]) 

351 datasetPtc.rawVars[ampName] = np.append(datasetPtc.rawVars[ampName], 

352 partialPtcDataset.rawVars[ampName][0]) 

353 datasetPtc.rowMeanVariance[ampName] = np.append(datasetPtc.rowMeanVariance[ampName], 

354 partialPtcDataset.rowMeanVariance[ampName][0]) 

355 datasetPtc.photoCharges[ampName] = np.append(datasetPtc.photoCharges[ampName], 

356 partialPtcDataset.photoCharges[ampName][0]) 

357 datasetPtc.histVars[ampName] = np.append(datasetPtc.histVars[ampName], 

358 partialPtcDataset.histVars[ampName][0]) 

359 datasetPtc.histChi2Dofs[ampName] = np.append(datasetPtc.histChi2Dofs[ampName], 

360 partialPtcDataset.histChi2Dofs[ampName][0]) 

361 datasetPtc.kspValues[ampName] = np.append(datasetPtc.kspValues[ampName], 

362 partialPtcDataset.kspValues[ampName][0]) 

363 datasetPtc.covariances[ampName] = np.append( 

364 datasetPtc.covariances[ampName].ravel(), 

365 partialPtcDataset.covariances[ampName].ravel() 

366 ).reshape( 

367 ( 

368 len(datasetPtc.rawExpTimes[ampName]), 

369 datasetPtc.covMatrixSide, 

370 datasetPtc.covMatrixSide, 

371 ) 

372 ) 

373 datasetPtc.covariancesSqrtWeights[ampName] = np.append( 

374 datasetPtc.covariancesSqrtWeights[ampName].ravel(), 

375 partialPtcDataset.covariancesSqrtWeights[ampName].ravel() 

376 ).reshape( 

377 ( 

378 len(datasetPtc.rawExpTimes[ampName]), 

379 datasetPtc.covMatrixSide, 

380 datasetPtc.covMatrixSide, 

381 ) 

382 ) 

383 

384 # Apply min/max masking. 

385 rawMean = partialPtcDataset.rawMeans[ampName][0] 

386 rawVar = partialPtcDataset.rawVars[ampName][0] 

387 expIdMask = partialPtcDataset.expIdMask[ampName][0] 

388 if (rawMean <= minMeanSignalDict[ampName]) or (rawMean >= maxMeanSignalDict[ampName]) \ 

389 or not np.isfinite(rawMean) or not np.isfinite(rawVar): 

390 expIdMask = False 

391 

392 kspValue = partialPtcDataset.kspValues[ampName][0] 

393 if not self.config.doLegacyTurnoffSelection and \ 

394 kspValue < self.config.ksTestMinPvalue: 

395 expIdMask = False 

396 

397 datasetPtc.expIdMask[ampName] = np.append(datasetPtc.expIdMask[ampName], expIdMask) 

398 

399 for key, value in partialPtcDataset.auxValues.items(): 

400 if key in datasetPtc.auxValues: 

401 datasetPtc.auxValues[key] = np.append(datasetPtc.auxValues[key], value) 

402 else: 

403 datasetPtc.auxValues[key] = value 

404 

405 # Sort arrays that are filled so far in the final dataset by 

406 # rawMeans index. 

407 # First compute the mean across all the amps to make sure that they are 

408 # all sorted the same way. 

409 detectorMeans = np.zeros(len(datasetPtc.inputExpIdPairs[ampNames[0]])) 

410 

411 for i in range(len(detectorMeans)): 

412 arr = np.array([datasetPtc.rawMeans[ampName][i] for ampName in ampNames]) 

413 good, = (np.isfinite(arr)).nonzero() 

414 if good.size == 0: 

415 detectorMeans[i] = np.nan 

416 else: 

417 detectorMeans[i] = np.mean(arr[good]) 

418 

419 index = np.argsort(detectorMeans) 

420 

421 for ampName in ampNames: 

422 datasetPtc.inputExpIdPairs[ampName] = np.array( 

423 datasetPtc.inputExpIdPairs[ampName] 

424 )[index].tolist() 

425 datasetPtc.rawExpTimes[ampName] = datasetPtc.rawExpTimes[ampName][index] 

426 datasetPtc.rawMeans[ampName] = datasetPtc.rawMeans[ampName][index] 

427 datasetPtc.rawVars[ampName] = datasetPtc.rawVars[ampName][index] 

428 datasetPtc.rowMeanVariance[ampName] = datasetPtc.rowMeanVariance[ampName][index] 

429 datasetPtc.photoCharges[ampName] = datasetPtc.photoCharges[ampName][index] 

430 datasetPtc.histVars[ampName] = datasetPtc.histVars[ampName][index] 

431 datasetPtc.histChi2Dofs[ampName] = datasetPtc.histChi2Dofs[ampName][index] 

432 datasetPtc.kspValues[ampName] = datasetPtc.kspValues[ampName][index] 

433 datasetPtc.expIdMask[ampName] = datasetPtc.expIdMask[ampName][index] 

434 datasetPtc.covariances[ampName] = datasetPtc.covariances[ampName][index] 

435 datasetPtc.covariancesSqrtWeights[ampName] = datasetPtc.covariancesSqrtWeights[ampName][index] 

436 for key, value in datasetPtc.auxValues.items(): 

437 datasetPtc.auxValues[key] = value[index] 

438 

439 if self.config.ptcFitType == "FULLCOVARIANCE": 

440 # Fit the measured covariances vs mean signal to 

441 # the Astier+19 full model (Eq. 20). Before that 

442 # do a preliminary fit to the variance (C_00) vs mean 

443 # signal (mu) curve using the EXPAPPROXIMATION model 

444 # (Eq. 16 in Astier+19) in order to 

445 # get the flat pairs that are masked. The 

446 # points at these fluxes will also be masked when 

447 # calculating the other elements of the covariance 

448 # matrix, C_ij, i!=j). 

449 

450 # Preliminary fit, usign a temp dataset to get the mask 

451 tempDatasetPtc = copy.copy(datasetPtc) 

452 tempDatasetPtc.ptcFitType = "EXPAPPROXIMATION" 

453 tempDatasetPtc = self.fitMeasurementsToModel(tempDatasetPtc) 

454 

455 # "FULLCOVARIANCE", using the mask obtained from the 

456 # previous fit. 

457 for ampName in datasetPtc.ampNames: 

458 datasetPtc.expIdMask[ampName] = tempDatasetPtc.expIdMask[ampName] 

459 datasetPtc.fitType = "FULLCOVARIANCE" 

460 datasetPtc = self.fitMeasurementsToModel(datasetPtc) 

461 # The other options are: self.config.ptcFitType in 

462 # ("EXPAPPROXIMATION", "POLYNOMIAL") 

463 else: 

464 # Fit the PTC to a polynomial or to Astier+19 exponential 

465 # approximation (Eq. 16). Fill up 

466 # PhotonTransferCurveDataset object. 

467 datasetPtc = self.fitMeasurementsToModel(datasetPtc) 

468 

469 if camera: 

470 detector = camera[detId] 

471 else: 

472 detector = None 

473 datasetPtc.updateMetadataFromExposures(inputCovariances) 

474 datasetPtc.updateMetadata(setDate=True, camera=camera, detector=detector) 

475 

476 return pipeBase.Struct( 

477 outputPtcDataset=datasetPtc, 

478 ) 

479 

480 def fitMeasurementsToModel(self, dataset): 

481 """Fit the measured covariances vs mean signal to a 

482 polynomial or one of the models in Astier+19 

483 (Eq. 16 or Eq.20). 

484 

485 Parameters 

486 ---------- 

487 dataset : `lsst.ip.isr.ptcDataset.PhotonTransferCurveDataset` 

488 The dataset containing information such as the means, 

489 (co)variances, and exposure times. 

490 

491 Returns 

492 ------- 

493 dataset : `lsst.ip.isr.ptcDataset.PhotonTransferCurveDataset` 

494 This is the same dataset as the input parameter, however, 

495 it has been modified to include information such as the 

496 fit vectors and the fit parameters. See the class 

497 `PhotonTransferCurveDatase`. 

498 """ 

499 fitType = dataset.ptcFitType 

500 if fitType in ["FULLCOVARIANCE", ]: 

501 # This model uses the full covariance matrix in the fit. 

502 # The PTC is technically defined as variance vs signal, 

503 # with variance = Cov_00 

504 dataset = self.fitDataFullCovariance(dataset) 

505 elif fitType in ["POLYNOMIAL", "EXPAPPROXIMATION"]: 

506 # The PTC is technically defined as variance vs signal 

507 dataset = self.fitPtc(dataset) 

508 else: 

509 raise RuntimeError( 

510 f"Fitting option {fitType} not one of " 

511 "'POLYNOMIAL', 'EXPAPPROXIMATION', or 'FULLCOVARIANCE'" 

512 ) 

513 

514 return dataset 

515 

516 def fitDataFullCovariance(self, dataset): 

517 """Fit measured flat covariances to the full model in 

518 Astier+19 (Eq. 20). 

519 

520 Parameters 

521 ---------- 

522 dataset : `lsst.ip.isr.ptcDataset.PhotonTransferCurveDataset` 

523 The dataset containing information such as the means, 

524 (co)variances, and exposure times. 

525 

526 Returns 

527 ------- 

528 dataset : `lsst.ip.isr.ptcDataset.PhotonTransferCurveDataset` 

529 This is the same dataset as the input parameter, however, 

530 it has been modified to include information such as the 

531 fit vectors and the fit parameters. See the class 

532 `PhotonTransferCurveDatase`. 

533 

534 Notes 

535 ----- 

536 The parameters of the full model for C_ij(mu) ("C_ij" and "mu" 

537 in ADU^2 and ADU, respectively) in Astier+19 (Eq. 20) are: 

538 

539 - "a" coefficients (r by r matrix), units: 1/e 

540 - "b" coefficients (r by r matrix), units: 1/e 

541 - noise matrix (r by r matrix), units: e^2 

542 - gain, units: e/ADU 

543 

544 "b" appears in Eq. 20 only through the "ab" combination, which 

545 is defined in this code as "c=ab". 

546 

547 Total number of parameters: #entries(a) + #entries(c) + #entries(noise) 

548 + 1. This is equivalent to r^2 + r^2 + r^2 + 1, where "r" is the 

549 maximum lag considered for the covariances calculation, and the 

550 extra "1" is the gain. If "b" is 0, then "c" is 0, and len(pInit) will 

551 have r^2 fewer entries. 

552 """ 

553 matrixSide = dataset.covMatrixSide 

554 matrixSideFit = dataset.covMatrixSideFullCovFit 

555 lenParams = matrixSideFit*matrixSideFit 

556 

557 for ampName in dataset.ampNames: 

558 lenInputTimes = len(dataset.rawExpTimes[ampName]) 

559 # Not used when ptcFitType is 'FULLCOVARIANCE' 

560 dataset.ptcFitPars[ampName] = np.array([np.nan]) 

561 dataset.ptcFitParsError[ampName] = np.array([np.nan]) 

562 dataset.ptcFitChiSq[ampName] = np.nan 

563 

564 if ampName in dataset.badAmps: 

565 # Bad amp 

566 # Entries need to have proper dimensions so read/write 

567 # with astropy.Table works. 

568 nanMatrixFit = np.full((matrixSideFit, matrixSideFit), np.nan) 

569 listNanMatrix = np.full((lenInputTimes, matrixSide, matrixSide), np.nan) 

570 listNanMatrixFit = np.full((lenInputTimes, matrixSideFit, matrixSideFit), np.nan) 

571 dataset.covariancesModel[ampName] = listNanMatrixFit 

572 dataset.covariancesSqrtWeights[ampName] = listNanMatrix 

573 dataset.aMatrix[ampName] = nanMatrixFit 

574 dataset.bMatrix[ampName] = nanMatrixFit 

575 dataset.covariancesModelNoB[ampName] = listNanMatrixFit 

576 dataset.aMatrixNoB[ampName] = nanMatrixFit 

577 dataset.noiseMatrix[ampName] = nanMatrixFit 

578 dataset.noiseMatrixNoB[ampName] = nanMatrixFit 

579 

580 dataset.expIdMask[ampName] = np.repeat(False, lenInputTimes) 

581 dataset.gain[ampName] = np.nan 

582 dataset.gainErr[ampName] = np.nan 

583 dataset.noise[ampName] = np.nan 

584 dataset.noiseErr[ampName] = np.nan 

585 dataset.finalVars[ampName] = np.repeat(np.nan, lenInputTimes) 

586 dataset.finalModelVars[ampName] = np.repeat(np.nan, lenInputTimes) 

587 dataset.finalMeans[ampName] = np.repeat(np.nan, lenInputTimes) 

588 continue 

589 

590 muAtAmp = dataset.rawMeans[ampName] 

591 maskAtAmp = dataset.expIdMask[ampName] 

592 if len(maskAtAmp) == 0: 

593 maskAtAmp = np.repeat(True, len(muAtAmp)) 

594 

595 muAtAmpMasked = muAtAmp[maskAtAmp] 

596 covAtAmp = dataset.covariances[ampName] 

597 covAtAmpMasked = np.nan_to_num(covAtAmp)[maskAtAmp] 

598 covSqrtWeightsAtAmp = dataset.covariancesSqrtWeights[ampName] 

599 covSqrtWeightsAtAmpMasked = np.nan_to_num(covSqrtWeightsAtAmp)[maskAtAmp] 

600 

601 # Subtract long-range covariances 

602 if self.config.doSubtractLongRangeCovariances: 

603 startLag = self.config.startLongRangeCovariances 

604 covAtAmpMasked, covSqrtWeightsAtAmpMasked = self.subtractDistantOffset( 

605 muAtAmpMasked, covAtAmpMasked, 

606 covSqrtWeightsAtAmpMasked, 

607 start=startLag, 

608 degree=self.config.polyDegLongRangeCovariances) 

609 

610 # In principle, we could fit to a lag smaller than the measured 

611 # covariances. 

612 r = self.config.maximumRangeCovariancesAstierFullCovFit 

613 covAtAmpForFitMasked = covAtAmpMasked[:, :r, :r] 

614 covSqrtWeightsAtAmpForFitMasked = covSqrtWeightsAtAmpMasked[:, :r, :r] 

615 

616 # Initial fit, to approximate parameters, with c=0 

617 a0, c0, noise0, gain0 = self.initialFitFullCovariance( 

618 muAtAmpMasked, 

619 covAtAmpForFitMasked, 

620 covSqrtWeightsAtAmpForFitMasked 

621 ) 

622 

623 # Fit full model (Eq. 20 of Astier+19) and same model with 

624 # b=0 (c=0 in this code) 

625 pInit = np.concatenate((a0.ravel(), c0.ravel(), noise0.ravel(), np.array(gain0)), axis=None) 

626 functionsDict = {'fullModel': self.funcFullCovarianceModel, 

627 'fullModelNoB': self.funcFullCovarianceModelNoB} 

628 fitResults = {'fullModel': {'a': [], 'c': [], 'noise': [], 'gain': [], 'paramsErr': []}, 

629 'fullModelNoB': {'a': [], 'c': [], 'noise': [], 'gain': [], 'paramsErr': []}} 

630 for key in functionsDict: 

631 params, paramsErr, _ = fitLeastSq(pInit, muAtAmpMasked, 

632 covAtAmpForFitMasked.ravel(), functionsDict[key], 

633 weightsY=covSqrtWeightsAtAmpForFitMasked.ravel()) 

634 a = params[:lenParams].reshape((matrixSideFit, matrixSideFit)) 

635 c = params[lenParams:2*lenParams].reshape((matrixSideFit, matrixSideFit)) 

636 noise = params[2*lenParams:3*lenParams].reshape((matrixSideFit, matrixSideFit)) 

637 gain = params[-1] 

638 

639 fitResults[key]['a'] = a 

640 fitResults[key]['c'] = c 

641 fitResults[key]['noise'] = noise 

642 fitResults[key]['gain'] = gain 

643 fitResults[key]['paramsErr'] = paramsErr 

644 

645 # Put the information in the PTC dataset 

646 

647 # Not used when ptcFitType is 'FULLCOVARIANCE' 

648 dataset.ptcFitPars[ampName] = np.array([np.nan]) 

649 dataset.ptcFitParsError[ampName] = np.array([np.nan]) 

650 dataset.ptcFitChiSq[ampName] = np.nan 

651 

652 # Save full covariances, covariances models, and their weights. 

653 # dataset.expIdMask is already full, but needs to be 

654 # converted to bool. 

655 dataset.expIdMask[ampName] = np.array(dataset.expIdMask[ampName], dtype=bool) 

656 dataset.covariances[ampName] = covAtAmp 

657 # We evaluate the covariance model everywhere, even the 

658 # masked amps. 

659 dataset.covariancesModel[ampName] = self.evalCovModel(muAtAmp, 

660 fitResults['fullModel']['a'], 

661 fitResults['fullModel']['c'], 

662 fitResults['fullModel']['noise'], 

663 fitResults['fullModel']['gain']) 

664 dataset.covariancesSqrtWeights[ampName] = covSqrtWeightsAtAmp 

665 dataset.aMatrix[ampName] = fitResults['fullModel']['a'] 

666 dataset.bMatrix[ampName] = fitResults['fullModel']['c']/fitResults['fullModel']['a'] 

667 dataset.covariancesModelNoB[ampName] = self.evalCovModel(muAtAmp, 

668 fitResults['fullModelNoB']['a'], 

669 fitResults['fullModelNoB']['c'], 

670 fitResults['fullModelNoB']['noise'], 

671 fitResults['fullModelNoB']['gain'], 

672 setBtoZero=True) 

673 dataset.aMatrixNoB[ampName] = fitResults['fullModelNoB']['a'] 

674 dataset.gain[ampName] = fitResults['fullModel']['gain'] 

675 dataset.gainErr[ampName] = fitResults['fullModel']['paramsErr'][-1] 

676 readoutNoise = fitResults['fullModel']['noise'][0][0] 

677 readoutNoiseSqrt = np.sqrt(np.fabs(readoutNoise)) 

678 dataset.noise[ampName] = readoutNoise 

679 readoutNoiseSigma = fitResults['fullModel']['paramsErr'][2*lenParams] 

680 dataset.noiseErr[ampName] = 0.5*(readoutNoiseSigma/np.fabs(readoutNoise))*readoutNoiseSqrt 

681 dataset.noiseMatrix[ampName] = fitResults['fullModel']['noise'] 

682 dataset.noiseMatrixNoB[ampName] = fitResults['fullModelNoB']['noise'] 

683 

684 dataset.finalVars[ampName] = covAtAmp[:, 0, 0] 

685 dataset.finalModelVars[ampName] = dataset.covariancesModel[ampName][:, 0, 0] 

686 dataset.finalMeans[ampName] = muAtAmp 

687 

688 return dataset 

689 

690 def initialFitFullCovariance(self, mu, cov, sqrtW): 

691 """ Performs a crude parabolic fit of the data in order to start 

692 the full fit close to the solution, setting b=0 (c=0) in Eq. 20 

693 of Astier+19. 

694 

695 Parameters 

696 ---------- 

697 mu : `numpy.array`, (N,) 

698 Signal `mu` (ADU) 

699 cov : `numpy.array`, (N, M, M) 

700 Covariance arrays of size `(M, M)` (with 

701 `M = config.maximumRangeCovariancesAstier`), 

702 indexed by mean signal `mu`. 

703 sqrtW : `numpy.array`, (N,) 

704 Covariance weights, defined as 1./sqrt(Variances) 

705 

706 Returns 

707 ------- 

708 a : `numpy.array`, (M, M) 

709 "a" parameter per flux in Eq. 20 of Astier+19. 

710 c : `numpy.array`, (M, M) 

711 "c"="ab" parameter per flux in Eq. 20 of Astier+19. 

712 noise : `numpy.array`, (M, M) 

713 "noise" parameter per flux in Eq. 20 of Astier+19. 

714 gain : `float` 

715 Amplifier gain (e/ADU) 

716 """ 

717 matrixSideFit = self.config.maximumRangeCovariancesAstierFullCovFit 

718 

719 # Initialize fit parameters 

720 a = np.zeros((matrixSideFit, matrixSideFit)) 

721 c = np.zeros((matrixSideFit, matrixSideFit)) 

722 noise = np.zeros((matrixSideFit, matrixSideFit)) 

723 gain = 1. 

724 

725 # iterate the fit to account for higher orders 

726 # the chi2 does not necessarily go down, so one could 

727 # stop when it increases 

728 oldChi2 = 1e30 

729 for _ in range(5): 

730 model = np.nan_to_num(self.evalCovModel(mu, a, c, noise, gain, setBtoZero=True)) 

731 # loop on lags 

732 for i in range(matrixSideFit): 

733 for j in range(matrixSideFit): 

734 # fit a parabola for a given lag 

735 parsFit = np.polyfit(mu, cov[:, i, j] - model[:, i, j], 

736 2, w=sqrtW[:, i, j]) 

737 # model equation (Eq. 20) in Astier+19, with c=a*b=0: 

738 a[i, j] += parsFit[0] 

739 noise[i, j] += parsFit[2] 

740 if i + j == 0: 

741 gain = 1./(1/gain+parsFit[1]) 

742 weightedRes = (model - cov)*sqrtW 

743 chi2 = (weightedRes.flatten()**2).sum() 

744 if chi2 > oldChi2: 

745 break 

746 oldChi2 = chi2 

747 

748 return a, c, noise, gain 

749 

750 def funcFullCovarianceModel(self, params, x): 

751 """Model to fit covariances from flat fields; Equation 20 of 

752 Astier+19. 

753 

754 Parameters 

755 ---------- 

756 params : `list` 

757 Parameters of the model: aMatrix, CMatrix, noiseMatrix, 

758 gain (e/ADU). 

759 x : `numpy.array`, (N,) 

760 Signal `mu` (ADU) 

761 

762 Returns 

763 ------- 

764 y : `numpy.array`, (N,) 

765 Covariance matrix. 

766 """ 

767 matrixSideFit = self.config.maximumRangeCovariancesAstierFullCovFit 

768 lenParams = matrixSideFit*matrixSideFit 

769 aMatrix = params[:lenParams].reshape((matrixSideFit, matrixSideFit)) 

770 cMatrix = params[lenParams:2*lenParams].reshape((matrixSideFit, matrixSideFit)) 

771 noiseMatrix = params[2*lenParams:3*lenParams].reshape((matrixSideFit, matrixSideFit)) 

772 gain = params[-1] 

773 

774 return self.evalCovModel(x, aMatrix, cMatrix, noiseMatrix, gain).flatten() 

775 

776 def funcFullCovarianceModelNoB(self, params, x): 

777 """Model to fit covariances from flat fields; Equation 20 of 

778 Astier+19, with b=0 (equivalent to c=a*b=0 in this code). 

779 

780 Parameters 

781 ---------- 

782 params : `list` 

783 Parameters of the model: aMatrix, CMatrix, noiseMatrix, 

784 gain (e/ADU). 

785 x : `numpy.array`, (N,) 

786 Signal mu (ADU) 

787 

788 Returns 

789 ------- 

790 y : `numpy.array`, (N,) 

791 Covariance matrix. 

792 """ 

793 matrixSideFit = self.config.maximumRangeCovariancesAstierFullCovFit 

794 lenParams = matrixSideFit*matrixSideFit 

795 aMatrix = params[:lenParams].reshape((matrixSideFit, matrixSideFit)) 

796 cMatrix = params[lenParams:2*lenParams].reshape((matrixSideFit, matrixSideFit)) 

797 noiseMatrix = params[2*lenParams:3*lenParams].reshape((matrixSideFit, matrixSideFit)) 

798 gain = params[-1] 

799 

800 return self.evalCovModel(x, aMatrix, cMatrix, noiseMatrix, gain, setBtoZero=True).flatten() 

801 

802 def evalCovModel(self, mu, aMatrix, cMatrix, noiseMatrix, gain, setBtoZero=False): 

803 """Computes full covariances model (Eq. 20 of Astier+19). 

804 

805 Parameters 

806 ---------- 

807 mu : `numpy.array`, (N,) 

808 List of mean signals. 

809 aMatrix : `numpy.array`, (M, M) 

810 "a" parameter per flux in Eq. 20 of Astier+19. 

811 cMatrix : `numpy.array`, (M, M) 

812 "c"="ab" parameter per flux in Eq. 20 of Astier+19. 

813 noiseMatrix : `numpy.array`, (M, M) 

814 "noise" parameter per flux in Eq. 20 of Astier+19. 

815 gain : `float` 

816 Amplifier gain (e/ADU) 

817 setBtoZero=False : `bool`, optional 

818 Set "b" parameter in full model (see Astier+19) to zero. 

819 

820 Returns 

821 ------- 

822 covModel : `numpy.array`, (N, M, M) 

823 Covariances model. 

824 

825 Notes 

826 ----- 

827 By default, computes the covModel for the mu's stored(self.mu). 

828 Returns cov[Nmu, M, M]. The variance for the PTC is 

829 cov[:, 0, 0]. mu and cov are in ADUs and ADUs squared. To use 

830 electrons for both, the gain should be set to 1. This routine 

831 implements the model in Astier+19 (1905.08677). 

832 The parameters of the full model for C_ij(mu) ("C_ij" and "mu" 

833 in ADU^2 and ADU, respectively) in Astier+19 (Eq. 20) are: 

834 

835 - "a" coefficients (M by M matrix), units: 1/e 

836 - "b" coefficients (M by M matrix), units: 1/e 

837 - noise matrix (M by M matrix), units: e^2 

838 - gain, units: e/ADU 

839 

840 "b" appears in Eq. 20 only through the "ab" combination, which 

841 is defined in this code as "c=ab". 

842 """ 

843 matrixSideFit = self.config.maximumRangeCovariancesAstierFullCovFit 

844 sa = (matrixSideFit, matrixSideFit) 

845 # pad a with zeros and symmetrize 

846 aEnlarged = np.zeros((int(sa[0]*1.5)+1, int(sa[1]*1.5)+1)) 

847 aEnlarged[0:sa[0], 0:sa[1]] = aMatrix 

848 aSym = symmetrize(aEnlarged) 

849 # pad c with zeros and symmetrize 

850 cEnlarged = np.zeros((int(sa[0]*1.5)+1, int(sa[1]*1.5)+1)) 

851 cEnlarged[0:sa[0], 0:sa[1]] = cMatrix 

852 cSym = symmetrize(cEnlarged) 

853 a2 = fftconvolve(aSym, aSym, mode='same') 

854 a3 = fftconvolve(a2, aSym, mode='same') 

855 ac = fftconvolve(aSym, cSym, mode='same') 

856 (xc, yc) = np.unravel_index(np.abs(aSym).argmax(), a2.shape) 

857 

858 a1 = aMatrix[np.newaxis, :, :] 

859 a2 = a2[np.newaxis, xc:xc + matrixSideFit, yc:yc + matrixSideFit] 

860 a3 = a3[np.newaxis, xc:xc + matrixSideFit, yc:yc + matrixSideFit] 

861 ac = ac[np.newaxis, xc:xc + matrixSideFit, yc:yc + matrixSideFit] 

862 c1 = cMatrix[np.newaxis, ::] 

863 

864 # assumes that mu is 1d 

865 bigMu = mu[:, np.newaxis, np.newaxis]*gain 

866 # c(=a*b in Astier+19) also has a contribution to the last 

867 # term, that is absent for now. 

868 if setBtoZero: 

869 c1 = np.zeros_like(c1) 

870 ac = np.zeros_like(ac) 

871 covModel = (bigMu/(gain*gain)*(a1*bigMu+2./3.*(bigMu*bigMu)*(a2 + c1) 

872 + (1./3.*a3 + 5./6.*ac)*(bigMu*bigMu*bigMu)) + noiseMatrix[np.newaxis, :, :]/gain**2) 

873 # add the Poisson term, and the read out noise (variance) 

874 covModel[:, 0, 0] += mu/gain 

875 

876 return covModel 

877 

878 def subtractDistantOffset(self, muAtAmpMasked, covAtAmpMasked, covSqrtWeightsAtAmpMasked, 

879 start, degree=1): 

880 """Subtract distant offset from the covariance matrices. 

881 

882 Parameters 

883 ---------- 

884 muAtAmpMasked : `numpy.array` 

885 Masked mean flux array for a particular amplifier. 

886 covAtAmpMasked : `numpy.array` 

887 Masked measured covariances for a particular amplifier. 

888 covSqrtWeightsAtAmpMasked : `numpy.array` 

889 Masked inverse covariance weights for a particular amplifier. 

890 start : int, optional 

891 The starting index to eliminate the core for the fit. 

892 degree : int, optional 

893 Degree of the polynomial fit. 

894 

895 Returns 

896 ------- 

897 covAtAmpMasked : `numpy.array` 

898 Subtracted measured covariances for a particular amplifier. 

899 covSqrtWeightsAtAmpMasked : `numpy.array` 

900 Masked inverse covariance weights for a particular amplifier. 

901 

902 Notes 

903 ----- 

904 Ported from https://gitlab.in2p3.fr/astier/bfptc by P. Astier. 

905 

906 This function subtracts a distant offset from the 

907 covariance matrices using polynomial fitting. The core 

908 of the matrices is eliminated for the fit. 

909 

910 The function modifies the internal state of the object, updating the 

911 covariance matrices and related attributes. 

912 """ 

913 for k in range(len(muAtAmpMasked)): 

914 # Make a copy because it will be altered 

915 w = np.copy(covSqrtWeightsAtAmpMasked[k, ...]) 

916 wShape = w.shape 

917 i, j = np.meshgrid(range(wShape[0]), range(wShape[1]), indexing='ij') 

918 

919 # Eliminate the core for the fit 

920 w[:start, :start] = 0 

921 

922 poly = Pol2D(i, j, covAtAmpMasked[k, ...], degree, w=w) 

923 back = poly.eval(i, j) 

924 

925 covAtAmpMasked[k, ...] -= back 

926 

927 return covAtAmpMasked, covSqrtWeightsAtAmpMasked 

928 

929 # EXPAPPROXIMATION and POLYNOMIAL fit methods 

930 @staticmethod 

931 def _initialParsForPolynomial(order): 

932 assert order >= 2 

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

934 pars[0] = 10 

935 pars[1] = 1 

936 pars[2:] = 0.0001 

937 return pars 

938 

939 @staticmethod 

940 def _boundsForPolynomial(initialPars, lowers=[], uppers=[]): 

941 if not len(lowers): 

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

943 if not len(uppers): 

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

945 lowers[1] = 0 # no negative gains 

946 return (lowers, uppers) 

947 

948 @staticmethod 

949 def _boundsForAstier(initialPars, lowers=[], uppers=[]): 

950 if not len(lowers): 

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

952 if not len(uppers): 

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

954 return (lowers, uppers) 

955 

956 @staticmethod 

957 def _getInitialGoodPoints(means, variances, minVarPivotSearch, consecutivePointsVarDecreases): 

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

959 

960 Parameters 

961 ---------- 

962 means : `numpy.array` 

963 Input array with mean signal values. 

964 variances : `numpy.array` 

965 Input array with variances at each mean value. 

966 minVarPivotSearch : `float` 

967 The variance (in ADU^2), above which, the point 

968 of decreasing variance should be sought. 

969 consecutivePointsVarDecreases : `int` 

970 Required number of consecutive points/fluxes 

971 in the PTC where the variance 

972 decreases in order to find a first 

973 estimate of the PTC turn-off. 

974 

975 Returns 

976 ------ 

977 goodPoints : `numpy.array` [`bool`] 

978 Boolean array to select good (`True`) and bad (`False`) 

979 points. 

980 

981 Notes 

982 ----- 

983 Eliminate points beyond which the variance decreases. 

984 """ 

985 goodPoints = np.ones_like(means, dtype=bool) 

986 # Variances are sorted and should monotonically increase 

987 pivotList = np.where(np.array(np.diff(variances)) < 0)[0] 

988 if len(pivotList) > 0: 

989 # For small values, sometimes the variance decreases slightly 

990 # Only look when var > self.config.minVarPivotSearch 

991 pivotList = [p for p in pivotList if variances[p] > minVarPivotSearch] 

992 # Require that the varince decreases during 

993 # consecutivePointsVarDecreases 

994 # consecutive points. This will give a first 

995 # estimate of the PTC turn-off, which 

996 # may be updated (reduced) further in the code. 

997 if len(pivotList) > 1: 

998 # enumerate(pivotList) creates tuples (index, value), for 

999 # each value in pivotList. The lambda function subtracts 

1000 # each value from the index. 

1001 # groupby groups elements by equal key value. 

1002 for k, g in groupby(enumerate(pivotList), lambda x: x[0]-x[1]): 

1003 group = (map(itemgetter(1), g)) 

1004 # Form groups of consecute values from pivotList 

1005 group = list(map(int, group)) 

1006 # values in pivotList are indices where np.diff(variances) 

1007 # is negative, i.e., where the variance starts decreasing. 

1008 # Find the first group of consecutive numbers when 

1009 # variance decreases. 

1010 if len(group) >= consecutivePointsVarDecreases: 

1011 pivotIndex = np.min(group) 

1012 goodPoints[pivotIndex+1:] = False 

1013 break 

1014 

1015 # Finally, we filter out any infinities or NaNs. 

1016 goodPoints[(~np.isfinite(means)) | (~np.isfinite(variances))] = False 

1017 

1018 return goodPoints 

1019 

1020 def _makeZeroSafe(self, array, substituteValue=1e-9): 

1021 """""" 

1022 array = np.array(array) 

1023 nBad = Counter(np.ravel(array))[0] 

1024 if nBad == 0: 

1025 return array 

1026 

1027 index, = np.where(array == 0) 

1028 if len(index): 

1029 msg = f"Found {nBad} zeros in array at elements {index}" 

1030 self.log.warning(msg) 

1031 

1032 array[index] = substituteValue 

1033 

1034 return array 

1035 

1036 def fitPtc(self, dataset): 

1037 """Fit the photon transfer curve to a polynomial or to the 

1038 Astier+19 approximation (Eq. 16). 

1039 

1040 Fit the photon transfer curve with either a polynomial of 

1041 the order specified in the task config (POLNOMIAL), 

1042 or using the exponential approximation in Astier+19 

1043 (Eq. 16, FULLCOVARIANCE). 

1044 

1045 Sigma clipping is performed iteratively for the fit, as 

1046 well as an initial clipping of data points that are more 

1047 than `config.initialNonLinearityExclusionThreshold` away 

1048 from lying on a straight line. This other step is necessary 

1049 because the photon transfer curve turns over catastrophically 

1050 at very high flux (because saturation 

1051 drops the variance to ~0) and these far outliers cause the 

1052 initial fit to fail, meaning the sigma cannot be calculated 

1053 to perform the sigma-clipping. 

1054 

1055 Parameters 

1056 ---------- 

1057 dataset : `lsst.ip.isr.ptcDataset.PhotonTransferCurveDataset` 

1058 The dataset containing the means, variances and 

1059 exposure times. 

1060 

1061 Returns 

1062 ------- 

1063 dataset : `lsst.ip.isr.ptcDataset.PhotonTransferCurveDataset` 

1064 This is the same dataset as the input parameter, however, 

1065 it has been modified to include information such as the 

1066 fit vectors and the fit parameters. See the class 

1067 `PhotonTransferCurveDatase`. 

1068 

1069 Raises 

1070 ------ 

1071 RuntimeError 

1072 Raised if dataset.ptcFitType is None or empty. 

1073 """ 

1074 if dataset.ptcFitType: 

1075 ptcFitType = dataset.ptcFitType 

1076 else: 

1077 raise RuntimeError("ptcFitType is None of empty in PTC dataset.") 

1078 # For FULLCOVARIANCE model fit 

1079 matrixSideFit = dataset.covMatrixSideFullCovFit 

1080 nanMatrixFit = np.empty((matrixSideFit, matrixSideFit)) 

1081 nanMatrixFit[:] = np.nan 

1082 

1083 for amp in dataset.ampNames: 

1084 lenInputTimes = len(dataset.rawExpTimes[amp]) 

1085 listNanMatrixFit = np.empty((lenInputTimes, matrixSideFit, matrixSideFit)) 

1086 listNanMatrixFit[:] = np.nan 

1087 

1088 dataset.covariancesModel[amp] = listNanMatrixFit 

1089 dataset.aMatrix[amp] = nanMatrixFit 

1090 dataset.bMatrix[amp] = nanMatrixFit 

1091 dataset.covariancesModelNoB[amp] = listNanMatrixFit 

1092 dataset.aMatrixNoB[amp] = nanMatrixFit 

1093 dataset.noiseMatrix[amp] = nanMatrixFit 

1094 dataset.noiseMatrixNoB[amp] = nanMatrixFit 

1095 

1096 def errFunc(p, x, y): 

1097 return ptcFunc(p, x) - y 

1098 

1099 sigmaCutPtcOutliers = self.config.sigmaCutPtcOutliers 

1100 maxIterationsPtcOutliers = self.config.maxIterationsPtcOutliers 

1101 

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

1103 meanVecOriginal = dataset.rawMeans[ampName].copy() 

1104 varVecOriginal = dataset.rawVars[ampName].copy() 

1105 varVecOriginal = self._makeZeroSafe(varVecOriginal) 

1106 

1107 if self.config.doLegacyTurnoffSelection: 

1108 # Discard points when the variance starts to decrease after two 

1109 # consecutive signal levels 

1110 goodPoints = self._getInitialGoodPoints(meanVecOriginal, varVecOriginal, 

1111 self.config.minVarPivotSearch, 

1112 self.config.consecutivePointsVarDecreases) 

1113 else: 

1114 goodPoints = dataset.expIdMask[ampName] 

1115 

1116 # Check if all points are bad from the 'cpExtractPtcTask' 

1117 initialExpIdMask = dataset.expIdMask[ampName] 

1118 

1119 if not (goodPoints.any() and initialExpIdMask.any()): 

1120 msg = (f"SERIOUS: All points in goodPoints: {goodPoints} or " 

1121 f"in initialExpIdMask: {initialExpIdMask} are bad." 

1122 f"Setting {ampName} to BAD.") 

1123 self.log.warning(msg) 

1124 # Fill entries with NaNs 

1125 self.fillBadAmp(dataset, ptcFitType, ampName) 

1126 continue 

1127 

1128 mask = goodPoints 

1129 

1130 if ptcFitType == 'EXPAPPROXIMATION': 

1131 ptcFunc = funcAstier 

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

1133 # lowers and uppers obtained from BOT data studies by 

1134 # C. Lage (UC Davis, 11/2020). 

1135 if self.config.binSize > 1: 

1136 bounds = self._boundsForAstier(parsIniPtc) 

1137 else: 

1138 bounds = self._boundsForAstier(parsIniPtc, lowers=[-1e-4, 0.1, -2000], 

1139 uppers=[1e-4, 10.0, 2000]) 

1140 if ptcFitType == 'POLYNOMIAL': 

1141 ptcFunc = funcPolynomial 

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

1143 bounds = self._boundsForPolynomial(parsIniPtc) 

1144 

1145 # We perform an initial (unweighted) fit of variance vs signal 

1146 # (after initial KS test or post-drop selection) to look for 

1147 # outliers, particularly at the high-flux end. The initial fit 

1148 # is performed only for points that are guaranteed to be below 

1149 # the PTC turnoff and then extrapolated to ensure that high 

1150 # flux points that have abnormal variance values can be properly 

1151 # rejected in this phase without biasing the initial fit. 

1152 # This algorithm was initially developed by Seth Digel for 

1153 # the EO Testing pipeline. 

1154 

1155 if self.config.scaleMaxSignalInitialPtcOutlierFit: 

1156 approxGain = np.nanmedian(meanVecOriginal/varVecOriginal) 

1157 maxADUInitialPtcOutlierFit = self.config.maxSignalInitialPtcOutlierFit/approxGain 

1158 self.log.info( 

1159 "Using approximate gain %.3f and ADU signal cutoff of %.1f for amplifier %s", 

1160 approxGain, 

1161 maxADUInitialPtcOutlierFit, 

1162 ampName, 

1163 ) 

1164 else: 

1165 maxADUInitialPtcOutlierFit = self.config.maxSignalInitialPtcOutlierFit 

1166 

1167 if maxIterationsPtcOutliers == 0: 

1168 # We are not doing any outlier rejection here, but we do want 

1169 # an initial fit. 

1170 res = least_squares( 

1171 errFunc, 

1172 parsIniPtc, 

1173 bounds=bounds, 

1174 args=(meanVecOriginal[mask], varVecOriginal[mask]), 

1175 ) 

1176 pars = res.x 

1177 newMask = mask.copy() 

1178 else: 

1179 newMask = (mask & (meanVecOriginal <= maxADUInitialPtcOutlierFit)) 

1180 

1181 count = 0 

1182 lastMask = mask.copy() 

1183 while count < maxIterationsPtcOutliers: 

1184 res = least_squares( 

1185 errFunc, 

1186 parsIniPtc, 

1187 bounds=bounds, 

1188 args=(meanVecOriginal[newMask], varVecOriginal[newMask]), 

1189 ) 

1190 pars = res.x 

1191 

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

1193 # The new mask includes points where the residuals are 

1194 # finite, are less than the cut, and include the original 

1195 # mask of known points that should not be used. 

1196 newMask = ( 

1197 np.isfinite(sigResids) 

1198 & (np.abs(np.nan_to_num(sigResids)) < sigmaCutPtcOutliers) 

1199 & mask 

1200 ) 

1201 if np.count_nonzero(newMask) == 0: 

1202 msg = (f"SERIOUS: All points after outlier rejection are bad. " 

1203 f"Setting {ampName} to BAD.") 

1204 self.log.warning(msg) 

1205 # Fill entries with NaNs 

1206 self.fillBadAmp(dataset, ptcFitType, ampName) 

1207 break 

1208 

1209 self.log.debug( 

1210 "Iteration %d: Removed %d points in total for %s.", 

1211 count, 

1212 np.count_nonzero(mask) - np.count_nonzero(newMask), 

1213 ampName, 

1214 ) 

1215 

1216 # If the mask hasn't changed then break out. 

1217 if np.all(newMask == lastMask): 

1218 self.log.debug("Convergence at iteration %d; breaking loop for %s.", count, ampName) 

1219 break 

1220 

1221 lastMask = newMask.copy() 

1222 

1223 count += 1 

1224 

1225 # Set the mask to the new mask 

1226 mask = newMask.copy() 

1227 

1228 if not mask.any(): 

1229 # We hae already filled the bad amp above, so continue. 

1230 continue 

1231 

1232 dataset.expIdMask[ampName] = mask 

1233 

1234 parsIniPtc = pars 

1235 meanVecFinal = meanVecOriginal[mask] 

1236 varVecFinal = varVecOriginal[mask] 

1237 

1238 # Save the maximum point after outlier detection as the 

1239 # PTC turnoff point. 

1240 dataset.ptcTurnoff[ampName] = meanVecFinal[-1] 

1241 

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

1243 self.log.info("Number of points discarded in PTC of amplifier %s:" 

1244 " %d out of %d", ampName, Counter(mask)[False], len(meanVecOriginal)) 

1245 

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

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

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

1249 self.log.warning(msg) 

1250 # Fill entries with NaNs 

1251 self.fillBadAmp(dataset, ptcFitType, ampName) 

1252 continue 

1253 # Fit the PTC. 

1254 # The variance of the variance is Var(v)=2*v^2/Npix. This is 

1255 # already calculated in `makeCovArray` of CpPtcExtract. 

1256 # dataset.covariancesSqrtWeights[ampName][:,0,0] 

1257 # has 1/sqrt(Var(v)). 

1258 weightsY = dataset.covariancesSqrtWeights[ampName][:, 0, 0][mask] 

1259 if self.config.doFitBootstrap: 

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

1261 varVecFinal, ptcFunc, 

1262 weightsY=weightsY) 

1263 else: 

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

1265 varVecFinal, ptcFunc, 

1266 weightsY=weightsY) 

1267 dataset.ptcFitPars[ampName] = parsFit 

1268 dataset.ptcFitParsError[ampName] = parsFitErr 

1269 dataset.ptcFitChiSq[ampName] = reducedChiSqPtc 

1270 

1271 dataset.finalVars[ampName] = varVecOriginal 

1272 dataset.finalVars[ampName][~mask] = np.nan 

1273 dataset.finalModelVars[ampName] = ptcFunc(parsFit, meanVecOriginal) 

1274 dataset.finalModelVars[ampName][~mask] = np.nan 

1275 dataset.finalMeans[ampName] = meanVecOriginal 

1276 dataset.finalMeans[ampName][~mask] = np.nan 

1277 

1278 if ptcFitType == 'EXPAPPROXIMATION': 

1279 ptcGain = parsFit[1] 

1280 ptcGainErr = parsFitErr[1] 

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

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

1283 if ptcFitType == 'POLYNOMIAL': 

1284 ptcGain = 1./parsFit[1] 

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

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

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

1288 dataset.gain[ampName] = ptcGain 

1289 dataset.gainErr[ampName] = ptcGainErr 

1290 dataset.noise[ampName] = ptcNoise 

1291 dataset.noiseErr[ampName] = ptcNoiseErr 

1292 

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

1294 dataset.ptcFitType = ptcFitType 

1295 if len(dataset.badAmps) == 0: 

1296 dataset.badAmps = [] 

1297 

1298 return dataset 

1299 

1300 def fillBadAmp(self, dataset, ptcFitType, ampName): 

1301 """Fill the dataset with NaNs if there are not enough 

1302 good points. 

1303 

1304 Parameters 

1305 ---------- 

1306 dataset : `lsst.ip.isr.ptcDataset.PhotonTransferCurveDataset` 

1307 The dataset containing the means, variances and 

1308 exposure times. 

1309 ptcFitType : {'POLYNOMIAL', 'EXPAPPROXIMATION'} 

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

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

1312 ampName : `str` 

1313 Amplifier name. 

1314 """ 

1315 dataset.badAmps.append(ampName) 

1316 dataset.expIdMask[ampName] = np.repeat(False, len(dataset.rawExpTimes[ampName])) 

1317 dataset.gain[ampName] = np.nan 

1318 dataset.gainErr[ampName] = np.nan 

1319 dataset.noise[ampName] = np.nan 

1320 dataset.noiseErr[ampName] = np.nan 

1321 dataset.ptcFitPars[ampName] = (np.repeat(np.nan, self.config.polynomialFitDegree + 1) if 

1322 ptcFitType in ["POLYNOMIAL", ] else np.repeat(np.nan, 3)) 

1323 dataset.ptcFitParsError[ampName] = (np.repeat(np.nan, self.config.polynomialFitDegree + 1) if 

1324 ptcFitType in ["POLYNOMIAL", ] else np.repeat(np.nan, 3)) 

1325 dataset.ptcFitChiSq[ampName] = np.nan 

1326 dataset.ptcTurnoff[ampName] = np.nan 

1327 dataset.finalVars[ampName] = np.repeat(np.nan, len(dataset.rawExpTimes[ampName])) 

1328 dataset.finalModelVars[ampName] = np.repeat(np.nan, len(dataset.rawExpTimes[ampName])) 

1329 dataset.finalMeans[ampName] = np.repeat(np.nan, len(dataset.rawExpTimes[ampName])) 

1330 

1331 return