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

305 statements  

« prev     ^ index     » next       coverage.py v7.5.1, created at 2024-05-15 02:24 -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# 

22 

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

24 

25import numpy as np 

26from scipy.stats import median_abs_deviation 

27 

28import lsst.afw.image as afwImage 

29import lsst.pipe.base as pipeBase 

30import lsst.pipe.base.connectionTypes as cT 

31import lsst.pex.config as pexConfig 

32 

33from lsstDebug import getDebugFrame 

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

35 

36from .utils import (funcPolynomial, irlsFit, AstierSplineLinearityFitter, 

37 extractCalibDate) 

38 

39 

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

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

42 

43 Parameters 

44 ---------- 

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

46 Dataset type to look up. 

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

48 Registry for the data repository being searched. 

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

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

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

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

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

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

55 would allow a real validity-range lookup. 

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

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

58 by this function (see notes below). 

59 

60 Returns 

61 ------- 

62 refs : `list` [ `DatasetRef` ] 

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

64 dataset, if one was found. 

65 

66 Raises 

67 ------ 

68 RuntimeError 

69 Raised if more than one PTC reference is found. 

70 """ 

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

72 findFirst=False)) 

73 if len(refs) >= 2: 

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

75 

76 return refs 

77 

78 

79class LinearitySolveConnections(pipeBase.PipelineTaskConnections, 

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

81 dummy = cT.Input( 

82 name="raw", 

83 doc="Dummy exposure.", 

84 storageClass='Exposure', 

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

86 multiple=True, 

87 deferLoad=True, 

88 ) 

89 

90 camera = cT.PrerequisiteInput( 

91 name="camera", 

92 doc="Camera Geometry definition.", 

93 storageClass="Camera", 

94 dimensions=("instrument", ), 

95 isCalibration=True, 

96 ) 

97 

98 inputPtc = cT.PrerequisiteInput( 

99 name="ptc", 

100 doc="Input PTC dataset.", 

101 storageClass="PhotonTransferCurveDataset", 

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

103 isCalibration=True, 

104 lookupFunction=ptcLookup, 

105 ) 

106 

107 inputPhotodiodeCorrection = cT.Input( 

108 name="pdCorrection", 

109 doc="Input photodiode correction.", 

110 storageClass="IsrCalib", 

111 dimensions=("instrument", ), 

112 isCalibration=True, 

113 ) 

114 

115 outputLinearizer = cT.Output( 

116 name="linearity", 

117 doc="Output linearity measurements.", 

118 storageClass="Linearizer", 

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

120 isCalibration=True, 

121 ) 

122 

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

124 if not config.applyPhotodiodeCorrection: 

125 del self.inputPhotodiodeCorrection 

126 

127 

128class LinearitySolveConfig(pipeBase.PipelineTaskConfig, 

129 pipelineConnections=LinearitySolveConnections): 

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

131 """ 

132 linearityType = pexConfig.ChoiceField( 

133 dtype=str, 

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

135 default="Squared", 

136 allowed={ 

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

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

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

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

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

142 } 

143 ) 

144 polynomialOrder = pexConfig.RangeField( 

145 dtype=int, 

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

147 default=3, 

148 min=2, 

149 ) 

150 splineKnots = pexConfig.Field( 

151 dtype=int, 

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

153 default=10, 

154 ) 

155 maxLookupTableAdu = pexConfig.Field( 

156 dtype=int, 

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

158 default=2**18, 

159 ) 

160 maxLinearAdu = pexConfig.Field( 

161 dtype=float, 

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

163 default=20000.0, 

164 ) 

165 minLinearAdu = pexConfig.Field( 

166 dtype=float, 

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

168 default=30.0, 

169 ) 

170 nSigmaClipLinear = pexConfig.Field( 

171 dtype=float, 

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

173 default=5.0, 

174 ) 

175 ignorePtcMask = pexConfig.Field( 

176 dtype=bool, 

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

178 default=False, 

179 ) 

180 usePhotodiode = pexConfig.Field( 

181 dtype=bool, 

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

183 default=False, 

184 ) 

185 applyPhotodiodeCorrection = pexConfig.Field( 

186 dtype=bool, 

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

188 default=False, 

189 ) 

190 minPhotodiodeCurrent = pexConfig.Field( 

191 dtype=float, 

192 doc="Minimum value to trust photodiode signals.", 

193 default=0.0, 

194 ) 

195 splineGroupingColumn = pexConfig.Field( 

196 dtype=str, 

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

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

199 "will be done.", 

200 default=None, 

201 optional=True, 

202 ) 

203 splineGroupingMinPoints = pexConfig.Field( 

204 dtype=int, 

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

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

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

208 default=100, 

209 ) 

210 splineFitMinIter = pexConfig.Field( 

211 dtype=int, 

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

213 default=3, 

214 ) 

215 splineFitMaxIter = pexConfig.Field( 

216 dtype=int, 

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

218 default=20, 

219 ) 

220 splineFitMaxRejectionPerIteration = pexConfig.Field( 

221 dtype=int, 

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

223 default=5, 

224 ) 

225 doSplineFitOffset = pexConfig.Field( 

226 dtype=bool, 

227 doc="Fit a scattered light offset in the spline fit.", 

228 default=True, 

229 ) 

230 doSplineFitWeights = pexConfig.Field( 

231 dtype=bool, 

232 doc="Fit linearity weight parameters in the spline fit.", 

233 default=False, 

234 ) 

235 splineFitWeightParsStart = pexConfig.ListField( 

236 dtype=float, 

237 doc="Starting parameters for weight fit, if doSplineFitWeights=True. " 

238 "Parameters are such that sigma = sqrt(par[0]**2. + par[1]**2./mu)." 

239 "If doSplineFitWeights=False then these are used as-is; otherwise " 

240 "they are used as the initial values for fitting these parameters.", 

241 length=2, 

242 default=[1.0, 0.0], 

243 ) 

244 doSplineFitTemperature = pexConfig.Field( 

245 dtype=bool, 

246 doc="Fit temperature coefficient in spline fit?", 

247 default=False, 

248 ) 

249 splineFitTemperatureColumn = pexConfig.Field( 

250 dtype=str, 

251 doc="Name of the temperature column to use when fitting temperature " 

252 "coefficients in spline fit; this must not be None if " 

253 "doSplineFitTemperature is True.", 

254 default=None, 

255 optional=True, 

256 ) 

257 

258 def validate(self): 

259 super().validate() 

260 

261 if self.doSplineFitTemperature and self.splineFitTemperatureColumn is None: 

262 raise ValueError("Must set splineFitTemperatureColumn if doSplineFitTemperature is True.") 

263 

264 

265class LinearitySolveTask(pipeBase.PipelineTask): 

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

267 """ 

268 

269 ConfigClass = LinearitySolveConfig 

270 _DefaultName = 'cpLinearitySolve' 

271 

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

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

274 

275 Parameters 

276 ---------- 

277 butlerQC : `lsst.daf.butler.QuantumContext` 

278 Butler to operate on. 

279 inputRefs : `lsst.pipe.base.InputQuantizedConnection` 

280 Input data refs to load. 

281 ouptutRefs : `lsst.pipe.base.OutputQuantizedConnection` 

282 Output data refs to persist. 

283 """ 

284 inputs = butlerQC.get(inputRefs) 

285 

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

287 inputs['inputDims'] = dict(inputRefs.inputPtc.dataId.required) 

288 

289 # Add calibration provenance info to header. 

290 kwargs = dict() 

291 reference = getattr(inputRefs, "inputPtc", None) 

292 

293 if reference is not None and hasattr(reference, "run"): 

294 runKey = "PTC_RUN" 

295 runValue = reference.run 

296 idKey = "PTC_UUID" 

297 idValue = str(reference.id) 

298 dateKey = "PTC_DATE" 

299 calib = inputs.get("inputPtc", None) 

300 dateValue = extractCalibDate(calib) 

301 

302 kwargs[runKey] = runValue 

303 kwargs[idKey] = idValue 

304 kwargs[dateKey] = dateValue 

305 

306 self.log.info("Using " + str(reference.run)) 

307 

308 outputs = self.run(**inputs) 

309 outputs.outputLinearizer.updateMetadata(setDate=False, **kwargs) 

310 

311 butlerQC.put(outputs, outputRefs) 

312 

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

314 inputPhotodiodeCorrection=None): 

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

316 object. 

317 

318 Parameters 

319 ---------- 

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

321 Pre-measured PTC dataset. 

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

323 The exposure used to select the appropriate PTC dataset. 

324 In almost all circumstances, one of the input exposures 

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

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

327 Pre-measured photodiode correction used in the case when 

328 applyPhotodiodeCorrection=True. 

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

330 Camera geometry. 

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

332 DataIds to use to populate the output calibration. 

333 

334 Returns 

335 ------- 

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

337 The results struct containing: 

338 

339 ``outputLinearizer`` 

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

341 ``outputProvenance`` 

342 Provenance data for the new calibration 

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

344 

345 Notes 

346 ----- 

347 This task currently fits only polynomial-defined corrections, 

348 where the correction coefficients are defined such that: 

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

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

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

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

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

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

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

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

357 """ 

358 if len(dummy) == 0: 

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

360 

361 detector = camera[inputDims['detector']] 

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

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

364 tableIndex = 0 

365 else: 

366 table = None 

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

368 

369 # Initialize the linearizer. 

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

371 linearizer.updateMetadataFromExposures([inputPtc]) 

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

373 abscissaCorrections = inputPhotodiodeCorrection.abscissaCorrections 

374 

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

376 if self.config.splineGroupingColumn is not None: 

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

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

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

380 groupingValues = inputPtc.auxValues[self.config.splineGroupingColumn] 

381 else: 

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

383 

384 if self.config.doSplineFitTemperature: 

385 if self.config.splineFitTemperatureColumn not in inputPtc.auxValues: 

386 raise ValueError("Config requests fitting temperature coefficient for " 

387 f"{self.config.splineFitTemperatureColumn} but this column " 

388 "is not available in inputPtc.auxValues.") 

389 temperatureValues = inputPtc.auxValues[self.config.splineFitTemperatureColumn] 

390 else: 

391 temperatureValues = None 

392 

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

394 fitOrder = self.config.splineKnots 

395 else: 

396 fitOrder = self.config.polynomialOrder 

397 

398 for i, amp in enumerate(detector): 

399 ampName = amp.getName() 

400 if ampName in inputPtc.badAmps: 

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

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

403 ampName, detector.getName()) 

404 continue 

405 

406 # Check for too few points. 

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

408 and self.config.splineGroupingColumn is not None \ 

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

410 raise RuntimeError( 

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

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

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

414 "config.splineGroupingMinPoints.") 

415 

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

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

418 ampName, detector.getName()) 

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

420 else: 

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

422 

423 if self.config.linearityType == "Spline" and temperatureValues is not None: 

424 mask &= np.isfinite(temperatureValues) 

425 

426 if self.config.usePhotodiode: 

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

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

429 # are masked. 

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

431 

432 # Make sure any photodiode measurements below the configured 

433 # minimum are masked. 

434 mask[modExpTimes < self.config.minPhotodiodeCurrent] = False 

435 

436 # Get the photodiode correction. 

437 if self.config.applyPhotodiodeCorrection: 

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

439 try: 

440 correction = abscissaCorrections[str(pair)] 

441 except KeyError: 

442 correction = 0.0 

443 modExpTimes[j] += correction 

444 

445 inputAbscissa = modExpTimes 

446 else: 

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

448 

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

450 

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

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

453 

454 if mask.sum() < 2: 

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

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

457 ampName, detector.getName()) 

458 continue 

459 

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

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

462 inputOrdinate[mask], funcPolynomial) 

463 

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

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

466 # Exclude low end outliers. 

467 # This is compared to the original values. 

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

469 

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

471 

472 if mask.sum() < 2: 

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

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

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

476 continue 

477 

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

479 

480 # Do fits 

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

482 polyFit = np.zeros(fitOrder + 1) 

483 polyFit[1] = 1.0 

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

485 inputOrdinate[mask], funcPolynomial) 

486 

487 # Truncate the polynomial fit to the squared term. 

488 k1 = polyFit[1] 

489 linearityCoeffs = np.array( 

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

491 )[2:] 

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

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

494 

495 modelOrdinate = funcPolynomial(polyFit, linearOrdinate) 

496 

497 self.debugFit( 

498 'polyFit', 

499 inputAbscissa[mask], 

500 inputOrdinate[mask], 

501 modelOrdinate[mask], 

502 None, 

503 ampName, 

504 ) 

505 

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

507 # The first term is the squared term. 

508 linearityCoeffs = linearityCoeffs[0: 1] 

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

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

511 # maxAduForLookupTableLinearizer DN 

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

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

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

515 signalUncorrected = funcPolynomial(polyFit, timeRange) 

516 lookupTableRow = signalIdeal - signalUncorrected # LinearizerLookupTable has correction 

517 

518 linearizer.tableData[tableIndex, :] = lookupTableRow 

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

520 tableIndex += 1 

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

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

523 # from Pierre Astier. 

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

525 # to allow for different linearity coefficients with different 

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

527 # fit with the residual of 

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

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

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

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

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

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

534 # CCOBCURR); and O is a constant offset to allow for light 

535 # leaks (and is only fit if doSplineFitOffset=True). In 

536 # addition, if config.doSplineFitTemperature is True then 

537 # the fit will adjust mu such that 

538 # mu = mu_input*(1 + alpha*(T - T_ref)) 

539 # and T_ref is taken as the median temperature of the run. 

540 

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

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

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

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

545 # fits deviations from linearity, rather than the linear 

546 # function itself which is degenerate with the gain. 

547 

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

549 

550 if temperatureValues is not None: 

551 temperatureValuesScaled = temperatureValues - np.median(temperatureValues[~mask]) 

552 else: 

553 temperatureValuesScaled = None 

554 

555 fitter = AstierSplineLinearityFitter( 

556 nodes, 

557 groupingValues, 

558 inputAbscissa, 

559 inputOrdinate, 

560 mask=mask, 

561 log=self.log, 

562 fit_offset=self.config.doSplineFitOffset, 

563 fit_weights=self.config.doSplineFitWeights, 

564 weight_pars_start=self.config.splineFitWeightParsStart, 

565 fit_temperature=self.config.doSplineFitTemperature, 

566 temperature_scaled=temperatureValuesScaled, 

567 ) 

568 p0 = fitter.estimate_p0() 

569 pars = fitter.fit( 

570 p0, 

571 min_iter=self.config.splineFitMinIter, 

572 max_iter=self.config.splineFitMaxIter, 

573 max_rejection_per_iteration=self.config.splineFitMaxRejectionPerIteration, 

574 n_sigma_clip=self.config.nSigmaClipLinear, 

575 ) 

576 

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

578 # exactly zero. 

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

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

581 "be consistent with zero.") 

582 pars[0] = 0.0 

583 

584 linearityChisq = fitter.compute_chisq_dof(pars) 

585 

586 linearityCoeffs = np.concatenate([nodes, pars[fitter.par_indices["values"]]]) 

587 linearFit = np.array([0.0, np.mean(pars[fitter.par_indices["groups"]])]) 

588 

589 # We must modify the inputOrdinate according to the 

590 # nuisance terms in the linearity fit for the residual 

591 # computation code to work properly. 

592 # The true mu (inputOrdinate) is given by 

593 # mu = mu_in * (1 + alpha*t_scale) 

594 if self.config.doSplineFitTemperature: 

595 inputOrdinate *= (1.0 

596 + pars[fitter.par_indices["temperature_coeff"]]*temperatureValuesScaled) 

597 # Divide by the relative scaling of the different groups. 

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

599 inputOrdinate[group_index] /= (pars[fitter.par_indices["groups"][j]] / linearFit[1]) 

600 # And remove the offset term. 

601 if self.config.doSplineFitOffset: 

602 inputOrdinate -= pars[fitter.par_indices["offset"]] 

603 

604 linearOrdinate = linearFit[1] * inputOrdinate 

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

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

607 # We additionally append the offset and weight_pars; 

608 # however these will be zero-length arrays if these were 

609 # not configured to be fit. 

610 polyFit = np.concatenate(( 

611 pars[fitter.par_indices["groups"]], 

612 pars[fitter.par_indices["offset"]], 

613 pars[fitter.par_indices["weight_pars"]], 

614 pars[fitter.par_indices["temperature_coeff"]], 

615 )) 

616 polyFitErr = np.zeros_like(polyFit) 

617 chiSq = linearityChisq 

618 

619 # Update mask based on what the fitter rejected. 

620 mask = fitter.mask 

621 else: 

622 polyFit = np.zeros(1) 

623 polyFitErr = np.zeros(1) 

624 chiSq = np.nan 

625 linearityCoeffs = np.zeros(1) 

626 

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

628 linearizer.linearityCoeffs[ampName] = linearityCoeffs 

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

630 linearizer.fitParams[ampName] = polyFit 

631 linearizer.fitParamsErr[ampName] = polyFitErr 

632 linearizer.fitChiSq[ampName] = chiSq 

633 linearizer.linearFit[ampName] = linearFit 

634 

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

636 image.array[:, :] = inputOrdinate 

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

638 linearizeFunction()( 

639 image, 

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

641 'table': linearizer.tableData, 

642 'log': linearizer.log} 

643 ) 

644 linearizeModel = image.array[0, :] 

645 

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

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

648 if mask.sum() < 2: 

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

650 "for residuals. Skipping!", ampName, detector.getName()) 

651 residuals = np.full_like(linearizeModel, np.nan) 

652 else: 

653 postLinearFit, _, _, _ = irlsFit( 

654 [0.0, 100.0], 

655 inputAbscissa[mask], 

656 linearizeModel[mask], 

657 funcPolynomial, 

658 ) 

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

660 # We set masked residuals to nan. 

661 residuals[~mask] = np.nan 

662 

663 linearizer.fitResiduals[ampName] = residuals 

664 

665 finite = np.isfinite(residuals) 

666 if finite.sum() == 0: 

667 sigmad = np.nan 

668 else: 

669 sigmad = median_abs_deviation(residuals[finite]/inputOrdinate[finite], scale="normal") 

670 linearizer.fitResidualsSigmaMad[ampName] = sigmad 

671 

672 self.debugFit( 

673 'solution', 

674 inputOrdinate[mask], 

675 linearOrdinate[mask], 

676 linearizeModel[mask], 

677 None, 

678 ampName, 

679 ) 

680 

681 self.fixupBadAmps(linearizer) 

682 

683 linearizer.hasLinearity = True 

684 linearizer.validate() 

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

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

687 provenance = IsrProvenance(calibType='linearizer') 

688 

689 return pipeBase.Struct( 

690 outputLinearizer=linearizer, 

691 outputProvenance=provenance, 

692 ) 

693 

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

695 # Need to fill linearizer with empty values 

696 # if the amp is non-functional 

697 ampName = amp.getName() 

698 nEntries = 1 

699 pEntries = 1 

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

701 # We discard the first 2 entries in the polynomial. 

702 nEntries = fitOrder + 1 - 2 

703 pEntries = fitOrder + 1 - 2 

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

705 nEntries = fitOrder * 2 

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

707 nEntries = 1 

708 pEntries = fitOrder + 1 

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

710 nEntries = 2 

711 pEntries = fitOrder + 1 

712 

713 linearizer.linearityType[ampName] = "None" 

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

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

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

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

718 linearizer.fitChiSq[ampName] = np.nan 

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

720 linearizer.fitResidualsSigmaMad[ampName] = np.nan 

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

722 return linearizer 

723 

724 def fixupBadAmps(self, linearizer): 

725 """Fix nan padding in bad amplifiers. 

726 

727 Parameters 

728 ---------- 

729 linearizer : `lsst.ip.isr.Linearizer` 

730 """ 

731 fitParamsMaxLen = 0 

732 for ampName in linearizer.ampNames: 

733 if (length := len(linearizer.fitParams[ampName])) > fitParamsMaxLen: 

734 fitParamsMaxLen = length 

735 

736 for ampName in linearizer.ampNames: 

737 if linearizer.linearityType[ampName] == "None": 

738 # Bad amplifier. 

739 linearizer.fitParams[ampName] = np.zeros(fitParamsMaxLen) 

740 linearizer.fitParamsErr[ampName] = np.zeros(fitParamsMaxLen) 

741 elif len(linearizer.fitParams[ampName]) != fitParamsMaxLen: 

742 raise RuntimeError("Linearity has mismatched fitParams; check code/data.") 

743 

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

745 """Debug method for linearity fitting. 

746 

747 Parameters 

748 ---------- 

749 stepname : `str` 

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

751 line of code. 

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

753 The values to use as the independent variable in the 

754 linearity fit. 

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

756 The values to use as the dependent variable in the 

757 linearity fit. 

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

759 The values to use as the linearized result. 

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

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

762 ``yVector`` to keep. 

763 ampName : `str` 

764 Amplifier name to lookup linearity correction values. 

765 """ 

766 frame = getDebugFrame(self._display, stepname) 

767 if frame: 

768 import matplotlib.pyplot as plt 

769 fig, axs = plt.subplots(2) 

770 

771 if mask is None: 

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

773 

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

775 if stepname == 'linearFit': 

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

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

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

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

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

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

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

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

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

785 elif stepname == 'solution': 

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

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

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

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

790 

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

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

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

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

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

796 

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

798 fig.tight_layout() 

799 fig.show() 

800 

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

802 while True: 

803 ans = input(prompt).lower() 

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

805 break 

806 elif ans in ("p", ): 

807 import pdb 

808 pdb.set_trace() 

809 elif ans in ("h", ): 

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

811 elif ans in ('x', ): 

812 exit() 

813 plt.close()