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

406 statements  

« prev     ^ index     » next       coverage.py v6.5.0, created at 2023-03-10 03:13 -0800

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, funcAstier, symmetrize) 

28 

29from scipy.signal import fftconvolve 

30from scipy.optimize import least_squares 

31from itertools import groupby 

32from operator import itemgetter 

33 

34import lsst.pipe.base.connectionTypes as cT 

35 

36from lsst.ip.isr import PhotonTransferCurveDataset 

37 

38from lsst.cp.pipe._lookupStaticCalibration import lookupStaticCalibration 

39 

40import copy 

41 

42 

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

44 

45 

46class PhotonTransferCurveSolveConnections(pipeBase.PipelineTaskConnections, 

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

48 inputCovariances = cT.Input( 

49 name="ptcCovariances", 

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

51 storageClass="PhotonTransferCurveDataset", 

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

53 isCalibration=True, 

54 multiple=True, 

55 ) 

56 camera = cT.PrerequisiteInput( 

57 name="camera", 

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

59 storageClass="Camera", 

60 dimensions=("instrument",), 

61 isCalibration=True, 

62 lookupFunction=lookupStaticCalibration, 

63 ) 

64 outputPtcDataset = cT.Output( 

65 name="ptcDatsetProposal", 

66 doc="Output proposed ptc dataset.", 

67 storageClass="PhotonTransferCurveDataset", 

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

69 multiple=False, 

70 isCalibration=True, 

71 ) 

72 

73 

74class PhotonTransferCurveSolveConfig(pipeBase.PipelineTaskConfig, 

75 pipelineConnections=PhotonTransferCurveSolveConnections): 

76 """Configuration for fitting measured covariances. 

77 """ 

78 

79 ptcFitType = pexConfig.ChoiceField( 

80 dtype=str, 

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

82 default="POLYNOMIAL", 

83 allowed={ 

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

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

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

87 } 

88 ) 

89 maximumRangeCovariancesAstier = pexConfig.Field( 

90 dtype=int, 

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

92 default=8, 

93 ) 

94 sigmaClipFullFitCovariancesAstier = pexConfig.Field( 

95 dtype=float, 

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

97 default=5.0, 

98 ) 

99 maxIterFullFitCovariancesAstier = pexConfig.Field( 

100 dtype=int, 

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

102 default=3, 

103 ) 

104 polynomialFitDegree = pexConfig.Field( 

105 dtype=int, 

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

107 default=3, 

108 ) 

109 sigmaCutPtcOutliers = pexConfig.Field( 

110 dtype=float, 

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

112 default=5.0, 

113 ) 

114 maxIterationsPtcOutliers = pexConfig.RangeField( 

115 dtype=int, 

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

117 default=2, 

118 min=0 

119 ) 

120 minVarPivotSearch = pexConfig.Field( 

121 dtype=float, 

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

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

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

125 " should be sought.", 

126 default=10000, 

127 ) 

128 consecutivePointsVarDecreases = pexConfig.RangeField( 

129 dtype=int, 

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

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

132 default=2, 

133 min=2 

134 ) 

135 doFitBootstrap = pexConfig.Field( 

136 dtype=bool, 

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

138 default=False, 

139 ) 

140 binSize = pexConfig.Field( 

141 dtype=int, 

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

143 default=1, 

144 ) 

145 

146 

147class PhotonTransferCurveSolveTask(pipeBase.PipelineTask): 

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

149 

150 The first task of the PTC measurement pipeline, 

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

152 before this task), produced a list of 

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

154 contains the mean signal and covariances of the 

155 difference image of the flat-field images taken at 

156 the same exposure time. The list also contains dummy 

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

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

159 match. 

160 

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

162 of individual PTC datasets produced 

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

164 dataset, discarding the dummy datset as appropiate. 

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

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

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

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

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

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

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

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

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

174 

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

176 of CCD sensors", arXiv:1905.08677 

177 """ 

178 

179 ConfigClass = PhotonTransferCurveSolveConfig 

180 _DefaultName = 'cpPhotonTransferCurveSolve' 

181 

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

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

184 

185 Parameters 

186 ---------- 

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

188 Butler to operate on. 

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

190 Input data refs to load. 

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

192 Output data refs to persist. 

193 """ 

194 inputs = butlerQC.get(inputRefs) 

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

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

197 butlerQC.put(outputs, outputRefs) 

198 

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

200 """Fit measured covariances to different models. 

201 

202 Parameters 

203 ---------- 

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

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

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

207 Input camera. 

208 detId : `int` 

209 Detector ID to locate the detector in the camera and 

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

211 metadata. 

212 Returns 

213 ------- 

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

215 The resultins structure contains: 

216 

217 ``outputPtcDatset`` 

218 Final PTC dataset, containing information such as the 

219 means, variances, and exposure times 

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

221 """ 

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

223 ampNames = [] 

224 for partialPtcDataset in inputCovariances: 

225 if partialPtcDataset.ptcFitType != 'DUMMY': 

226 ampNames = partialPtcDataset.ampNames 

227 break 

228 

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

230 datasetPtc = PhotonTransferCurveDataset(ampNames=ampNames, 

231 ptcFitType=self.config.ptcFitType, 

232 covMatrixSide=self.config.maximumRangeCovariancesAstier) 

233 for partialPtcDataset in inputCovariances: 

234 # Ignore dummy datasets 

235 if partialPtcDataset.ptcFitType == 'DUMMY': 

236 continue 

237 for ampName in ampNames: 

238 datasetPtc.inputExpIdPairs[ampName].append(partialPtcDataset.inputExpIdPairs[ampName]) 

239 if type(partialPtcDataset.rawExpTimes[ampName]) is list: 

240 datasetPtc.rawExpTimes[ampName].append(partialPtcDataset.rawExpTimes[ampName][0]) 

241 else: 

242 datasetPtc.rawExpTimes[ampName].append(partialPtcDataset.rawExpTimes[ampName]) 

243 if type(partialPtcDataset.rawMeans[ampName]) is list: 

244 datasetPtc.rawMeans[ampName].append(partialPtcDataset.rawMeans[ampName][0]) 

245 else: 

246 datasetPtc.rawMeans[ampName].append(partialPtcDataset.rawMeans[ampName]) 

247 if type(partialPtcDataset.rawVars[ampName]) is list: 

248 datasetPtc.rawVars[ampName].append(partialPtcDataset.rawVars[ampName][0]) 

249 else: 

250 datasetPtc.rawVars[ampName].append(partialPtcDataset.rawVars[ampName]) 

251 if type(partialPtcDataset.expIdMask[ampName]) is list: 

252 datasetPtc.expIdMask[ampName].append(partialPtcDataset.expIdMask[ampName][0]) 

253 else: 

254 datasetPtc.expIdMask[ampName].append(partialPtcDataset.expIdMask[ampName]) 

255 datasetPtc.covariances[ampName].append(np.array(partialPtcDataset.covariances[ampName][0])) 

256 datasetPtc.covariancesSqrtWeights[ampName].append( 

257 np.array(partialPtcDataset.covariancesSqrtWeights[ampName][0])) 

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

259 # rawMeans index 

260 for ampName in ampNames: 

261 index = np.argsort(np.ravel(np.array(datasetPtc.rawMeans[ampName]))) 

262 datasetPtc.inputExpIdPairs[ampName] = np.array(datasetPtc.inputExpIdPairs[ampName])[index] 

263 datasetPtc.rawExpTimes[ampName] = np.array(datasetPtc.rawExpTimes[ampName])[index] 

264 datasetPtc.rawMeans[ampName] = np.array(datasetPtc.rawMeans[ampName])[index] 

265 datasetPtc.rawVars[ampName] = np.array(datasetPtc.rawVars[ampName])[index] 

266 datasetPtc.expIdMask[ampName] = np.array(datasetPtc.expIdMask[ampName])[index] 

267 datasetPtc.covariances[ampName] = np.array(datasetPtc.covariances[ampName])[index] 

268 datasetPtc.covariancesSqrtWeights[ampName] = np.array( 

269 datasetPtc.covariancesSqrtWeights[ampName])[index] 

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

271 # Fit the measured covariances vs mean signal to 

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

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

274 # signal (mu) curve using the EXPAPPROXIMATION model 

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

276 # get the flat pairs that are masked. The 

277 # points at these fluxes will also be masked when 

278 # calculating the other elements of the covariance 

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

280 

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

282 tempDatasetPtc = copy.copy(datasetPtc) 

283 tempDatasetPtc.ptcFitType = "EXPAPPROXIMATION" 

284 tempDatasetPtc = self.fitMeasurementsToModel(tempDatasetPtc) 

285 

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

287 # previous fit. 

288 for ampName in datasetPtc.ampNames: 

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

290 datasetPtc.fitType = "FULLCOVARIANCE" 

291 datasetPtc = self.fitMeasurementsToModel(datasetPtc) 

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

293 # ("EXPAPPROXIMATION", "POLYNOMIAL") 

294 else: 

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

296 # approximation (Eq. 16). Fill up 

297 # PhotonTransferCurveDataset object. 

298 datasetPtc = self.fitMeasurementsToModel(datasetPtc) 

299 

300 if camera: 

301 detector = camera[detId] 

302 else: 

303 detector = None 

304 datasetPtc.updateMetadataFromExposures(inputCovariances) 

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

306 

307 return pipeBase.Struct( 

308 outputPtcDataset=datasetPtc, 

309 ) 

310 

311 def fitMeasurementsToModel(self, dataset): 

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

313 polynomial or one of the models in Astier+19 

314 (Eq. 16 or Eq.20). 

315 

316 Parameters 

317 ---------- 

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

319 The dataset containing information such as the means, 

320 (co)variances, and exposure times. 

321 

322 Returns 

323 ------- 

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

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

326 it has been modified to include information such as the 

327 fit vectors and the fit parameters. See the class 

328 `PhotonTransferCurveDatase`. 

329 """ 

330 fitType = dataset.ptcFitType 

331 if fitType in ["FULLCOVARIANCE", ]: 

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

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

334 # with variance = Cov_00 

335 dataset = self.fitDataFullCovariance(dataset) 

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

337 # The PTC is technically defined as variance vs signal 

338 dataset = self.fitPtc(dataset) 

339 else: 

340 raise RuntimeError( 

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

342 "'POLYNOMIAL', 'EXPAPPROXIMATION', or 'FULLCOVARIANCE'" 

343 ) 

344 

345 return dataset 

346 

347 def fitDataFullCovariance(self, dataset): 

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

349 Astier+19 (Eq. 20). 

350 

351 Parameters 

352 ---------- 

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

354 The dataset containing information such as the means, 

355 (co)variances, and exposure times. 

356 

357 Returns 

358 ------- 

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

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

361 it has been modified to include information such as the 

362 fit vectors and the fit parameters. See the class 

363 `PhotonTransferCurveDatase`. 

364 

365 Notes 

366 ----- 

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

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

369 

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

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

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

373 - gain, units: e/ADU 

374 

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

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

377 

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

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

380 maximum lag considered for the covariances calculation, and the 

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

382 have r^2 fewer entries. 

383 """ 

384 matrixSide = self.config.maximumRangeCovariancesAstier 

385 lenParams = matrixSide*matrixSide 

386 

387 for ampName in dataset.ampNames: 

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

389 # Not used when ptcFitType is 'FULLCOVARIANCE' 

390 dataset.ptcFitPars[ampName] = [np.nan] 

391 dataset.ptcFitParsError[ampName] = [np.nan] 

392 dataset.ptcFitChiSq[ampName] = np.nan 

393 

394 if ampName in dataset.badAmps: 

395 # Bad amp 

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

397 # with astropy.Table works. 

398 nanMatrix = np.full((matrixSide, matrixSide), np.nan) 

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

400 dataset.covariancesModel[ampName] = listNanMatrix 

401 dataset.covariancesSqrtWeights[ampName] = listNanMatrix 

402 dataset.aMatrix[ampName] = nanMatrix 

403 dataset.bMatrix[ampName] = nanMatrix 

404 dataset.covariancesModelNoB[ampName] = listNanMatrix 

405 dataset.aMatrixNoB[ampName] = nanMatrix 

406 

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

408 dataset.gain[ampName] = np.nan 

409 dataset.gainErr[ampName] = np.nan 

410 dataset.noise[ampName] = np.nan 

411 dataset.noiseErr[ampName] = np.nan 

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

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

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

415 continue 

416 

417 muAtAmp = dataset.rawMeans[ampName] 

418 maskAtAmp = dataset.expIdMask[ampName] 

419 if len(maskAtAmp) == 0: 

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

421 

422 muAtAmp = muAtAmp[maskAtAmp] 

423 covAtAmp = np.nan_to_num(dataset.covariances[ampName])[maskAtAmp] 

424 covSqrtWeightsAtAmp = np.nan_to_num(dataset.covariancesSqrtWeights[ampName])[maskAtAmp] 

425 

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

427 a0, c0, noise0, gain0 = self.initialFitFullCovariance(muAtAmp, covAtAmp, covSqrtWeightsAtAmp) 

428 

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

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

431 pInit = np.concatenate((a0.flatten(), c0.flatten(), noise0.flatten(), np.array(gain0)), axis=None) 

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

433 'fullModelNoB': self.funcFullCovarianceModelNoB} 

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

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

436 for key in functionsDict: 

437 params, paramsErr, _ = fitLeastSq(pInit, muAtAmp, 

438 covAtAmp.flatten(), functionsDict[key], 

439 weightsY=covSqrtWeightsAtAmp.flatten()) 

440 a = params[:lenParams].reshape((matrixSide, matrixSide)) 

441 c = params[lenParams:2*lenParams].reshape((matrixSide, matrixSide)) 

442 noise = params[2*lenParams:3*lenParams].reshape((matrixSide, matrixSide)) 

443 gain = params[-1] 

444 

445 fitResults[key]['a'] = a 

446 fitResults[key]['c'] = c 

447 fitResults[key]['noise'] = noise 

448 fitResults[key]['gain'] = gain 

449 fitResults[key]['paramsErr'] = paramsErr 

450 

451 # Put the information in the PTC dataset 

452 

453 # Not used when ptcFitType is 'FULLCOVARIANCE' 

454 dataset.ptcFitPars[ampName] = [np.nan] 

455 dataset.ptcFitParsError[ampName] = [np.nan] 

456 dataset.ptcFitChiSq[ampName] = np.nan 

457 

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

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

460 # converted to bool. 

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

462 dataset.covariances[ampName] = covAtAmp 

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

464 fitResults['fullModel']['a'], 

465 fitResults['fullModel']['c'], 

466 fitResults['fullModel']['noise'], 

467 fitResults['fullModel']['gain']) 

468 dataset.covariancesSqrtWeights[ampName] = covSqrtWeightsAtAmp 

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

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

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

472 fitResults['fullModelNoB']['a'], 

473 fitResults['fullModelNoB']['c'], 

474 fitResults['fullModelNoB']['noise'], 

475 fitResults['fullModelNoB']['gain'], 

476 setBtoZero=True) 

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

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

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

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

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

482 dataset.noise[ampName] = readoutNoise 

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

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

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

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

487 dataset.finalMeans[ampName] = muAtAmp 

488 

489 return dataset 

490 

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

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

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

494 of Astier+19. 

495 

496 Parameters 

497 ---------- 

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

499 Signal `mu` (ADU) 

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

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

502 `M = config.maximumRangeCovariancesAstier`), 

503 indexed by mean signal `mu`. 

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

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

506 

507 Returns 

508 ------- 

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

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

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

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

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

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

515 gain : `float` 

516 Amplifier gain (e/ADU) 

517 """ 

518 matrixSide = self.config.maximumRangeCovariancesAstier 

519 

520 # Initialize fit parameters 

521 a = np.zeros((matrixSide, matrixSide)) 

522 c = np.zeros((matrixSide, matrixSide)) 

523 noise = np.zeros((matrixSide, matrixSide)) 

524 gain = 1. 

525 

526 # iterate the fit to account for higher orders 

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

528 # stop when it increases 

529 oldChi2 = 1e30 

530 for _ in range(5): 

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

532 # loop on lags 

533 for i in range(matrixSide): 

534 for j in range(matrixSide): 

535 # fit a parabola for a given lag 

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

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

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

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

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

541 if(i + j == 0): 

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

543 weightedRes = (model - cov)*sqrtW 

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

545 if chi2 > oldChi2: 

546 break 

547 oldChi2 = chi2 

548 

549 return a, c, noise, gain 

550 

551 def funcFullCovarianceModel(self, params, x): 

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

553 Astier+19. 

554 

555 Parameters 

556 ---------- 

557 params : `list` 

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

559 gain (e/ADU). 

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

561 Signal `mu` (ADU) 

562 

563 Returns 

564 ------- 

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

566 Covariance matrix. 

567 """ 

568 matrixSide = self.config.maximumRangeCovariancesAstier 

569 lenParams = matrixSide*matrixSide 

570 aMatrix = params[:lenParams].reshape((matrixSide, matrixSide)) 

571 cMatrix = params[lenParams:2*lenParams].reshape((matrixSide, matrixSide)) 

572 noiseMatrix = params[2*lenParams:3*lenParams].reshape((matrixSide, matrixSide)) 

573 gain = params[-1] 

574 

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

576 

577 def funcFullCovarianceModelNoB(self, params, x): 

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

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

580 

581 Parameters 

582 ---------- 

583 params : `list` 

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

585 gain (e/ADU). 

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

587 Signal mu (ADU) 

588 

589 Returns 

590 ------- 

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

592 Covariance matrix. 

593 """ 

594 matrixSide = self.config.maximumRangeCovariancesAstier 

595 lenParams = matrixSide*matrixSide 

596 aMatrix = params[:lenParams].reshape((matrixSide, matrixSide)) 

597 cMatrix = params[lenParams:2*lenParams].reshape((matrixSide, matrixSide)) 

598 noiseMatrix = params[2*lenParams:3*lenParams].reshape((matrixSide, matrixSide)) 

599 gain = params[-1] 

600 

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

602 

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

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

605 

606 Parameters 

607 ---------- 

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

609 List of mean signals. 

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

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

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

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

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

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

616 gain : `float` 

617 Amplifier gain (e/ADU) 

618 setBtoZero=False : `bool`, optional 

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

620 

621 Returns 

622 ------- 

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

624 Covariances model. 

625 

626 Notes 

627 ----- 

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

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

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

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

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

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

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

635 

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

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

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

639 - gain, units: e/ADU 

640 

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

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

643 """ 

644 matrixSide = self.config.maximumRangeCovariancesAstier 

645 sa = (matrixSide, matrixSide) 

646 # pad a with zeros and symmetrize 

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

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

649 aSym = symmetrize(aEnlarged) 

650 # pad c with zeros and symmetrize 

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

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

653 cSym = symmetrize(cEnlarged) 

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

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

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

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

658 

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

660 a2 = a2[np.newaxis, xc:xc + matrixSide, yc:yc + matrixSide] 

661 a3 = a3[np.newaxis, xc:xc + matrixSide, yc:yc + matrixSide] 

662 ac = ac[np.newaxis, xc:xc + matrixSide, yc:yc + matrixSide] 

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

664 

665 # assumes that mu is 1d 

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

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

668 # term, that is absent for now. 

669 if setBtoZero: 

670 c1 = np.zeros_like(c1) 

671 ac = np.zeros_like(ac) 

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

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

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

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

676 

677 return covModel 

678 

679 # EXPAPPROXIMATION and POLYNOMIAL fit methods 

680 @staticmethod 

681 def _initialParsForPolynomial(order): 

682 assert(order >= 2) 

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

684 pars[0] = 10 

685 pars[1] = 1 

686 pars[2:] = 0.0001 

687 return pars 

688 

689 @staticmethod 

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

691 if not len(lowers): 

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

693 if not len(uppers): 

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

695 lowers[1] = 0 # no negative gains 

696 return (lowers, uppers) 

697 

698 @staticmethod 

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

700 if not len(lowers): 

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

702 if not len(uppers): 

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

704 return (lowers, uppers) 

705 

706 @staticmethod 

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

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

709 

710 Parameters 

711 ---------- 

712 means : `numpy.array` 

713 Input array with mean signal values. 

714 variances : `numpy.array` 

715 Input array with variances at each mean value. 

716 minVarPivotSearch : `float` 

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

718 of decreasing variance should be sought. 

719 consecutivePointsVarDecreases : `int` 

720 Required number of consecutive points/fluxes 

721 in the PTC where the variance 

722 decreases in order to find a first 

723 estimate of the PTC turn-off. 

724 

725 Returns 

726 ------ 

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

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

729 points. 

730 

731 Notes 

732 ----- 

733 Eliminate points beyond which the variance decreases. 

734 """ 

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

736 # Variances are sorted and should monotonically increase 

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

738 if len(pivotList) > 0: 

739 # For small values, sometimes the variance decreases slightly 

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

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

742 # Require that the varince decreases during 

743 # consecutivePointsVarDecreases 

744 # consecutive points. This will give a first 

745 # estimate of the PTC turn-off, which 

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

747 if len(pivotList) > 1: 

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

749 # each value in pivotList. The lambda function subtracts 

750 # each value from the index. 

751 # groupby groups elements by equal key value. 

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

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

754 # Form groups of consecute values from pivotList 

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

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

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

758 # Find the first group of consecutive numbers when 

759 # variance decreases. 

760 if len(group) >= consecutivePointsVarDecreases: 

761 pivotIndex = np.min(group) 

762 goodPoints[pivotIndex+1:] = False 

763 break 

764 

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

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

767 

768 return goodPoints 

769 

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

771 """""" 

772 array = np.array(array) 

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

774 if nBad == 0: 

775 return array 

776 

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

778 if len(index): 

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

780 self.log.warning(msg) 

781 

782 array[index] = substituteValue 

783 

784 return array 

785 

786 def fitPtc(self, dataset): 

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

788 Astier+19 approximation (Eq. 16). 

789 

790 Fit the photon transfer curve with either a polynomial of 

791 the order specified in the task config, or using the 

792 exponential approximation in Astier+19 (Eq. 16). 

793 

794 Sigma clipping is performed iteratively for the fit, as 

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

796 than `config.initialNonLinearityExclusionThreshold` away 

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

798 because the photon transfer curve turns over catastrophically 

799 at very high flux (because saturation 

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

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

802 to perform the sigma-clipping. 

803 

804 Parameters 

805 ---------- 

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

807 The dataset containing the means, variances and 

808 exposure times. 

809 

810 Returns 

811 ------- 

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

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

814 it has been modified to include information such as the 

815 fit vectors and the fit parameters. See the class 

816 `PhotonTransferCurveDatase`. 

817 

818 Raises 

819 ------ 

820 RuntimeError 

821 Raised if dataset.ptcFitType is None or empty. 

822 """ 

823 if dataset.ptcFitType: 

824 ptcFitType = dataset.ptcFitType 

825 else: 

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

827 matrixSide = self.config.maximumRangeCovariancesAstier 

828 nanMatrix = np.empty((matrixSide, matrixSide)) 

829 nanMatrix[:] = np.nan 

830 

831 for amp in dataset.ampNames: 

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

833 listNanMatrix = np.empty((lenInputTimes, matrixSide, matrixSide)) 

834 listNanMatrix[:] = np.nan 

835 

836 dataset.covariancesModel[amp] = listNanMatrix 

837 dataset.aMatrix[amp] = nanMatrix 

838 dataset.bMatrix[amp] = nanMatrix 

839 dataset.covariancesModelNoB[amp] = listNanMatrix 

840 dataset.aMatrixNoB[amp] = nanMatrix 

841 

842 def errFunc(p, x, y): 

843 return ptcFunc(p, x) - y 

844 

845 sigmaCutPtcOutliers = self.config.sigmaCutPtcOutliers 

846 maxIterationsPtcOutliers = self.config.maxIterationsPtcOutliers 

847 

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

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

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

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

852 varVecOriginal = self._makeZeroSafe(varVecOriginal) 

853 

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

855 # consecutive signal levels 

856 goodPoints = self._getInitialGoodPoints(meanVecOriginal, varVecOriginal, 

857 self.config.minVarPivotSearch, 

858 self.config.consecutivePointsVarDecreases) 

859 

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

861 initialExpIdMask = np.ravel(np.array(dataset.expIdMask[ampName])) 

862 

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

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

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

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

867 self.log.warning(msg) 

868 # Fill entries with NaNs 

869 self.fillBadAmp(dataset, ptcFitType, ampName) 

870 continue 

871 

872 # Save the point where the variance starts decreasing as the 

873 # PTC turnoff point 

874 ptcTurnoff = meanVecOriginal[goodPoints][-1] 

875 dataset.ptcTurnoff[ampName] = ptcTurnoff 

876 

877 mask = goodPoints 

878 

879 if ptcFitType == 'EXPAPPROXIMATION': 

880 ptcFunc = funcAstier 

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

882 # lowers and uppers obtained from BOT data studies by 

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

884 if self.config.binSize > 1: 

885 bounds = self._boundsForAstier(parsIniPtc) 

886 else: 

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

888 uppers=[1e-4, 2.5, 2000]) 

889 if ptcFitType == 'POLYNOMIAL': 

890 ptcFunc = funcPolynomial 

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

892 bounds = self._boundsForPolynomial(parsIniPtc) 

893 

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

895 # This further process of outlier rejection be skipped 

896 # if self.config.maxIterationsPtcOutliers = 0. 

897 # We already did some initial outlier rejection above in 

898 # self._getInitialGoodPoints. 

899 count = 1 

900 newMask = np.ones_like(meanVecOriginal, dtype=bool) 

901 pars = parsIniPtc 

902 while count <= maxIterationsPtcOutliers: 

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

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

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

906 meanTempVec = meanVecOriginal[mask] 

907 varTempVec = varVecOriginal[mask] 

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

909 pars = res.x 

910 

911 # change this to the original from the temp because 

912 # the masks are ANDed meaning once a point is masked 

913 # it's always masked, and the masks must always be the 

914 # same length for broadcasting 

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

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

917 mask = mask & newMask 

918 if not (mask.any() and newMask.any()): 

919 msg = (f"SERIOUS: All points in either mask: {mask} or newMask: {newMask} are bad. " 

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

921 self.log.warning(msg) 

922 # Fill entries with NaNs 

923 self.fillBadAmp(dataset, ptcFitType, ampName) 

924 break 

925 nDroppedTotal = Counter(mask)[False] 

926 self.log.debug("Iteration %d: discarded %d points in total for %s", 

927 count, nDroppedTotal, ampName) 

928 count += 1 

929 # objects should never shrink 

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

931 if not (mask.any() and newMask.any()): 

932 continue 

933 dataset.expIdMask[ampName] = np.array(dataset.expIdMask[ampName]) 

934 # store the final mask 

935 if len(dataset.expIdMask[ampName]): 

936 dataset.expIdMask[ampName] &= mask # bitwise_and if there is already a mask 

937 else: 

938 dataset.expIdMask[ampName] = mask 

939 # In case there was a previous mask stored 

940 mask = dataset.expIdMask[ampName] 

941 parsIniPtc = pars 

942 meanVecFinal = meanVecOriginal[mask] 

943 varVecFinal = varVecOriginal[mask] 

944 

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

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

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

948 

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

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

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

952 self.log.warning(msg) 

953 # Fill entries with NaNs 

954 self.fillBadAmp(dataset, ptcFitType, ampName) 

955 continue 

956 # Fit the PTC. 

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

958 # already calculated in `makeCovArray` of CpPtcExtract. 

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

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

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

962 if self.config.doFitBootstrap: 

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

964 varVecFinal, ptcFunc, 

965 weightsY=weightsY) 

966 else: 

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

968 varVecFinal, ptcFunc, 

969 weightsY=weightsY) 

970 dataset.ptcFitPars[ampName] = parsFit 

971 dataset.ptcFitParsError[ampName] = parsFitErr 

972 dataset.ptcFitChiSq[ampName] = reducedChiSqPtc 

973 # Masked variances (measured and modeled) and means. Need 

974 # to pad the array so astropy.Table does not crash (the 

975 # mask may vary per amp). 

976 padLength = len(dataset.rawExpTimes[ampName]) - len(varVecFinal) 

977 dataset.finalVars[ampName] = np.pad(varVecFinal, (0, padLength), 'constant', 

978 constant_values=np.nan) 

979 dataset.finalModelVars[ampName] = np.pad(ptcFunc(parsFit, meanVecFinal), (0, padLength), 

980 'constant', constant_values=np.nan) 

981 dataset.finalMeans[ampName] = np.pad(meanVecFinal, (0, padLength), 'constant', 

982 constant_values=np.nan) 

983 if ptcFitType == 'EXPAPPROXIMATION': 

984 ptcGain = parsFit[1] 

985 ptcGainErr = parsFitErr[1] 

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

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

988 if ptcFitType == 'POLYNOMIAL': 

989 ptcGain = 1./parsFit[1] 

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

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

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

993 dataset.gain[ampName] = ptcGain 

994 dataset.gainErr[ampName] = ptcGainErr 

995 dataset.noise[ampName] = ptcNoise 

996 dataset.noiseErr[ampName] = ptcNoiseErr 

997 

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

999 dataset.ptcFitType = ptcFitType 

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

1001 dataset.badAmps = np.repeat(np.nan, len(list(dataset.rawExpTimes.values())[0])) 

1002 

1003 return dataset 

1004 

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

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

1007 good points. 

1008 

1009 Parameters 

1010 ---------- 

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

1012 The dataset containing the means, variances and 

1013 exposure times. 

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

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

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

1017 ampName : `str` 

1018 Amplifier name. 

1019 """ 

1020 dataset.badAmps.append(ampName) 

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

1022 dataset.gain[ampName] = np.nan 

1023 dataset.gainErr[ampName] = np.nan 

1024 dataset.noise[ampName] = np.nan 

1025 dataset.noiseErr[ampName] = np.nan 

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

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

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

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

1030 dataset.ptcFitChiSq[ampName] = np.nan 

1031 dataset.ptcTurnoff[ampName] = np.nan 

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

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

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

1035 

1036 return