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

445 statements  

« prev     ^ index     » next       coverage.py v7.2.7, created at 2023-06-25 03:07 -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 dataset.noiseMatrix[ampName] = nanMatrix 

508 dataset.noiseMatrixNoB[ampName] = nanMatrix 

509 

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

511 dataset.gain[ampName] = np.nan 

512 dataset.gainErr[ampName] = np.nan 

513 dataset.noise[ampName] = np.nan 

514 dataset.noiseErr[ampName] = np.nan 

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

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

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

518 continue 

519 

520 muAtAmp = dataset.rawMeans[ampName] 

521 maskAtAmp = dataset.expIdMask[ampName] 

522 if len(maskAtAmp) == 0: 

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

524 

525 muAtAmpMasked = muAtAmp[maskAtAmp] 

526 covAtAmp = dataset.covariances[ampName] 

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

528 covSqrtWeightsAtAmp = dataset.covariancesSqrtWeights[ampName] 

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

530 

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

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

533 muAtAmpMasked, 

534 covAtAmpMasked, 

535 covSqrtWeightsAtAmpMasked 

536 ) 

537 

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

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

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

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

542 'fullModelNoB': self.funcFullCovarianceModelNoB} 

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

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

545 for key in functionsDict: 

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

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

548 weightsY=covSqrtWeightsAtAmpMasked.ravel()) 

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

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

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

552 gain = params[-1] 

553 

554 fitResults[key]['a'] = a 

555 fitResults[key]['c'] = c 

556 fitResults[key]['noise'] = noise 

557 fitResults[key]['gain'] = gain 

558 fitResults[key]['paramsErr'] = paramsErr 

559 

560 # Put the information in the PTC dataset 

561 

562 # Not used when ptcFitType is 'FULLCOVARIANCE' 

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

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

565 dataset.ptcFitChiSq[ampName] = np.nan 

566 

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

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

569 # converted to bool. 

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

571 dataset.covariances[ampName] = covAtAmp 

572 # We evaluate the covariance model everywhere, even the 

573 # masked amps. 

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

575 fitResults['fullModel']['a'], 

576 fitResults['fullModel']['c'], 

577 fitResults['fullModel']['noise'], 

578 fitResults['fullModel']['gain']) 

579 dataset.covariancesSqrtWeights[ampName] = covSqrtWeightsAtAmp 

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

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

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

583 fitResults['fullModelNoB']['a'], 

584 fitResults['fullModelNoB']['c'], 

585 fitResults['fullModelNoB']['noise'], 

586 fitResults['fullModelNoB']['gain'], 

587 setBtoZero=True) 

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

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

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

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

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

593 dataset.noise[ampName] = readoutNoise 

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

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

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

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

598 

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

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

601 dataset.finalMeans[ampName] = muAtAmp 

602 

603 return dataset 

604 

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

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

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

608 of Astier+19. 

609 

610 Parameters 

611 ---------- 

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

613 Signal `mu` (ADU) 

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

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

616 `M = config.maximumRangeCovariancesAstier`), 

617 indexed by mean signal `mu`. 

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

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

620 

621 Returns 

622 ------- 

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

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

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

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

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

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

629 gain : `float` 

630 Amplifier gain (e/ADU) 

631 """ 

632 matrixSide = self.config.maximumRangeCovariancesAstier 

633 

634 # Initialize fit parameters 

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

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

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

638 gain = 1. 

639 

640 # iterate the fit to account for higher orders 

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

642 # stop when it increases 

643 oldChi2 = 1e30 

644 for _ in range(5): 

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

646 # loop on lags 

647 for i in range(matrixSide): 

648 for j in range(matrixSide): 

649 # fit a parabola for a given lag 

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

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

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

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

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

655 if(i + j == 0): 

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

657 weightedRes = (model - cov)*sqrtW 

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

659 if chi2 > oldChi2: 

660 break 

661 oldChi2 = chi2 

662 

663 return a, c, noise, gain 

664 

665 def funcFullCovarianceModel(self, params, x): 

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

667 Astier+19. 

668 

669 Parameters 

670 ---------- 

671 params : `list` 

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

673 gain (e/ADU). 

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

675 Signal `mu` (ADU) 

676 

677 Returns 

678 ------- 

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

680 Covariance matrix. 

681 """ 

682 matrixSide = self.config.maximumRangeCovariancesAstier 

683 lenParams = matrixSide*matrixSide 

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

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

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

687 gain = params[-1] 

688 

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

690 

691 def funcFullCovarianceModelNoB(self, params, x): 

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

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

694 

695 Parameters 

696 ---------- 

697 params : `list` 

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

699 gain (e/ADU). 

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

701 Signal mu (ADU) 

702 

703 Returns 

704 ------- 

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

706 Covariance matrix. 

707 """ 

708 matrixSide = self.config.maximumRangeCovariancesAstier 

709 lenParams = matrixSide*matrixSide 

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

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

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

713 gain = params[-1] 

714 

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

716 

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

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

719 

720 Parameters 

721 ---------- 

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

723 List of mean signals. 

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

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

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

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

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

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

730 gain : `float` 

731 Amplifier gain (e/ADU) 

732 setBtoZero=False : `bool`, optional 

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

734 

735 Returns 

736 ------- 

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

738 Covariances model. 

739 

740 Notes 

741 ----- 

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

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

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

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

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

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

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

749 

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

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

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

753 - gain, units: e/ADU 

754 

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

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

757 """ 

758 matrixSide = self.config.maximumRangeCovariancesAstier 

759 sa = (matrixSide, matrixSide) 

760 # pad a with zeros and symmetrize 

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

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

763 aSym = symmetrize(aEnlarged) 

764 # pad c with zeros and symmetrize 

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

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

767 cSym = symmetrize(cEnlarged) 

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

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

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

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

772 

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

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

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

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

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

778 

779 # assumes that mu is 1d 

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

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

782 # term, that is absent for now. 

783 if setBtoZero: 

784 c1 = np.zeros_like(c1) 

785 ac = np.zeros_like(ac) 

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

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

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

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

790 

791 return covModel 

792 

793 # EXPAPPROXIMATION and POLYNOMIAL fit methods 

794 @staticmethod 

795 def _initialParsForPolynomial(order): 

796 assert(order >= 2) 

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

798 pars[0] = 10 

799 pars[1] = 1 

800 pars[2:] = 0.0001 

801 return pars 

802 

803 @staticmethod 

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

805 if not len(lowers): 

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

807 if not len(uppers): 

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

809 lowers[1] = 0 # no negative gains 

810 return (lowers, uppers) 

811 

812 @staticmethod 

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

814 if not len(lowers): 

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

816 if not len(uppers): 

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

818 return (lowers, uppers) 

819 

820 @staticmethod 

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

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

823 

824 Parameters 

825 ---------- 

826 means : `numpy.array` 

827 Input array with mean signal values. 

828 variances : `numpy.array` 

829 Input array with variances at each mean value. 

830 minVarPivotSearch : `float` 

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

832 of decreasing variance should be sought. 

833 consecutivePointsVarDecreases : `int` 

834 Required number of consecutive points/fluxes 

835 in the PTC where the variance 

836 decreases in order to find a first 

837 estimate of the PTC turn-off. 

838 

839 Returns 

840 ------ 

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

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

843 points. 

844 

845 Notes 

846 ----- 

847 Eliminate points beyond which the variance decreases. 

848 """ 

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

850 # Variances are sorted and should monotonically increase 

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

852 if len(pivotList) > 0: 

853 # For small values, sometimes the variance decreases slightly 

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

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

856 # Require that the varince decreases during 

857 # consecutivePointsVarDecreases 

858 # consecutive points. This will give a first 

859 # estimate of the PTC turn-off, which 

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

861 if len(pivotList) > 1: 

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

863 # each value in pivotList. The lambda function subtracts 

864 # each value from the index. 

865 # groupby groups elements by equal key value. 

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

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

868 # Form groups of consecute values from pivotList 

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

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

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

872 # Find the first group of consecutive numbers when 

873 # variance decreases. 

874 if len(group) >= consecutivePointsVarDecreases: 

875 pivotIndex = np.min(group) 

876 goodPoints[pivotIndex+1:] = False 

877 break 

878 

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

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

881 

882 return goodPoints 

883 

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

885 """""" 

886 array = np.array(array) 

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

888 if nBad == 0: 

889 return array 

890 

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

892 if len(index): 

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

894 self.log.warning(msg) 

895 

896 array[index] = substituteValue 

897 

898 return array 

899 

900 def fitPtc(self, dataset): 

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

902 Astier+19 approximation (Eq. 16). 

903 

904 Fit the photon transfer curve with either a polynomial of 

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

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

907 

908 Sigma clipping is performed iteratively for the fit, as 

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

910 than `config.initialNonLinearityExclusionThreshold` away 

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

912 because the photon transfer curve turns over catastrophically 

913 at very high flux (because saturation 

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

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

916 to perform the sigma-clipping. 

917 

918 Parameters 

919 ---------- 

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

921 The dataset containing the means, variances and 

922 exposure times. 

923 

924 Returns 

925 ------- 

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

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

928 it has been modified to include information such as the 

929 fit vectors and the fit parameters. See the class 

930 `PhotonTransferCurveDatase`. 

931 

932 Raises 

933 ------ 

934 RuntimeError 

935 Raised if dataset.ptcFitType is None or empty. 

936 """ 

937 if dataset.ptcFitType: 

938 ptcFitType = dataset.ptcFitType 

939 else: 

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

941 matrixSide = self.config.maximumRangeCovariancesAstier 

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

943 nanMatrix[:] = np.nan 

944 

945 for amp in dataset.ampNames: 

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

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

948 listNanMatrix[:] = np.nan 

949 

950 dataset.covariancesModel[amp] = listNanMatrix 

951 dataset.aMatrix[amp] = nanMatrix 

952 dataset.bMatrix[amp] = nanMatrix 

953 dataset.covariancesModelNoB[amp] = listNanMatrix 

954 dataset.aMatrixNoB[amp] = nanMatrix 

955 dataset.noiseMatrix[amp] = nanMatrix 

956 dataset.noiseMatrixNoB[amp] = nanMatrix 

957 

958 def errFunc(p, x, y): 

959 return ptcFunc(p, x) - y 

960 

961 sigmaCutPtcOutliers = self.config.sigmaCutPtcOutliers 

962 maxIterationsPtcOutliers = self.config.maxIterationsPtcOutliers 

963 

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

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

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

967 varVecOriginal = self._makeZeroSafe(varVecOriginal) 

968 

969 if self.config.doLegacyTurnoffSelection: 

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

971 # consecutive signal levels 

972 goodPoints = self._getInitialGoodPoints(meanVecOriginal, varVecOriginal, 

973 self.config.minVarPivotSearch, 

974 self.config.consecutivePointsVarDecreases) 

975 else: 

976 goodPoints = dataset.expIdMask[ampName] 

977 

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

979 initialExpIdMask = dataset.expIdMask[ampName] 

980 

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

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

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

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

985 self.log.warning(msg) 

986 # Fill entries with NaNs 

987 self.fillBadAmp(dataset, ptcFitType, ampName) 

988 continue 

989 

990 mask = goodPoints 

991 

992 if ptcFitType == 'EXPAPPROXIMATION': 

993 ptcFunc = funcAstier 

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

995 # lowers and uppers obtained from BOT data studies by 

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

997 if self.config.binSize > 1: 

998 bounds = self._boundsForAstier(parsIniPtc) 

999 else: 

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

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

1002 if ptcFitType == 'POLYNOMIAL': 

1003 ptcFunc = funcPolynomial 

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

1005 bounds = self._boundsForPolynomial(parsIniPtc) 

1006 

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

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

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

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

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

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

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

1014 # This algorithm was initially developed by Seth Digel for 

1015 # the EO Testing pipeline. 

1016 

1017 if maxIterationsPtcOutliers == 0: 

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

1019 # an initial fit. 

1020 res = least_squares( 

1021 errFunc, 

1022 parsIniPtc, 

1023 bounds=bounds, 

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

1025 ) 

1026 pars = res.x 

1027 newMask = mask.copy() 

1028 else: 

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

1030 

1031 count = 0 

1032 lastMask = mask.copy() 

1033 while count < maxIterationsPtcOutliers: 

1034 res = least_squares( 

1035 errFunc, 

1036 parsIniPtc, 

1037 bounds=bounds, 

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

1039 ) 

1040 pars = res.x 

1041 

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

1043 # The new mask includes points where the residuals are 

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

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

1046 newMask = ( 

1047 np.isfinite(sigResids) 

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

1049 & mask 

1050 ) 

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

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

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

1054 self.log.warning(msg) 

1055 # Fill entries with NaNs 

1056 self.fillBadAmp(dataset, ptcFitType, ampName) 

1057 break 

1058 

1059 self.log.debug( 

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

1061 count, 

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

1063 ampName, 

1064 ) 

1065 

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

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

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

1069 break 

1070 

1071 lastMask = newMask.copy() 

1072 

1073 count += 1 

1074 

1075 # Set the mask to the new mask 

1076 mask = newMask.copy() 

1077 

1078 if not mask.any(): 

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

1080 continue 

1081 

1082 dataset.expIdMask[ampName] = mask 

1083 

1084 parsIniPtc = pars 

1085 meanVecFinal = meanVecOriginal[mask] 

1086 varVecFinal = varVecOriginal[mask] 

1087 

1088 # Save the maximum point after outlier detection as the 

1089 # PTC turnoff point. 

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

1091 

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

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

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

1095 

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

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

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

1099 self.log.warning(msg) 

1100 # Fill entries with NaNs 

1101 self.fillBadAmp(dataset, ptcFitType, ampName) 

1102 continue 

1103 # Fit the PTC. 

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

1105 # already calculated in `makeCovArray` of CpPtcExtract. 

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

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

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

1109 if self.config.doFitBootstrap: 

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

1111 varVecFinal, ptcFunc, 

1112 weightsY=weightsY) 

1113 else: 

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

1115 varVecFinal, ptcFunc, 

1116 weightsY=weightsY) 

1117 dataset.ptcFitPars[ampName] = parsFit 

1118 dataset.ptcFitParsError[ampName] = parsFitErr 

1119 dataset.ptcFitChiSq[ampName] = reducedChiSqPtc 

1120 

1121 dataset.finalVars[ampName] = varVecOriginal 

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

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

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

1125 dataset.finalMeans[ampName] = meanVecOriginal 

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

1127 

1128 if ptcFitType == 'EXPAPPROXIMATION': 

1129 ptcGain = parsFit[1] 

1130 ptcGainErr = parsFitErr[1] 

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

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

1133 if ptcFitType == 'POLYNOMIAL': 

1134 ptcGain = 1./parsFit[1] 

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

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

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

1138 dataset.gain[ampName] = ptcGain 

1139 dataset.gainErr[ampName] = ptcGainErr 

1140 dataset.noise[ampName] = ptcNoise 

1141 dataset.noiseErr[ampName] = ptcNoiseErr 

1142 

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

1144 dataset.ptcFitType = ptcFitType 

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

1146 dataset.badAmps = [] 

1147 

1148 return dataset 

1149 

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

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

1152 good points. 

1153 

1154 Parameters 

1155 ---------- 

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

1157 The dataset containing the means, variances and 

1158 exposure times. 

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

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

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

1162 ampName : `str` 

1163 Amplifier name. 

1164 """ 

1165 dataset.badAmps.append(ampName) 

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

1167 dataset.gain[ampName] = np.nan 

1168 dataset.gainErr[ampName] = np.nan 

1169 dataset.noise[ampName] = np.nan 

1170 dataset.noiseErr[ampName] = np.nan 

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

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

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

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

1175 dataset.ptcFitChiSq[ampName] = np.nan 

1176 dataset.ptcTurnoff[ampName] = np.nan 

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

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

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

1180 

1181 return