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

439 statements  

« prev     ^ index     » next       coverage.py v7.2.7, created at 2023-06-11 02:49 -0700

1# This file is part of cp_pipe. 

2# 

3# Developed for the LSST Data Management System. 

4# This product includes software developed by the LSST Project 

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

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

7# for details of code ownership. 

8# 

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

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

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

12# (at your option) any later version. 

13# 

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

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

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

17# GNU General Public License for more details. 

18# 

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

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

21# 

22import numpy as np 

23from 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 minMeanSignal = pexConfig.DictField( 

90 keytype=str, 

91 itemtype=float, 

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

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

94 " {'ALL_AMPS': value}", 

95 default={'ALL_AMPS': 0.0}, 

96 ) 

97 maxMeanSignal = pexConfig.DictField( 

98 keytype=str, 

99 itemtype=float, 

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

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

102 " {'ALL_AMPS': value}", 

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

104 ) 

105 maximumRangeCovariancesAstier = pexConfig.Field( 

106 dtype=int, 

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

108 default=8, 

109 ) 

110 sigmaClipFullFitCovariancesAstier = pexConfig.Field( 

111 dtype=float, 

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

113 default=5.0, 

114 ) 

115 maxIterFullFitCovariancesAstier = pexConfig.Field( 

116 dtype=int, 

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

118 default=3, 

119 ) 

120 polynomialFitDegree = pexConfig.Field( 

121 dtype=int, 

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

123 default=3, 

124 ) 

125 doLegacyTurnoffSelection = pexConfig.Field( 

126 dtype=bool, 

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

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

129 default=False, 

130 ) 

131 sigmaCutPtcOutliers = pexConfig.Field( 

132 dtype=float, 

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

134 default=5.0, 

135 ) 

136 maxIterationsPtcOutliers = pexConfig.RangeField( 

137 dtype=int, 

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

139 default=2, 

140 min=0 

141 ) 

142 maxSignalInitialPtcOutlierFit = pexConfig.Field( 

143 dtype=float, 

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

145 "the PTC turnoff to ensure accurate outlier rejection.", 

146 default=30_000., 

147 ) 

148 minVarPivotSearch = pexConfig.Field( 

149 dtype=float, 

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

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

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

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

154 default=10000, 

155 ) 

156 consecutivePointsVarDecreases = pexConfig.RangeField( 

157 dtype=int, 

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

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

160 "Only used if doLegacyTurnoffSelection is True.", 

161 default=2, 

162 min=2 

163 ) 

164 ksTestMinPvalue = pexConfig.Field( 

165 dtype=float, 

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

167 "Only used if doLegacyTurnoffSelection is False.", 

168 default=0.01, 

169 ) 

170 doFitBootstrap = pexConfig.Field( 

171 dtype=bool, 

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

173 default=False, 

174 ) 

175 binSize = pexConfig.Field( 

176 dtype=int, 

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

178 default=1, 

179 ) 

180 

181 

182class PhotonTransferCurveSolveTask(pipeBase.PipelineTask): 

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

184 

185 The first task of the PTC measurement pipeline, 

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

187 before this task), produced a list of 

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

189 contains the mean signal and covariances of the 

190 difference image of the flat-field images taken at 

191 the same exposure time. The list also contains dummy 

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

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

194 match. 

195 

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

197 of individual PTC datasets produced 

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

199 dataset, discarding the dummy datset as appropiate. 

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

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

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

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

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

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

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

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

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

209 

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

211 of CCD sensors", arXiv:1905.08677 

212 """ 

213 

214 ConfigClass = PhotonTransferCurveSolveConfig 

215 _DefaultName = 'cpPhotonTransferCurveSolve' 

216 

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

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

219 

220 Parameters 

221 ---------- 

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

223 Butler to operate on. 

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

225 Input data refs to load. 

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

227 Output data refs to persist. 

228 """ 

229 inputs = butlerQC.get(inputRefs) 

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

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

232 butlerQC.put(outputs, outputRefs) 

233 

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

235 """Fit measured covariances to different models. 

236 

237 Parameters 

238 ---------- 

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

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

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

242 Input camera. 

243 detId : `int` 

244 Detector ID to locate the detector in the camera and 

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

246 metadata. 

247 Returns 

248 ------- 

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

250 The resultins structure contains: 

251 

252 ``outputPtcDatset`` 

253 Final PTC dataset, containing information such as the 

254 means, variances, and exposure times 

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

256 """ 

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

258 ampNames = [] 

259 for partialPtcDataset in inputCovariances: 

260 if partialPtcDataset.ptcFitType != 'DUMMY': 

261 ampNames = partialPtcDataset.ampNames 

262 break 

263 

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

265 # specified in the config. 

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

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

268 for ampName in ampNames: 

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

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

271 elif ampName in self.config.maxMeanSignal: 

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

273 

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

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

276 elif ampName in self.config.minMeanSignal: 

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

278 

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

280 datasetPtc = PhotonTransferCurveDataset(ampNames=ampNames, 

281 ptcFitType=self.config.ptcFitType, 

282 covMatrixSide=self.config.maximumRangeCovariancesAstier) 

283 for partialPtcDataset in inputCovariances: 

284 # Ignore dummy datasets 

285 if partialPtcDataset.ptcFitType == 'DUMMY': 

286 continue 

287 for ampName in ampNames: 

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

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

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

291 # (and only) element of the list. 

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

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

294 partialPtcDataset.rawExpTimes[ampName][0]) 

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

296 partialPtcDataset.rawMeans[ampName][0]) 

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

298 partialPtcDataset.rawVars[ampName][0]) 

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

300 partialPtcDataset.histVars[ampName][0]) 

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

302 partialPtcDataset.histChi2Dofs[ampName][0]) 

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

304 partialPtcDataset.kspValues[ampName][0]) 

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

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

307 partialPtcDataset.covariances[ampName].ravel() 

308 ).reshape( 

309 ( 

310 len(datasetPtc.rawExpTimes[ampName]), 

311 datasetPtc.covMatrixSide, 

312 datasetPtc.covMatrixSide, 

313 ) 

314 ) 

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

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

317 partialPtcDataset.covariancesSqrtWeights[ampName].ravel() 

318 ).reshape( 

319 ( 

320 len(datasetPtc.rawExpTimes[ampName]), 

321 datasetPtc.covMatrixSide, 

322 datasetPtc.covMatrixSide, 

323 ) 

324 ) 

325 

326 # Apply min/max masking. 

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

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

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

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

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

332 expIdMask = False 

333 

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

335 if not self.config.doLegacyTurnoffSelection and \ 

336 kspValue < self.config.ksTestMinPvalue: 

337 expIdMask = False 

338 

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

340 

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

342 # rawMeans index. 

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

344 # all sorted the same way. 

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

346 

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

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

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

350 if good.size == 0: 

351 detectorMeans[i] = np.nan 

352 else: 

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

354 

355 index = np.argsort(detectorMeans) 

356 

357 for ampName in ampNames: 

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

359 datasetPtc.inputExpIdPairs[ampName] 

360 )[index].tolist() 

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

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

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

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

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

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

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

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

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

370 

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

372 # Fit the measured covariances vs mean signal to 

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

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

375 # signal (mu) curve using the EXPAPPROXIMATION model 

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

377 # get the flat pairs that are masked. The 

378 # points at these fluxes will also be masked when 

379 # calculating the other elements of the covariance 

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

381 

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

383 tempDatasetPtc = copy.copy(datasetPtc) 

384 tempDatasetPtc.ptcFitType = "EXPAPPROXIMATION" 

385 tempDatasetPtc = self.fitMeasurementsToModel(tempDatasetPtc) 

386 

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

388 # previous fit. 

389 for ampName in datasetPtc.ampNames: 

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

391 datasetPtc.fitType = "FULLCOVARIANCE" 

392 datasetPtc = self.fitMeasurementsToModel(datasetPtc) 

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

394 # ("EXPAPPROXIMATION", "POLYNOMIAL") 

395 else: 

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

397 # approximation (Eq. 16). Fill up 

398 # PhotonTransferCurveDataset object. 

399 datasetPtc = self.fitMeasurementsToModel(datasetPtc) 

400 

401 if camera: 

402 detector = camera[detId] 

403 else: 

404 detector = None 

405 datasetPtc.updateMetadataFromExposures(inputCovariances) 

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

407 

408 return pipeBase.Struct( 

409 outputPtcDataset=datasetPtc, 

410 ) 

411 

412 def fitMeasurementsToModel(self, dataset): 

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

414 polynomial or one of the models in Astier+19 

415 (Eq. 16 or Eq.20). 

416 

417 Parameters 

418 ---------- 

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

420 The dataset containing information such as the means, 

421 (co)variances, and exposure times. 

422 

423 Returns 

424 ------- 

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

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

427 it has been modified to include information such as the 

428 fit vectors and the fit parameters. See the class 

429 `PhotonTransferCurveDatase`. 

430 """ 

431 fitType = dataset.ptcFitType 

432 if fitType in ["FULLCOVARIANCE", ]: 

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

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

435 # with variance = Cov_00 

436 dataset = self.fitDataFullCovariance(dataset) 

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

438 # The PTC is technically defined as variance vs signal 

439 dataset = self.fitPtc(dataset) 

440 else: 

441 raise RuntimeError( 

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

443 "'POLYNOMIAL', 'EXPAPPROXIMATION', or 'FULLCOVARIANCE'" 

444 ) 

445 

446 return dataset 

447 

448 def fitDataFullCovariance(self, dataset): 

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

450 Astier+19 (Eq. 20). 

451 

452 Parameters 

453 ---------- 

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

455 The dataset containing information such as the means, 

456 (co)variances, and exposure times. 

457 

458 Returns 

459 ------- 

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

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

462 it has been modified to include information such as the 

463 fit vectors and the fit parameters. See the class 

464 `PhotonTransferCurveDatase`. 

465 

466 Notes 

467 ----- 

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

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

470 

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

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

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

474 - gain, units: e/ADU 

475 

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

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

478 

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

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

481 maximum lag considered for the covariances calculation, and the 

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

483 have r^2 fewer entries. 

484 """ 

485 matrixSide = self.config.maximumRangeCovariancesAstier 

486 lenParams = matrixSide*matrixSide 

487 

488 for ampName in dataset.ampNames: 

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

490 # Not used when ptcFitType is 'FULLCOVARIANCE' 

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

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

493 dataset.ptcFitChiSq[ampName] = np.nan 

494 

495 if ampName in dataset.badAmps: 

496 # Bad amp 

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

498 # with astropy.Table works. 

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

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

501 dataset.covariancesModel[ampName] = listNanMatrix 

502 dataset.covariancesSqrtWeights[ampName] = listNanMatrix 

503 dataset.aMatrix[ampName] = nanMatrix 

504 dataset.bMatrix[ampName] = nanMatrix 

505 dataset.covariancesModelNoB[ampName] = listNanMatrix 

506 dataset.aMatrixNoB[ampName] = nanMatrix 

507 

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

509 dataset.gain[ampName] = np.nan 

510 dataset.gainErr[ampName] = np.nan 

511 dataset.noise[ampName] = np.nan 

512 dataset.noiseErr[ampName] = np.nan 

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

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

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

516 continue 

517 

518 muAtAmp = dataset.rawMeans[ampName] 

519 maskAtAmp = dataset.expIdMask[ampName] 

520 if len(maskAtAmp) == 0: 

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

522 

523 muAtAmpMasked = muAtAmp[maskAtAmp] 

524 covAtAmp = dataset.covariances[ampName] 

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

526 covSqrtWeightsAtAmp = dataset.covariancesSqrtWeights[ampName] 

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

528 

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

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

531 muAtAmpMasked, 

532 covAtAmpMasked, 

533 covSqrtWeightsAtAmpMasked 

534 ) 

535 

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

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

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

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

540 'fullModelNoB': self.funcFullCovarianceModelNoB} 

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

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

543 for key in functionsDict: 

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

545 covAtAmpMasked.ravel(), functionsDict[key], 

546 weightsY=covSqrtWeightsAtAmpMasked.ravel()) 

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

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

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

550 gain = params[-1] 

551 

552 fitResults[key]['a'] = a 

553 fitResults[key]['c'] = c 

554 fitResults[key]['noise'] = noise 

555 fitResults[key]['gain'] = gain 

556 fitResults[key]['paramsErr'] = paramsErr 

557 

558 # Put the information in the PTC dataset 

559 

560 # Not used when ptcFitType is 'FULLCOVARIANCE' 

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

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

563 dataset.ptcFitChiSq[ampName] = np.nan 

564 

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

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

567 # converted to bool. 

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

569 dataset.covariances[ampName] = covAtAmp 

570 # We evaluate the covariance model everywhere, even the 

571 # masked amps. 

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

573 fitResults['fullModel']['a'], 

574 fitResults['fullModel']['c'], 

575 fitResults['fullModel']['noise'], 

576 fitResults['fullModel']['gain']) 

577 dataset.covariancesSqrtWeights[ampName] = covSqrtWeightsAtAmp 

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

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

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

581 fitResults['fullModelNoB']['a'], 

582 fitResults['fullModelNoB']['c'], 

583 fitResults['fullModelNoB']['noise'], 

584 fitResults['fullModelNoB']['gain'], 

585 setBtoZero=True) 

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

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

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

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

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

591 dataset.noise[ampName] = readoutNoise 

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

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

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

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

596 dataset.finalMeans[ampName] = muAtAmp 

597 

598 return dataset 

599 

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

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

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

603 of Astier+19. 

604 

605 Parameters 

606 ---------- 

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

608 Signal `mu` (ADU) 

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

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

611 `M = config.maximumRangeCovariancesAstier`), 

612 indexed by mean signal `mu`. 

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

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

615 

616 Returns 

617 ------- 

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

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

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

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

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

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

624 gain : `float` 

625 Amplifier gain (e/ADU) 

626 """ 

627 matrixSide = self.config.maximumRangeCovariancesAstier 

628 

629 # Initialize fit parameters 

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

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

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

633 gain = 1. 

634 

635 # iterate the fit to account for higher orders 

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

637 # stop when it increases 

638 oldChi2 = 1e30 

639 for _ in range(5): 

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

641 # loop on lags 

642 for i in range(matrixSide): 

643 for j in range(matrixSide): 

644 # fit a parabola for a given lag 

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

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

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

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

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

650 if(i + j == 0): 

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

652 weightedRes = (model - cov)*sqrtW 

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

654 if chi2 > oldChi2: 

655 break 

656 oldChi2 = chi2 

657 

658 return a, c, noise, gain 

659 

660 def funcFullCovarianceModel(self, params, x): 

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

662 Astier+19. 

663 

664 Parameters 

665 ---------- 

666 params : `list` 

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

668 gain (e/ADU). 

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

670 Signal `mu` (ADU) 

671 

672 Returns 

673 ------- 

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

675 Covariance matrix. 

676 """ 

677 matrixSide = self.config.maximumRangeCovariancesAstier 

678 lenParams = matrixSide*matrixSide 

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

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

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

682 gain = params[-1] 

683 

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

685 

686 def funcFullCovarianceModelNoB(self, params, x): 

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

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

689 

690 Parameters 

691 ---------- 

692 params : `list` 

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

694 gain (e/ADU). 

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

696 Signal mu (ADU) 

697 

698 Returns 

699 ------- 

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

701 Covariance matrix. 

702 """ 

703 matrixSide = self.config.maximumRangeCovariancesAstier 

704 lenParams = matrixSide*matrixSide 

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

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

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

708 gain = params[-1] 

709 

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

711 

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

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

714 

715 Parameters 

716 ---------- 

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

718 List of mean signals. 

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

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

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

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

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

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

725 gain : `float` 

726 Amplifier gain (e/ADU) 

727 setBtoZero=False : `bool`, optional 

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

729 

730 Returns 

731 ------- 

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

733 Covariances model. 

734 

735 Notes 

736 ----- 

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

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

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

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

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

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

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

744 

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

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

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

748 - gain, units: e/ADU 

749 

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

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

752 """ 

753 matrixSide = self.config.maximumRangeCovariancesAstier 

754 sa = (matrixSide, matrixSide) 

755 # pad a with zeros and symmetrize 

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

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

758 aSym = symmetrize(aEnlarged) 

759 # pad c with zeros and symmetrize 

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

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

762 cSym = symmetrize(cEnlarged) 

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

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

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

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

767 

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

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

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

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

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

773 

774 # assumes that mu is 1d 

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

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

777 # term, that is absent for now. 

778 if setBtoZero: 

779 c1 = np.zeros_like(c1) 

780 ac = np.zeros_like(ac) 

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

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

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

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

785 

786 return covModel 

787 

788 # EXPAPPROXIMATION and POLYNOMIAL fit methods 

789 @staticmethod 

790 def _initialParsForPolynomial(order): 

791 assert(order >= 2) 

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

793 pars[0] = 10 

794 pars[1] = 1 

795 pars[2:] = 0.0001 

796 return pars 

797 

798 @staticmethod 

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

800 if not len(lowers): 

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

802 if not len(uppers): 

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

804 lowers[1] = 0 # no negative gains 

805 return (lowers, uppers) 

806 

807 @staticmethod 

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

809 if not len(lowers): 

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

811 if not len(uppers): 

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

813 return (lowers, uppers) 

814 

815 @staticmethod 

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

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

818 

819 Parameters 

820 ---------- 

821 means : `numpy.array` 

822 Input array with mean signal values. 

823 variances : `numpy.array` 

824 Input array with variances at each mean value. 

825 minVarPivotSearch : `float` 

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

827 of decreasing variance should be sought. 

828 consecutivePointsVarDecreases : `int` 

829 Required number of consecutive points/fluxes 

830 in the PTC where the variance 

831 decreases in order to find a first 

832 estimate of the PTC turn-off. 

833 

834 Returns 

835 ------ 

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

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

838 points. 

839 

840 Notes 

841 ----- 

842 Eliminate points beyond which the variance decreases. 

843 """ 

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

845 # Variances are sorted and should monotonically increase 

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

847 if len(pivotList) > 0: 

848 # For small values, sometimes the variance decreases slightly 

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

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

851 # Require that the varince decreases during 

852 # consecutivePointsVarDecreases 

853 # consecutive points. This will give a first 

854 # estimate of the PTC turn-off, which 

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

856 if len(pivotList) > 1: 

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

858 # each value in pivotList. The lambda function subtracts 

859 # each value from the index. 

860 # groupby groups elements by equal key value. 

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

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

863 # Form groups of consecute values from pivotList 

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

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

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

867 # Find the first group of consecutive numbers when 

868 # variance decreases. 

869 if len(group) >= consecutivePointsVarDecreases: 

870 pivotIndex = np.min(group) 

871 goodPoints[pivotIndex+1:] = False 

872 break 

873 

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

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

876 

877 return goodPoints 

878 

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

880 """""" 

881 array = np.array(array) 

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

883 if nBad == 0: 

884 return array 

885 

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

887 if len(index): 

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

889 self.log.warning(msg) 

890 

891 array[index] = substituteValue 

892 

893 return array 

894 

895 def fitPtc(self, dataset): 

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

897 Astier+19 approximation (Eq. 16). 

898 

899 Fit the photon transfer curve with either a polynomial of 

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

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

902 

903 Sigma clipping is performed iteratively for the fit, as 

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

905 than `config.initialNonLinearityExclusionThreshold` away 

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

907 because the photon transfer curve turns over catastrophically 

908 at very high flux (because saturation 

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

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

911 to perform the sigma-clipping. 

912 

913 Parameters 

914 ---------- 

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

916 The dataset containing the means, variances and 

917 exposure times. 

918 

919 Returns 

920 ------- 

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

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

923 it has been modified to include information such as the 

924 fit vectors and the fit parameters. See the class 

925 `PhotonTransferCurveDatase`. 

926 

927 Raises 

928 ------ 

929 RuntimeError 

930 Raised if dataset.ptcFitType is None or empty. 

931 """ 

932 if dataset.ptcFitType: 

933 ptcFitType = dataset.ptcFitType 

934 else: 

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

936 matrixSide = self.config.maximumRangeCovariancesAstier 

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

938 nanMatrix[:] = np.nan 

939 

940 for amp in dataset.ampNames: 

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

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

943 listNanMatrix[:] = np.nan 

944 

945 dataset.covariancesModel[amp] = listNanMatrix 

946 dataset.aMatrix[amp] = nanMatrix 

947 dataset.bMatrix[amp] = nanMatrix 

948 dataset.covariancesModelNoB[amp] = listNanMatrix 

949 dataset.aMatrixNoB[amp] = nanMatrix 

950 

951 def errFunc(p, x, y): 

952 return ptcFunc(p, x) - y 

953 

954 sigmaCutPtcOutliers = self.config.sigmaCutPtcOutliers 

955 maxIterationsPtcOutliers = self.config.maxIterationsPtcOutliers 

956 

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

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

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

960 varVecOriginal = self._makeZeroSafe(varVecOriginal) 

961 

962 if self.config.doLegacyTurnoffSelection: 

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

964 # consecutive signal levels 

965 goodPoints = self._getInitialGoodPoints(meanVecOriginal, varVecOriginal, 

966 self.config.minVarPivotSearch, 

967 self.config.consecutivePointsVarDecreases) 

968 else: 

969 goodPoints = dataset.expIdMask[ampName] 

970 

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

972 initialExpIdMask = dataset.expIdMask[ampName] 

973 

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

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

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

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

978 self.log.warning(msg) 

979 # Fill entries with NaNs 

980 self.fillBadAmp(dataset, ptcFitType, ampName) 

981 continue 

982 

983 mask = goodPoints 

984 

985 if ptcFitType == 'EXPAPPROXIMATION': 

986 ptcFunc = funcAstier 

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

988 # lowers and uppers obtained from BOT data studies by 

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

990 if self.config.binSize > 1: 

991 bounds = self._boundsForAstier(parsIniPtc) 

992 else: 

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

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

995 if ptcFitType == 'POLYNOMIAL': 

996 ptcFunc = funcPolynomial 

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

998 bounds = self._boundsForPolynomial(parsIniPtc) 

999 

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

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

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

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

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

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

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

1007 # This algorithm was initially developed by Seth Digel for 

1008 # the EO Testing pipeline. 

1009 

1010 if maxIterationsPtcOutliers == 0: 

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

1012 # an initial fit. 

1013 res = least_squares( 

1014 errFunc, 

1015 parsIniPtc, 

1016 bounds=bounds, 

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

1018 ) 

1019 pars = res.x 

1020 newMask = mask.copy() 

1021 else: 

1022 newMask = (mask & (meanVecOriginal <= self.config.maxSignalInitialPtcOutlierFit)) 

1023 

1024 count = 0 

1025 lastMask = mask.copy() 

1026 while count < maxIterationsPtcOutliers: 

1027 res = least_squares( 

1028 errFunc, 

1029 parsIniPtc, 

1030 bounds=bounds, 

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

1032 ) 

1033 pars = res.x 

1034 

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

1036 # The new mask includes points where the residuals are 

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

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

1039 newMask = ( 

1040 np.isfinite(sigResids) 

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

1042 & mask 

1043 ) 

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

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

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

1047 self.log.warning(msg) 

1048 # Fill entries with NaNs 

1049 self.fillBadAmp(dataset, ptcFitType, ampName) 

1050 break 

1051 

1052 self.log.debug( 

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

1054 count, 

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

1056 ampName, 

1057 ) 

1058 

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

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

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

1062 break 

1063 

1064 lastMask = newMask.copy() 

1065 

1066 count += 1 

1067 

1068 # Set the mask to the new mask 

1069 mask = newMask.copy() 

1070 

1071 if not mask.any(): 

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

1073 continue 

1074 

1075 dataset.expIdMask[ampName] = mask 

1076 

1077 parsIniPtc = pars 

1078 meanVecFinal = meanVecOriginal[mask] 

1079 varVecFinal = varVecOriginal[mask] 

1080 

1081 # Save the maximum point after outlier detection as the 

1082 # PTC turnoff point. 

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

1084 

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

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

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

1088 

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

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

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

1092 self.log.warning(msg) 

1093 # Fill entries with NaNs 

1094 self.fillBadAmp(dataset, ptcFitType, ampName) 

1095 continue 

1096 # Fit the PTC. 

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

1098 # already calculated in `makeCovArray` of CpPtcExtract. 

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

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

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

1102 if self.config.doFitBootstrap: 

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

1104 varVecFinal, ptcFunc, 

1105 weightsY=weightsY) 

1106 else: 

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

1108 varVecFinal, ptcFunc, 

1109 weightsY=weightsY) 

1110 dataset.ptcFitPars[ampName] = parsFit 

1111 dataset.ptcFitParsError[ampName] = parsFitErr 

1112 dataset.ptcFitChiSq[ampName] = reducedChiSqPtc 

1113 

1114 dataset.finalVars[ampName] = varVecOriginal 

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

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

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

1118 dataset.finalMeans[ampName] = meanVecOriginal 

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

1120 

1121 if ptcFitType == 'EXPAPPROXIMATION': 

1122 ptcGain = parsFit[1] 

1123 ptcGainErr = parsFitErr[1] 

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

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

1126 if ptcFitType == 'POLYNOMIAL': 

1127 ptcGain = 1./parsFit[1] 

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

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

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

1131 dataset.gain[ampName] = ptcGain 

1132 dataset.gainErr[ampName] = ptcGainErr 

1133 dataset.noise[ampName] = ptcNoise 

1134 dataset.noiseErr[ampName] = ptcNoiseErr 

1135 

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

1137 dataset.ptcFitType = ptcFitType 

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

1139 dataset.badAmps = [] 

1140 

1141 return dataset 

1142 

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

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

1145 good points. 

1146 

1147 Parameters 

1148 ---------- 

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

1150 The dataset containing the means, variances and 

1151 exposure times. 

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

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

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

1155 ampName : `str` 

1156 Amplifier name. 

1157 """ 

1158 dataset.badAmps.append(ampName) 

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

1160 dataset.gain[ampName] = np.nan 

1161 dataset.gainErr[ampName] = np.nan 

1162 dataset.noise[ampName] = np.nan 

1163 dataset.noiseErr[ampName] = np.nan 

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

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

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

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

1168 dataset.ptcFitChiSq[ampName] = np.nan 

1169 dataset.ptcTurnoff[ampName] = np.nan 

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

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

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

1173 

1174 return