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

Shortcuts on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

290 statements  

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) 

28 

29from scipy.optimize import least_squares 

30 

31import lsst.pipe.base.connectionTypes as cT 

32 

33from .astierCovPtcUtils import fitDataFullCovariance 

34 

35from lsst.ip.isr import PhotonTransferCurveDataset 

36 

37from lsst.cp.pipe._lookupStaticCalibration import lookupStaticCalibration 

38 

39import copy 

40 

41 

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

43 

44 

45class PhotonTransferCurveSolveConnections(pipeBase.PipelineTaskConnections, 

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

47 inputCovariances = cT.Input( 

48 name="ptcCovariances", 

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

50 storageClass="PhotonTransferCurveDataset", 

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

52 multiple=True, 

53 ) 

54 camera = cT.PrerequisiteInput( 

55 name="camera", 

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

57 storageClass="Camera", 

58 dimensions=("instrument",), 

59 isCalibration=True, 

60 lookupFunction=lookupStaticCalibration, 

61 ) 

62 outputPtcDataset = cT.Output( 

63 name="ptcDatsetProposal", 

64 doc="Output proposed ptc dataset.", 

65 storageClass="PhotonTransferCurveDataset", 

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

67 multiple=False, 

68 isCalibration=True, 

69 ) 

70 

71 

72class PhotonTransferCurveSolveConfig(pipeBase.PipelineTaskConfig, 

73 pipelineConnections=PhotonTransferCurveSolveConnections): 

74 """Configuration for fitting measured covariances. 

75 """ 

76 

77 ptcFitType = pexConfig.ChoiceField( 

78 dtype=str, 

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

80 default="POLYNOMIAL", 

81 allowed={ 

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

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

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

85 } 

86 ) 

87 maximumRangeCovariancesAstier = pexConfig.Field( 

88 dtype=int, 

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

90 default=8, 

91 ) 

92 sigmaClipFullFitCovariancesAstier = pexConfig.Field( 

93 dtype=float, 

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

95 default=5.0, 

96 ) 

97 maxIterFullFitCovariancesAstier = pexConfig.Field( 

98 dtype=int, 

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

100 default=3, 

101 ) 

102 polynomialFitDegree = pexConfig.Field( 

103 dtype=int, 

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

105 default=3, 

106 ) 

107 sigmaCutPtcOutliers = pexConfig.Field( 

108 dtype=float, 

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

110 default=5.0, 

111 ) 

112 maxIterationsPtcOutliers = pexConfig.Field( 

113 dtype=int, 

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

115 default=2, 

116 ) 

117 minVarPivotSearch = pexConfig.Field( 

118 dtype=float, 

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

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

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

122 " should be sought.", 

123 default=10000, 

124 ) 

125 doFitBootstrap = pexConfig.Field( 

126 dtype=bool, 

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

128 default=False, 

129 ) 

130 

131 

132class PhotonTransferCurveSolveTask(pipeBase.PipelineTask, 

133 pipeBase.CmdLineTask): 

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

135 

136 This task assembles the list of individual PTC datasets produced 

137 by ``PhotonTransferCurveSolveTask`` into one single final PTC 

138 dataset. The task fits the measured (co)variances to a polynomial 

139 model or to the models described in equations 16 and 20 of 

140 Astier+19 (referred to as ``POLYNOMIAL``, ``EXPAPPROXIMATION``, 

141 and ``FULLCOVARIANCE`` in the configuration options of the task, 

142 respectively). Parameters of interest such as tghe gain and noise 

143 are derived from the fits. 

144 

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

146 of CCD sensors", arXiv:1905.08677 

147 """ 

148 

149 ConfigClass = PhotonTransferCurveSolveConfig 

150 _DefaultName = 'cpPhotonTransferCurveSolve' 

151 

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

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

154 

155 Parameters 

156 ---------- 

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

158 Butler to operate on. 

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

160 Input data refs to load. 

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

162 Output data refs to persist. 

163 """ 

164 inputs = butlerQC.get(inputRefs) 

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

166 butlerQC.put(outputs, outputRefs) 

167 

168 def run(self, inputCovariances, camera=None, inputExpList=None): 

169 """Fit measure covariances to different models. 

170 

171 Parameters 

172 ---------- 

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

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

175 

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

177 Input camera. 

178 

179 inputExpList : `list` [`~lsst.afw.image.ExposureF`], optional 

180 List of exposures. 

181 

182 Returns 

183 ------- 

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

185 The results struct containing: 

186 

187 ``outputPtcDatset`` 

188 Final PTC dataset, containing information such as the 

189 means, variances, and exposure times 

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

191 """ 

192 # Assemble partial PTC datasets into a single dataset. 

193 ampNames = np.unique(inputCovariances[0].ampNames) 

194 datasetPtc = PhotonTransferCurveDataset(ampNames, self.config.ptcFitType, 

195 self.config.maximumRangeCovariancesAstier) 

196 for partialPtcDataset in inputCovariances: 

197 if partialPtcDataset.ptcFitType == 'DUMMY': 

198 continue 

199 for ampName in ampNames: 

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

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

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

203 else: 

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

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

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

207 else: 

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

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

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

211 else: 

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

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

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

215 else: 

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

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

218 datasetPtc.covariancesSqrtWeights[ampName].append( 

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

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

221 # rawMeans index 

222 for ampName in ampNames: 

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

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

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

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

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

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

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

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

231 datasetPtc.covariancesSqrtWeights[ampName])[index] 

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

233 # Calculate covariances and fit them, including the PTC, 

234 # to Astier+19 full model (Eq. 20) First, fit get the flat 

235 # pairs that are masked, fitting C_00 vs mu to the 

236 # EXPAPPROXIMATION model (Eq. 16 in Astier+19). The 

237 # points at these fluxes will also be masked when 

238 # calculating the other covariances, C_ij) 

239 tempDatasetPtc = copy.copy(datasetPtc) 

240 tempDatasetPtc.ptcFitType = "EXPAPPROXIMATION" 

241 tempDatasetPtc = self.fitPtc(tempDatasetPtc) 

242 for ampName in datasetPtc.ampNames: 

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

244 datasetPtc.fitType = "FULLCOVARIANCE" 

245 datasetPtc = self.fitCovariancesAstier(datasetPtc) 

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

247 # ("EXPAPPROXIMATION", "POLYNOMIAL") 

248 else: 

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

250 # approximation (Eq. 16). Fill up 

251 # PhotonTransferCurveDataset object. 

252 datasetPtc = self.fitPtc(datasetPtc) 

253 if inputExpList is not None: 

254 # It should be a list of exposures, to get the detector. 

255 detector = inputExpList[0].getDetector() 

256 else: 

257 detector = None 

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

259 

260 return pipeBase.Struct( 

261 outputPtcDataset=datasetPtc, 

262 ) 

263 

264 def fitCovariancesAstier(self, dataset): 

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

266 

267 Parameters 

268 ---------- 

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

270 The dataset containing information such as the means, 

271 (co)variances, and exposure times. 

272 

273 Returns 

274 ------- 

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

276 This is the same dataset as the input paramter, however, 

277 it has been modified to include information such as the 

278 fit vectors and the fit parameters. See the class 

279 `PhotonTransferCurveDatase`. 

280 """ 

281 covFits, covFitsNoB = fitDataFullCovariance(dataset) 

282 dataset = self.getOutputPtcDataCovAstier(dataset, covFits, covFitsNoB) 

283 

284 return dataset 

285 

286 def getOutputPtcDataCovAstier(self, dataset, covFits, covFitsNoB): 

287 """Get output data for PhotonTransferCurveCovAstierDataset from CovFit 

288 objects. 

289 

290 Parameters 

291 ---------- 

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

293 The dataset containing information such as the means, 

294 variances and exposure times. 

295 covFits : `dict` 

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

297 covFitsNoB : `dict` 

298 Dictionary of CovFit objects, with amp names as keys, and 

299 'b=0' in Eq. 20 of Astier+19. 

300 

301 Returns 

302 ------- 

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

304 This is the same dataset as the input paramter, however, 

305 it has been modified to include extra information such as 

306 the mask 1D array, gains, reoudout noise, measured signal, 

307 measured variance, modeled variance, a, and b coefficient 

308 matrices (see Astier+19) per amplifier. See the class 

309 `PhotonTransferCurveDatase`. 

310 """ 

311 assert(len(covFits) == len(covFitsNoB)) 

312 

313 for i, amp in enumerate(dataset.ampNames): 

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

315 # Not used when ptcFitType is 'FULLCOVARIANCE' 

316 dataset.ptcFitPars[amp] = [np.nan] 

317 dataset.ptcFitParsError[amp] = [np.nan] 

318 dataset.ptcFitChiSq[amp] = np.nan 

319 if amp in covFits: 

320 fit = covFits[amp] 

321 fitNoB = covFitsNoB[amp] 

322 # Save full covariances, covariances models, and their weights 

323 # dataset.expIdMask is already full 

324 dataset.covariances[amp] = fit.cov 

325 dataset.covariancesModel[amp] = fit.evalCovModel() 

326 dataset.covariancesSqrtWeights[amp] = fit.sqrtW 

327 dataset.aMatrix[amp] = fit.getA() 

328 dataset.bMatrix[amp] = fit.getB() 

329 dataset.covariancesModelNoB[amp] = fitNoB.evalCovModel() 

330 dataset.aMatrixNoB[amp] = fitNoB.getA() 

331 

332 (meanVecFinal, varVecFinal, varVecModel, 

333 wc, varMask) = fit.getFitData(0, 0, divideByMu=False) 

334 gain = fit.getGain() 

335 

336 dataset.gain[amp] = gain 

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

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

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

340 dataset.finalVars[amp] = varVecFinal 

341 dataset.finalModelVars[amp] = varVecModel 

342 dataset.finalMeans[amp] = meanVecFinal 

343 

344 else: 

345 # Bad amp 

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

347 # with astropy.Table works. 

348 matrixSide = self.config.maximumRangeCovariancesAstier 

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

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

351 

352 dataset.covariances[amp] = listNanMatrix 

353 dataset.covariancesModel[amp] = listNanMatrix 

354 dataset.covariancesSqrtWeights[amp] = listNanMatrix 

355 dataset.aMatrix[amp] = nanMatrix 

356 dataset.bMatrix[amp] = nanMatrix 

357 dataset.covariancesModelNoB[amp] = listNanMatrix 

358 dataset.aMatrixNoB[amp] = nanMatrix 

359 

360 dataset.expIdMask[amp] = np.repeat(np.nan, lenInputTimes) 

361 dataset.gain[amp] = np.nan 

362 dataset.gainErr[amp] = np.nan 

363 dataset.noise[amp] = np.nan 

364 dataset.noiseErr[amp] = np.nan 

365 dataset.finalVars[amp] = np.repeat(np.nan, lenInputTimes) 

366 dataset.finalModelVars[amp] = np.repeat(np.nan, lenInputTimes) 

367 dataset.finalMeans[amp] = np.repeat(np.nan, lenInputTimes) 

368 

369 return dataset 

370 

371 @staticmethod 

372 def _initialParsForPolynomial(order): 

373 assert(order >= 2) 

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

375 pars[0] = 10 

376 pars[1] = 1 

377 pars[2:] = 0.0001 

378 return pars 

379 

380 @staticmethod 

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

382 if not len(lowers): 

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

384 if not len(uppers): 

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

386 lowers[1] = 0 # no negative gains 

387 return (lowers, uppers) 

388 

389 @staticmethod 

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

391 if not len(lowers): 

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

393 if not len(uppers): 

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

395 return (lowers, uppers) 

396 

397 @staticmethod 

398 def _getInitialGoodPoints(means, variances, minVarPivotSearch): 

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

400 

401 Parameters 

402 ---------- 

403 means : `numpy.array` 

404 Input array with mean signal values. 

405 variances : `numpy.array` 

406 Input array with variances at each mean value. 

407 minVarPivotSearch : `float` 

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

409 of decreasing variance should be sought. 

410 

411 Returns 

412 ------ 

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

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

415 points. 

416 

417 Notes 

418 ----- 

419 Eliminate points beyond which the variance decreases 

420 """ 

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

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

423 if len(pivotList) > 0: 

424 # For small values, sometimes the variance decreases slightly 

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

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

427 if len(pivotList) > 1: 

428 # Require that the decrease in variance happen for two 

429 # consecutive signal levels 

430 pivotIndex = np.min(np.where(np.diff(pivotList) == 1)[0]) 

431 pivot = pivotList[pivotIndex] 

432 goodPoints[pivot+1:] = False 

433 

434 return goodPoints 

435 

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

437 """""" 

438 array = np.array(array) 

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

440 if nBad == 0: 

441 return array 

442 

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

444 if len(index): 

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

446 self.log.warning(msg) 

447 

448 array[index] = substituteValue 

449 

450 return array 

451 

452 def fitPtc(self, dataset): 

453 """Fit the photon transfer curve to a polynomial or to Astier+19 

454 approximation. 

455 

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

457 specified in the task config, or using the exponential approximation 

458 in Astier+19 (Eq. 16). 

459 

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

461 initial clipping of data points that are more than 

462 config.initialNonLinearityExclusionThreshold away from lying on a 

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

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

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

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

467 sigma-clipping. 

468 

469 Parameters 

470 ---------- 

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

472 The dataset containing the means, variances and exposure times. 

473 

474 Returns 

475 ------- 

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

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

478 it has been modified to include information such as the 

479 fit vectors and the fit parameters. See the class 

480 `PhotonTransferCurveDatase`. 

481 

482 Raises 

483 ------ 

484 RuntimeError: 

485 Raises if dataset.ptcFitType is None or empty. 

486 """ 

487 if dataset.ptcFitType: 

488 ptcFitType = dataset.ptcFitType 

489 else: 

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

491 matrixSide = self.config.maximumRangeCovariancesAstier 

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

493 nanMatrix[:] = np.nan 

494 

495 for amp in dataset.ampNames: 

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

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

498 listNanMatrix[:] = np.nan 

499 

500 dataset.covariancesModel[amp] = listNanMatrix 

501 dataset.aMatrix[amp] = nanMatrix 

502 dataset.bMatrix[amp] = nanMatrix 

503 dataset.covariancesModelNoB[amp] = listNanMatrix 

504 dataset.aMatrixNoB[amp] = nanMatrix 

505 

506 def errFunc(p, x, y): 

507 return ptcFunc(p, x) - y 

508 

509 sigmaCutPtcOutliers = self.config.sigmaCutPtcOutliers 

510 maxIterationsPtcOutliers = self.config.maxIterationsPtcOutliers 

511 

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

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

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

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

516 varVecOriginal = self._makeZeroSafe(varVecOriginal) 

517 

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

519 # consecutive signal levels 

520 goodPoints = self._getInitialGoodPoints(meanVecOriginal, varVecOriginal, 

521 self.config.minVarPivotSearch) 

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

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

524 

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

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

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

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

529 self.log.warning(msg) 

530 # Fill entries with NaNs 

531 self.fillBadAmp(dataset, ptcFitType, ampName) 

532 continue 

533 

534 mask = goodPoints 

535 

536 if ptcFitType == 'EXPAPPROXIMATION': 

537 ptcFunc = funcAstier 

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

539 # lowers and uppers obtained from BOT data studies by 

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

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

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

543 if ptcFitType == 'POLYNOMIAL': 

544 ptcFunc = funcPolynomial 

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

546 bounds = self._boundsForPolynomial(parsIniPtc) 

547 

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

549 count = 1 

550 while count <= maxIterationsPtcOutliers: 

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

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

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

554 meanTempVec = meanVecOriginal[mask] 

555 varTempVec = varVecOriginal[mask] 

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

557 pars = res.x 

558 

559 # change this to the original from the temp because 

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

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

562 # same length for broadcasting 

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

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

565 mask = mask & newMask 

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

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

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

569 self.log.warning(msg) 

570 # Fill entries with NaNs 

571 self.fillBadAmp(dataset, ptcFitType, ampName) 

572 break 

573 nDroppedTotal = Counter(mask)[False] 

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

575 count, nDroppedTotal, ampName) 

576 count += 1 

577 # objects should never shrink 

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

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

580 continue 

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

582 # store the final mask 

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

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

585 else: 

586 dataset.expIdMask[ampName] = mask 

587 parsIniPtc = pars 

588 meanVecFinal = meanVecOriginal[mask] 

589 varVecFinal = varVecOriginal[mask] 

590 

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

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

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

594 

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

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

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

598 self.log.warning(msg) 

599 # Fill entries with NaNs 

600 self.fillBadAmp(dataset, ptcFitType, ampName) 

601 continue 

602 # Fit the PTC 

603 if self.config.doFitBootstrap: 

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

605 varVecFinal, ptcFunc, 

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

607 else: 

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

609 varVecFinal, ptcFunc, 

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

611 dataset.ptcFitPars[ampName] = parsFit 

612 dataset.ptcFitParsError[ampName] = parsFitErr 

613 dataset.ptcFitChiSq[ampName] = reducedChiSqPtc 

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

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

616 # mask may vary per amp). 

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

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

619 constant_values=np.nan) 

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

621 'constant', constant_values=np.nan) 

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

623 constant_values=np.nan) 

624 if ptcFitType == 'EXPAPPROXIMATION': 

625 ptcGain = parsFit[1] 

626 ptcGainErr = parsFitErr[1] 

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

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

629 if ptcFitType == 'POLYNOMIAL': 

630 ptcGain = 1./parsFit[1] 

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

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

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

634 dataset.gain[ampName] = ptcGain 

635 dataset.gainErr[ampName] = ptcGainErr 

636 dataset.noise[ampName] = ptcNoise 

637 dataset.noiseErr[ampName] = ptcNoiseErr 

638 

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

640 dataset.ptcFitType = ptcFitType 

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

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

643 

644 return dataset 

645 

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

647 """Fill the dataset with NaNs if there are not enough good points. 

648 

649 Parameters 

650 ---------- 

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

652 The dataset containing the means, variances and exposure times. 

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

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

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

656 ampName : `str` 

657 Amplifier name. 

658 """ 

659 dataset.badAmps.append(ampName) 

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

661 dataset.gain[ampName] = np.nan 

662 dataset.gainErr[ampName] = np.nan 

663 dataset.noise[ampName] = np.nan 

664 dataset.noiseErr[ampName] = np.nan 

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

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

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

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

669 dataset.ptcFitChiSq[ampName] = np.nan 

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

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

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

673 

674 return