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

245 statements  

« prev     ^ index     » next       coverage.py v7.3.0, created at 2023-08-29 08:58 +0000

1# This file is part of cp_pipe. 

2# 

3# Developed for the LSST Data Management System. 

4# This product includes software developed by the LSST Project 

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

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

7# for details of code ownership. 

8# 

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

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

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

12# (at your option) any later version. 

13# 

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

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

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

17# GNU General Public License for more details. 

18# 

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

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

21# 

22 

23__all__ = ["LinearitySolveTask", "LinearitySolveConfig"] 

24 

25import numpy as np 

26import lsst.afw.image as afwImage 

27import lsst.pipe.base as pipeBase 

28import lsst.pipe.base.connectionTypes as cT 

29import lsst.pex.config as pexConfig 

30 

31from lsstDebug import getDebugFrame 

32from lsst.ip.isr import (Linearizer, IsrProvenance) 

33 

34from .utils import funcPolynomial, irlsFit, AstierSplineLinearityFitter 

35from ._lookupStaticCalibration import lookupStaticCalibration 

36 

37 

38def ptcLookup(datasetType, registry, quantumDataId, collections): 

39 """Butler lookup function to allow PTC to be found. 

40 

41 Parameters 

42 ---------- 

43 datasetType : `lsst.daf.butler.DatasetType` 

44 Dataset type to look up. 

45 registry : `lsst.daf.butler.Registry` 

46 Registry for the data repository being searched. 

47 quantumDataId : `lsst.daf.butler.DataCoordinate` 

48 Data ID for the quantum of the task this dataset will be passed to. 

49 This must include an "instrument" key, and should also include any 

50 keys that are present in ``datasetType.dimensions``. If it has an 

51 ``exposure`` or ``visit`` key, that's a sign that this function is 

52 not actually needed, as those come with the temporal information that 

53 would allow a real validity-range lookup. 

54 collections : `lsst.daf.butler.registry.CollectionSearch` 

55 Collections passed by the user when generating a QuantumGraph. Ignored 

56 by this function (see notes below). 

57 

58 Returns 

59 ------- 

60 refs : `list` [ `DatasetRef` ] 

61 A zero- or single-element list containing the matching 

62 dataset, if one was found. 

63 

64 Raises 

65 ------ 

66 RuntimeError 

67 Raised if more than one PTC reference is found. 

68 """ 

69 refs = list(registry.queryDatasets(datasetType, dataId=quantumDataId, collections=collections, 

70 findFirst=False)) 

71 if len(refs) >= 2: 

72 RuntimeError("Too many PTC connections found. Incorrect collections supplied?") 

73 

74 return refs 

75 

76 

77class LinearitySolveConnections(pipeBase.PipelineTaskConnections, 

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

79 dummy = cT.Input( 

80 name="raw", 

81 doc="Dummy exposure.", 

82 storageClass='Exposure', 

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

84 multiple=True, 

85 deferLoad=True, 

86 ) 

87 

88 camera = cT.PrerequisiteInput( 

89 name="camera", 

90 doc="Camera Geometry definition.", 

91 storageClass="Camera", 

92 dimensions=("instrument", ), 

93 isCalibration=True, 

94 lookupFunction=lookupStaticCalibration, 

95 ) 

96 

97 inputPtc = cT.PrerequisiteInput( 

98 name="ptc", 

99 doc="Input PTC dataset.", 

100 storageClass="PhotonTransferCurveDataset", 

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

102 isCalibration=True, 

103 lookupFunction=ptcLookup, 

104 ) 

105 

106 inputPhotodiodeCorrection = cT.Input( 

107 name="pdCorrection", 

108 doc="Input photodiode correction.", 

109 storageClass="IsrCalib", 

110 dimensions=("instrument", ), 

111 isCalibration=True, 

112 ) 

113 

114 outputLinearizer = cT.Output( 

115 name="linearity", 

116 doc="Output linearity measurements.", 

117 storageClass="Linearizer", 

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

119 isCalibration=True, 

120 ) 

121 

122 def __init__(self, *, config=None): 

123 if not config.applyPhotodiodeCorrection: 

124 del self.inputPhotodiodeCorrection 

125 

126 

127class LinearitySolveConfig(pipeBase.PipelineTaskConfig, 

128 pipelineConnections=LinearitySolveConnections): 

129 """Configuration for solving the linearity from PTC dataset. 

130 """ 

131 linearityType = pexConfig.ChoiceField( 

132 dtype=str, 

133 doc="Type of linearizer to construct.", 

134 default="Squared", 

135 allowed={ 

136 "LookupTable": "Create a lookup table solution.", 

137 "Polynomial": "Create an arbitrary polynomial solution.", 

138 "Squared": "Create a single order squared solution.", 

139 "Spline": "Create a spline based solution.", 

140 "None": "Create a dummy solution.", 

141 } 

142 ) 

143 polynomialOrder = pexConfig.RangeField( 

144 dtype=int, 

145 doc="Degree of polynomial to fit. Must be at least 2.", 

146 default=3, 

147 min=2, 

148 ) 

149 splineKnots = pexConfig.Field( 

150 dtype=int, 

151 doc="Number of spline knots to use in fit.", 

152 default=10, 

153 ) 

154 maxLookupTableAdu = pexConfig.Field( 

155 dtype=int, 

156 doc="Maximum DN value for a LookupTable linearizer.", 

157 default=2**18, 

158 ) 

159 maxLinearAdu = pexConfig.Field( 

160 dtype=float, 

161 doc="Maximum DN value to use to estimate linear term.", 

162 default=20000.0, 

163 ) 

164 minLinearAdu = pexConfig.Field( 

165 dtype=float, 

166 doc="Minimum DN value to use to estimate linear term.", 

167 default=30.0, 

168 ) 

169 nSigmaClipLinear = pexConfig.Field( 

170 dtype=float, 

171 doc="Maximum deviation from linear solution for Poissonian noise.", 

172 default=5.0, 

173 ) 

174 ignorePtcMask = pexConfig.Field( 

175 dtype=bool, 

176 doc="Ignore the expIdMask set by the PTC solver?", 

177 default=False, 

178 ) 

179 usePhotodiode = pexConfig.Field( 

180 dtype=bool, 

181 doc="Use the photodiode info instead of the raw expTimes?", 

182 default=False, 

183 ) 

184 photodiodeIntegrationMethod = pexConfig.ChoiceField( 

185 dtype=str, 

186 doc="Integration method for photodiode monitoring data.", 

187 default="DIRECT_SUM", 

188 allowed={ 

189 "DIRECT_SUM": ("Use numpy's trapz integrator on all photodiode " 

190 "readout entries"), 

191 "TRIMMED_SUM": ("Use numpy's trapz integrator, clipping the " 

192 "leading and trailing entries, which are " 

193 "nominally at zero baseline level."), 

194 "CHARGE_SUM": ("Treat the current values as integrated charge " 

195 "over the sampling interval and simply sum " 

196 "the values, after subtracting a baseline level."), 

197 }, 

198 # TODO: remove on DM-40065. 

199 deprecated="This config has been moved to cpExtractPtcTask, and will be removed after v26.", 

200 ) 

201 photodiodeCurrentScale = pexConfig.Field( 

202 dtype=float, 

203 doc="Scale factor to apply to photodiode current values for the " 

204 "``CHARGE_SUM`` integration method.", 

205 default=-1.0, 

206 # TODO: remove on DM-40065. 

207 deprecated="This config has been moved to cpExtractPtcTask, and will be removed after v26.", 

208 ) 

209 applyPhotodiodeCorrection = pexConfig.Field( 

210 dtype=bool, 

211 doc="Calculate and apply a correction to the photodiode readings?", 

212 default=False, 

213 ) 

214 splineGroupingColumn = pexConfig.Field( 

215 dtype=str, 

216 doc="Column to use for grouping together points for Spline mode, to allow " 

217 "for different proportionality constants. If not set, no grouping " 

218 "will be done.", 

219 default=None, 

220 optional=True, 

221 ) 

222 splineGroupingMinPoints = pexConfig.Field( 

223 dtype=int, 

224 doc="Minimum number of linearity points to allow grouping together points " 

225 "for Spline mode with splineGroupingColumn. This configuration is here " 

226 "to prevent misuse of the Spline code to avoid over-fitting.", 

227 default=100, 

228 ) 

229 splineFitMinIter = pexConfig.Field( 

230 dtype=int, 

231 doc="Minimum number of iterations for spline fit.", 

232 default=3, 

233 ) 

234 splineFitMaxIter = pexConfig.Field( 

235 dtype=int, 

236 doc="Maximum number of iterations for spline fit.", 

237 default=20, 

238 ) 

239 splineFitMaxRejectionPerIteration = pexConfig.Field( 

240 dtype=int, 

241 doc="Maximum number of rejections per iteration for spline fit.", 

242 default=5, 

243 ) 

244 

245 

246class LinearitySolveTask(pipeBase.PipelineTask): 

247 """Fit the linearity from the PTC dataset. 

248 """ 

249 

250 ConfigClass = LinearitySolveConfig 

251 _DefaultName = 'cpLinearitySolve' 

252 

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

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

255 

256 Parameters 

257 ---------- 

258 butlerQC : `lsst.daf.butler.butlerQuantumContext.ButlerQuantumContext` 

259 Butler to operate on. 

260 inputRefs : `lsst.pipe.base.connections.InputQuantizedConnection` 

261 Input data refs to load. 

262 ouptutRefs : `lsst.pipe.base.connections.OutputQuantizedConnection` 

263 Output data refs to persist. 

264 """ 

265 inputs = butlerQC.get(inputRefs) 

266 

267 # Use the dimensions to set calib/provenance information. 

268 inputs['inputDims'] = inputRefs.inputPtc.dataId.byName() 

269 

270 outputs = self.run(**inputs) 

271 butlerQC.put(outputs, outputRefs) 

272 

273 def run(self, inputPtc, dummy, camera, inputDims, 

274 inputPhotodiodeCorrection=None): 

275 """Fit non-linearity to PTC data, returning the correct Linearizer 

276 object. 

277 

278 Parameters 

279 ---------- 

280 inputPtc : `lsst.ip.isr.PtcDataset` 

281 Pre-measured PTC dataset. 

282 dummy : `lsst.afw.image.Exposure` 

283 The exposure used to select the appropriate PTC dataset. 

284 In almost all circumstances, one of the input exposures 

285 used to generate the PTC dataset is the best option. 

286 inputPhotodiodeCorrection : `lsst.ip.isr.PhotodiodeCorrection` 

287 Pre-measured photodiode correction used in the case when 

288 applyPhotodiodeCorrection=True. 

289 camera : `lsst.afw.cameraGeom.Camera` 

290 Camera geometry. 

291 inputDims : `lsst.daf.butler.DataCoordinate` or `dict` 

292 DataIds to use to populate the output calibration. 

293 

294 Returns 

295 ------- 

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

297 The results struct containing: 

298 

299 ``outputLinearizer`` 

300 Final linearizer calibration (`lsst.ip.isr.Linearizer`). 

301 ``outputProvenance`` 

302 Provenance data for the new calibration 

303 (`lsst.ip.isr.IsrProvenance`). 

304 

305 Notes 

306 ----- 

307 This task currently fits only polynomial-defined corrections, 

308 where the correction coefficients are defined such that: 

309 :math:`corrImage = uncorrImage + \\sum_i c_i uncorrImage^(2 + i)` 

310 These :math:`c_i` are defined in terms of the direct polynomial fit: 

311 :math:`meanVector ~ P(x=timeVector) = \\sum_j k_j x^j` 

312 such that :math:`c_(j-2) = -k_j/(k_1^j)` in units of DN^(1-j) (c.f., 

313 Eq. 37 of 2003.05978). The `config.polynomialOrder` or 

314 `config.splineKnots` define the maximum order of :math:`x^j` to fit. 

315 As :math:`k_0` and :math:`k_1` are degenerate with bias level and gain, 

316 they are not included in the non-linearity correction. 

317 """ 

318 if len(dummy) == 0: 

319 self.log.warning("No dummy exposure found.") 

320 

321 detector = camera[inputDims['detector']] 

322 if self.config.linearityType == 'LookupTable': 

323 table = np.zeros((len(detector), self.config.maxLookupTableAdu), dtype=np.float32) 

324 tableIndex = 0 

325 else: 

326 table = None 

327 tableIndex = None # This will fail if we increment it. 

328 

329 # Initialize the linearizer. 

330 linearizer = Linearizer(detector=detector, table=table, log=self.log) 

331 linearizer.updateMetadataFromExposures([inputPtc]) 

332 if self.config.usePhotodiode and self.config.applyPhotodiodeCorrection: 

333 abscissaCorrections = inputPhotodiodeCorrection.abscissaCorrections 

334 

335 if self.config.linearityType == 'Spline': 

336 if self.config.splineGroupingColumn is not None: 

337 if self.config.splineGroupingColumn not in inputPtc.auxValues: 

338 raise ValueError(f"Config requests grouping by {self.config.splineGroupingColumn}, " 

339 "but this column is not available in inputPtc.auxValues.") 

340 groupingValue = inputPtc.auxValues[self.config.splineGroupingColumn] 

341 else: 

342 groupingValue = np.ones(len(inputPtc.rawMeans[inputPtc.ampNames[0]]), dtype=int) 

343 # We set this to have a value to fill the bad amps. 

344 fitOrder = self.config.splineKnots 

345 else: 

346 fitOrder = self.config.polynomialOrder 

347 

348 for i, amp in enumerate(detector): 

349 ampName = amp.getName() 

350 if ampName in inputPtc.badAmps: 

351 linearizer = self.fillBadAmp(linearizer, fitOrder, inputPtc, amp) 

352 self.log.warning("Amp %s in detector %s has no usable PTC information. Skipping!", 

353 ampName, detector.getName()) 

354 continue 

355 

356 # Check for too few points. 

357 if self.config.linearityType == "Spline" \ 

358 and self.config.splineGroupingColumn is not None \ 

359 and len(inputPtc.inputExpIdPairs[ampName]) < self.config.splineGroupingMinPoints: 

360 raise RuntimeError( 

361 "The input PTC has too few points to reliably run with PD grouping. " 

362 "The recommended course of action is to set splineGroupingColumn to None. " 

363 "If you really know what you are doing, you may reduce " 

364 "config.splineGroupingMinPoints.") 

365 

366 if (len(inputPtc.expIdMask[ampName]) == 0) or self.config.ignorePtcMask: 

367 self.log.warning("Mask not found for %s in detector %s in fit. Using all points.", 

368 ampName, detector.getName()) 

369 mask = np.ones(len(inputPtc.expIdMask[ampName]), dtype=bool) 

370 else: 

371 mask = inputPtc.expIdMask[ampName].copy() 

372 

373 if self.config.usePhotodiode: 

374 modExpTimes = inputPtc.photoCharges[ampName].copy() 

375 # Make sure any exposure pairs that do not have photodiode data 

376 # are masked. 

377 mask[~np.isfinite(modExpTimes)] = False 

378 

379 # Get the photodiode correction. 

380 if self.config.applyPhotodiodeCorrection: 

381 for j, pair in enumerate(inputPtc.inputExpIdPairs[ampName]): 

382 try: 

383 correction = abscissaCorrections[str(pair)] 

384 except KeyError: 

385 correction = 0.0 

386 modExpTimes[j] += correction 

387 

388 inputAbscissa = modExpTimes 

389 else: 

390 inputAbscissa = inputPtc.rawExpTimes[ampName].copy() 

391 

392 inputOrdinate = inputPtc.rawMeans[ampName].copy() 

393 

394 mask &= (inputOrdinate < self.config.maxLinearAdu) 

395 mask &= (inputOrdinate > self.config.minLinearAdu) 

396 

397 if mask.sum() < 2: 

398 linearizer = self.fillBadAmp(linearizer, fitOrder, inputPtc, amp) 

399 self.log.warning("Amp %s in detector %s has not enough points for fit. Skipping!", 

400 ampName, detector.getName()) 

401 continue 

402 

403 if self.config.linearityType != 'Spline': 

404 linearFit, linearFitErr, chiSq, weights = irlsFit([0.0, 100.0], inputAbscissa[mask], 

405 inputOrdinate[mask], funcPolynomial) 

406 

407 # Convert this proxy-to-flux fit into an expected linear flux 

408 linearOrdinate = linearFit[0] + linearFit[1] * inputAbscissa 

409 # Exclude low end outliers. 

410 # This is compared to the original values. 

411 threshold = self.config.nSigmaClipLinear * np.sqrt(abs(inputOrdinate)) 

412 

413 mask[np.abs(inputOrdinate - linearOrdinate) >= threshold] = False 

414 

415 if mask.sum() < 2: 

416 linearizer = self.fillBadAmp(linearizer, fitOrder, inputPtc, amp) 

417 self.log.warning("Amp %s in detector %s has not enough points in linear ordinate. " 

418 "Skipping!", ampName, detector.getName()) 

419 continue 

420 

421 self.debugFit('linearFit', inputAbscissa, inputOrdinate, linearOrdinate, mask, ampName) 

422 

423 # Do fits 

424 if self.config.linearityType in ['Polynomial', 'Squared', 'LookupTable']: 

425 polyFit = np.zeros(fitOrder + 1) 

426 polyFit[1] = 1.0 

427 polyFit, polyFitErr, chiSq, weights = irlsFit(polyFit, linearOrdinate[mask], 

428 inputOrdinate[mask], funcPolynomial) 

429 

430 # Truncate the polynomial fit to the squared term. 

431 k1 = polyFit[1] 

432 linearityCoeffs = np.array( 

433 [-coeff/(k1**order) for order, coeff in enumerate(polyFit)] 

434 )[2:] 

435 significant = np.where(np.abs(linearityCoeffs) > 1e-10) 

436 self.log.info("Significant polynomial fits: %s", significant) 

437 

438 modelOrdinate = funcPolynomial(polyFit, linearOrdinate) 

439 

440 self.debugFit( 

441 'polyFit', 

442 inputAbscissa[mask], 

443 inputOrdinate[mask], 

444 modelOrdinate[mask], 

445 None, 

446 ampName, 

447 ) 

448 

449 if self.config.linearityType == 'Squared': 

450 # The first term is the squared term. 

451 linearityCoeffs = linearityCoeffs[0: 1] 

452 elif self.config.linearityType == 'LookupTable': 

453 # Use linear part to get time at which signal is 

454 # maxAduForLookupTableLinearizer DN 

455 tMax = (self.config.maxLookupTableAdu - polyFit[0])/polyFit[1] 

456 timeRange = np.linspace(0, tMax, self.config.maxLookupTableAdu) 

457 signalIdeal = polyFit[0] + polyFit[1]*timeRange 

458 signalUncorrected = funcPolynomial(polyFit, timeRange) 

459 lookupTableRow = signalIdeal - signalUncorrected # LinearizerLookupTable has correction 

460 

461 linearizer.tableData[tableIndex, :] = lookupTableRow 

462 linearityCoeffs = np.array([tableIndex, 0]) 

463 tableIndex += 1 

464 elif self.config.linearityType in ['Spline']: 

465 # This is a spline fit with photodiode data based on a model 

466 # from Pierre Astier. 

467 # This model fits a spline with (optional) nuisance parameters 

468 # to allow for different linearity coefficients with different 

469 # photodiode settings. The minimization is a least-squares 

470 # fit with the residual of 

471 # Sum[(S(mu_i) + mu_i)/(k_j * D_i) - 1]**2, where S(mu_i) is 

472 # an Akima Spline function of mu_i, the observed flat-pair 

473 # mean; D_j is the photo-diode measurement corresponding to 

474 # that flat-pair; and k_j is a constant of proportionality 

475 # which is over index j as it is allowed to 

476 # be different based on different photodiode settings (e.g. 

477 # CCOBCURR). 

478 

479 # The fit has additional constraints to ensure that the spline 

480 # goes through the (0, 0) point, as well as a normalization 

481 # condition so that the average of the spline over the full 

482 # range is 0. The normalization ensures that the spline only 

483 # fits deviations from linearity, rather than the linear 

484 # function itself which is degenerate with the gain. 

485 

486 nodes = np.linspace(0.0, np.max(inputOrdinate[mask]), self.config.splineKnots) 

487 

488 fitter = AstierSplineLinearityFitter( 

489 nodes, 

490 groupingValue, 

491 inputAbscissa, 

492 inputOrdinate, 

493 mask=mask, 

494 log=self.log, 

495 ) 

496 p0 = fitter.estimate_p0() 

497 pars = fitter.fit( 

498 p0, 

499 min_iter=self.config.splineFitMinIter, 

500 max_iter=self.config.splineFitMaxIter, 

501 max_rejection_per_iteration=self.config.splineFitMaxRejectionPerIteration, 

502 n_sigma_clip=self.config.nSigmaClipLinear, 

503 ) 

504 

505 # Confirm that the first parameter is 0, and set it to 

506 # exactly zero. 

507 if not np.isclose(pars[0], 0): 

508 raise RuntimeError("Programmer error! First spline parameter must " 

509 "be consistent with zero.") 

510 pars[0] = 0.0 

511 

512 linearityCoeffs = np.concatenate([nodes, pars[0: len(nodes)]]) 

513 linearFit = np.array([0.0, np.mean(pars[len(nodes):])]) 

514 

515 # We modify the inputAbscissa according to the linearity fits 

516 # here, for proper residual computation. 

517 for j, group_index in enumerate(fitter.group_indices): 

518 inputOrdinate[group_index] /= (pars[len(nodes) + j] / linearFit[1]) 

519 

520 linearOrdinate = linearFit[1] * inputOrdinate 

521 # For the spline fit, reuse the "polyFit -> fitParams" 

522 # field to record the linear coefficients for the groups. 

523 polyFit = pars[len(nodes):] 

524 polyFitErr = np.zeros_like(polyFit) 

525 chiSq = np.nan 

526 

527 # Update mask based on what the fitter rejected. 

528 mask = fitter.mask 

529 else: 

530 polyFit = np.zeros(1) 

531 polyFitErr = np.zeros(1) 

532 chiSq = np.nan 

533 linearityCoeffs = np.zeros(1) 

534 

535 linearizer.linearityType[ampName] = self.config.linearityType 

536 linearizer.linearityCoeffs[ampName] = linearityCoeffs 

537 linearizer.linearityBBox[ampName] = amp.getBBox() 

538 linearizer.fitParams[ampName] = polyFit 

539 linearizer.fitParamsErr[ampName] = polyFitErr 

540 linearizer.fitChiSq[ampName] = chiSq 

541 linearizer.linearFit[ampName] = linearFit 

542 

543 image = afwImage.ImageF(len(inputOrdinate), 1) 

544 image.array[:, :] = inputOrdinate 

545 linearizeFunction = linearizer.getLinearityTypeByName(linearizer.linearityType[ampName]) 

546 linearizeFunction()( 

547 image, 

548 **{'coeffs': linearizer.linearityCoeffs[ampName], 

549 'table': linearizer.tableData, 

550 'log': linearizer.log} 

551 ) 

552 linearizeModel = image.array[0, :] 

553 

554 # The residuals that we record are the final residuals compared to 

555 # a linear model, after everything has been (properly?) linearized. 

556 postLinearFit, _, _, _ = irlsFit( 

557 [0.0, 100.0], 

558 inputAbscissa[mask], 

559 linearizeModel[mask], 

560 funcPolynomial, 

561 ) 

562 residuals = linearizeModel - (postLinearFit[0] + postLinearFit[1] * inputAbscissa) 

563 # We set masked residuals to nan. 

564 residuals[~mask] = np.nan 

565 

566 linearizer.fitResiduals[ampName] = residuals 

567 

568 self.debugFit( 

569 'solution', 

570 inputOrdinate[mask], 

571 linearOrdinate[mask], 

572 linearizeModel[mask], 

573 None, 

574 ampName, 

575 ) 

576 

577 linearizer.hasLinearity = True 

578 linearizer.validate() 

579 linearizer.updateMetadata(camera=camera, detector=detector, filterName='NONE') 

580 linearizer.updateMetadata(setDate=True, setCalibId=True) 

581 provenance = IsrProvenance(calibType='linearizer') 

582 

583 return pipeBase.Struct( 

584 outputLinearizer=linearizer, 

585 outputProvenance=provenance, 

586 ) 

587 

588 def fillBadAmp(self, linearizer, fitOrder, inputPtc, amp): 

589 # Need to fill linearizer with empty values 

590 # if the amp is non-functional 

591 ampName = amp.getName() 

592 nEntries = 1 

593 pEntries = 1 

594 if self.config.linearityType in ['Polynomial']: 

595 nEntries = fitOrder + 1 

596 pEntries = fitOrder + 1 

597 elif self.config.linearityType in ['Spline']: 

598 nEntries = fitOrder * 2 

599 elif self.config.linearityType in ['Squared', 'None']: 

600 nEntries = 1 

601 pEntries = fitOrder + 1 

602 elif self.config.linearityType in ['LookupTable']: 

603 nEntries = 2 

604 pEntries = fitOrder + 1 

605 

606 linearizer.linearityType[ampName] = "None" 

607 linearizer.linearityCoeffs[ampName] = np.zeros(nEntries) 

608 linearizer.linearityBBox[ampName] = amp.getBBox() 

609 linearizer.fitParams[ampName] = np.zeros(pEntries) 

610 linearizer.fitParamsErr[ampName] = np.zeros(pEntries) 

611 linearizer.fitChiSq[ampName] = np.nan 

612 linearizer.fitResiduals[ampName] = np.zeros(len(inputPtc.expIdMask[ampName])) 

613 linearizer.linearFit[ampName] = np.zeros(2) 

614 return linearizer 

615 

616 def debugFit(self, stepname, xVector, yVector, yModel, mask, ampName): 

617 """Debug method for linearity fitting. 

618 

619 Parameters 

620 ---------- 

621 stepname : `str` 

622 A label to use to check if we care to debug at a given 

623 line of code. 

624 xVector : `numpy.array`, (N,) 

625 The values to use as the independent variable in the 

626 linearity fit. 

627 yVector : `numpy.array`, (N,) 

628 The values to use as the dependent variable in the 

629 linearity fit. 

630 yModel : `numpy.array`, (N,) 

631 The values to use as the linearized result. 

632 mask : `numpy.array` [`bool`], (N,) , optional 

633 A mask to indicate which entries of ``xVector`` and 

634 ``yVector`` to keep. 

635 ampName : `str` 

636 Amplifier name to lookup linearity correction values. 

637 """ 

638 frame = getDebugFrame(self._display, stepname) 

639 if frame: 

640 import matplotlib.pyplot as plt 

641 fig, axs = plt.subplots(2) 

642 

643 if mask is None: 

644 mask = np.ones_like(xVector, dtype=bool) 

645 

646 fig.suptitle(f"{stepname} {ampName} {self.config.linearityType}") 

647 if stepname == 'linearFit': 

648 axs[0].set_xlabel("Input Abscissa (time or mondiode)") 

649 axs[0].set_ylabel("Input Ordinate (flux)") 

650 axs[1].set_xlabel("Linear Ordinate (linear flux)") 

651 axs[1].set_ylabel("Flux Difference: (input - linear)") 

652 elif stepname in ('polyFit', 'splineFit'): 

653 axs[0].set_xlabel("Linear Abscissa (linear flux)") 

654 axs[0].set_ylabel("Input Ordinate (flux)") 

655 axs[1].set_xlabel("Linear Ordinate (linear flux)") 

656 axs[1].set_ylabel("Flux Difference: (input - full model fit)") 

657 elif stepname == 'solution': 

658 axs[0].set_xlabel("Input Abscissa (time or mondiode)") 

659 axs[0].set_ylabel("Linear Ordinate (linear flux)") 

660 axs[1].set_xlabel("Model flux (linear flux)") 

661 axs[1].set_ylabel("Flux Difference: (linear - model)") 

662 

663 axs[0].set_yscale('log') 

664 axs[0].set_xscale('log') 

665 axs[0].scatter(xVector, yVector) 

666 axs[0].scatter(xVector[~mask], yVector[~mask], c='red', marker='x') 

667 axs[1].set_xscale('log') 

668 

669 axs[1].scatter(yModel, yVector[mask] - yModel) 

670 fig.tight_layout() 

671 fig.show() 

672 

673 prompt = "Press Enter or c to continue [chpx]..." 

674 while True: 

675 ans = input(prompt).lower() 

676 if ans in ("", " ", "c",): 

677 break 

678 elif ans in ("p", ): 

679 import pdb 

680 pdb.set_trace() 

681 elif ans in ("h", ): 

682 print("[h]elp [c]ontinue [p]db") 

683 elif ans in ('x', ): 

684 exit() 

685 plt.close()