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

394 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__ = ['ddict2dict', 'CovFastFourierTransform'] 

24 

25 

26import galsim 

27import logging 

28import numpy as np 

29import itertools 

30import numpy.polynomial.polynomial as poly 

31 

32from scipy.optimize import leastsq 

33from scipy.stats import median_abs_deviation, norm 

34 

35from lsst.ip.isr import isrMock 

36import lsst.afw.image 

37import lsst.afw.math 

38 

39 

40def sigmaClipCorrection(nSigClip): 

41 """Correct measured sigma to account for clipping. 

42 

43 If we clip our input data and then measure sigma, then the 

44 measured sigma is smaller than the true value because real 

45 points beyond the clip threshold have been removed. This is a 

46 small (1.5% at nSigClip=3) effect when nSigClip >~ 3, but the 

47 default parameters for measure crosstalk use nSigClip=2.0. 

48 This causes the measured sigma to be about 15% smaller than 

49 real. This formula corrects the issue, for the symmetric case 

50 (upper clip threshold equal to lower clip threshold). 

51 

52 Parameters 

53 ---------- 

54 nSigClip : `float` 

55 Number of sigma the measurement was clipped by. 

56 

57 Returns 

58 ------- 

59 scaleFactor : `float` 

60 Scale factor to increase the measured sigma by. 

61 """ 

62 varFactor = 1.0 - (2 * nSigClip * norm.pdf(nSigClip)) / (norm.cdf(nSigClip) - norm.cdf(-nSigClip)) 

63 return 1.0 / np.sqrt(varFactor) 

64 

65 

66def calculateWeightedReducedChi2(measured, model, weightsMeasured, nData, nParsModel): 

67 """Calculate weighted reduced chi2. 

68 

69 Parameters 

70 ---------- 

71 measured : `list` 

72 List with measured data. 

73 model : `list` 

74 List with modeled data. 

75 weightsMeasured : `list` 

76 List with weights for the measured data. 

77 nData : `int` 

78 Number of data points. 

79 nParsModel : `int` 

80 Number of parameters in the model. 

81 

82 Returns 

83 ------- 

84 redWeightedChi2 : `float` 

85 Reduced weighted chi2. 

86 """ 

87 wRes = (measured - model)*weightsMeasured 

88 return ((wRes*wRes).sum())/(nData-nParsModel) 

89 

90 

91def makeMockFlats(expTime, gain=1.0, readNoiseElectrons=5, fluxElectrons=1000, 

92 randomSeedFlat1=1984, randomSeedFlat2=666, powerLawBfParams=[], 

93 expId1=0, expId2=1): 

94 """Create a pair or mock flats with isrMock. 

95 

96 Parameters 

97 ---------- 

98 expTime : `float` 

99 Exposure time of the flats. 

100 gain : `float`, optional 

101 Gain, in e/ADU. 

102 readNoiseElectrons : `float`, optional 

103 Read noise rms, in electrons. 

104 fluxElectrons : `float`, optional 

105 Flux of flats, in electrons per second. 

106 randomSeedFlat1 : `int`, optional 

107 Random seed for the normal distrubutions for the mean signal 

108 and noise (flat1). 

109 randomSeedFlat2 : `int`, optional 

110 Random seed for the normal distrubutions for the mean signal 

111 and noise (flat2). 

112 powerLawBfParams : `list`, optional 

113 Parameters for `galsim.cdmodel.PowerLawCD` to simulate the 

114 brightter-fatter effect. 

115 expId1 : `int`, optional 

116 Exposure ID for first flat. 

117 expId2 : `int`, optional 

118 Exposure ID for second flat. 

119 

120 Returns 

121 ------- 

122 flatExp1 : `lsst.afw.image.exposure.ExposureF` 

123 First exposure of flat field pair. 

124 flatExp2 : `lsst.afw.image.exposure.ExposureF` 

125 Second exposure of flat field pair. 

126 

127 Notes 

128 ----- 

129 The parameters of `galsim.cdmodel.PowerLawCD` are `n, r0, t0, rx, 

130 tx, r, t, alpha`. For more information about their meaning, see 

131 the Galsim documentation 

132 https://galsim-developers.github.io/GalSim/_build/html/_modules/galsim/cdmodel.html # noqa: W505 

133 and Gruen+15 (1501.02802). 

134 

135 Example: galsim.cdmodel.PowerLawCD(8, 1.1e-7, 1.1e-7, 1.0e-8, 

136 1.0e-8, 1.0e-9, 1.0e-9, 2.0) 

137 """ 

138 flatFlux = fluxElectrons # e/s 

139 flatMean = flatFlux*expTime # e 

140 readNoise = readNoiseElectrons # e 

141 

142 mockImageConfig = isrMock.IsrMock.ConfigClass() 

143 

144 mockImageConfig.flatDrop = 0.99999 

145 mockImageConfig.isTrimmed = True 

146 

147 flatExp1 = isrMock.FlatMock(config=mockImageConfig).run() 

148 flatExp2 = flatExp1.clone() 

149 (shapeY, shapeX) = flatExp1.getDimensions() 

150 flatWidth = np.sqrt(flatMean) 

151 

152 rng1 = np.random.RandomState(randomSeedFlat1) 

153 flatData1 = rng1.normal(flatMean, flatWidth, (shapeX, shapeY)) + rng1.normal(0.0, readNoise, 

154 (shapeX, shapeY)) 

155 rng2 = np.random.RandomState(randomSeedFlat2) 

156 flatData2 = rng2.normal(flatMean, flatWidth, (shapeX, shapeY)) + rng2.normal(0.0, readNoise, 

157 (shapeX, shapeY)) 

158 # Simulate BF with power law model in galsim 

159 if len(powerLawBfParams): 

160 if not len(powerLawBfParams) == 8: 

161 raise RuntimeError("Wrong number of parameters for `galsim.cdmodel.PowerLawCD`. " 

162 f"Expected 8; passed {len(powerLawBfParams)}.") 

163 cd = galsim.cdmodel.PowerLawCD(*powerLawBfParams) 

164 tempFlatData1 = galsim.Image(flatData1) 

165 temp2FlatData1 = cd.applyForward(tempFlatData1) 

166 

167 tempFlatData2 = galsim.Image(flatData2) 

168 temp2FlatData2 = cd.applyForward(tempFlatData2) 

169 

170 flatExp1.image.array[:] = temp2FlatData1.array/gain # ADU 

171 flatExp2.image.array[:] = temp2FlatData2.array/gain # ADU 

172 else: 

173 flatExp1.image.array[:] = flatData1/gain # ADU 

174 flatExp2.image.array[:] = flatData2/gain # ADU 

175 

176 visitInfoExp1 = lsst.afw.image.VisitInfo(exposureTime=expTime) 

177 visitInfoExp2 = lsst.afw.image.VisitInfo(exposureTime=expTime) 

178 

179 flatExp1.info.id = expId1 

180 flatExp1.getInfo().setVisitInfo(visitInfoExp1) 

181 flatExp2.info.id = expId2 

182 flatExp2.getInfo().setVisitInfo(visitInfoExp2) 

183 

184 return flatExp1, flatExp2 

185 

186 

187def irlsFit(initialParams, dataX, dataY, function, weightsY=None, weightType='Cauchy', scaleResidual=True): 

188 """Iteratively reweighted least squares fit. 

189 

190 This uses the `lsst.cp.pipe.utils.fitLeastSq`, but applies weights 

191 based on the Cauchy distribution by default. Other weight options 

192 are implemented. See e.g. Holland and Welsch, 1977, 

193 doi:10.1080/03610927708827533 

194 

195 Parameters 

196 ---------- 

197 initialParams : `list` [`float`] 

198 Starting parameters. 

199 dataX : `numpy.array`, (N,) 

200 Abscissa data. 

201 dataY : `numpy.array`, (N,) 

202 Ordinate data. 

203 function : callable 

204 Function to fit. 

205 weightsY : `numpy.array`, (N,) 

206 Weights to apply to the data. 

207 weightType : `str`, optional 

208 Type of weighting to use. One of Cauchy, Anderson, bisquare, 

209 box, Welsch, Huber, logistic, or Fair. 

210 scaleResidual : `bool`, optional 

211 If true, the residual is scaled by the sqrt of the Y values. 

212 

213 Returns 

214 ------- 

215 polyFit : `list` [`float`] 

216 Final best fit parameters. 

217 polyFitErr : `list` [`float`] 

218 Final errors on fit parameters. 

219 chiSq : `float` 

220 Reduced chi squared. 

221 weightsY : `list` [`float`] 

222 Final weights used for each point. 

223 

224 Raises 

225 ------ 

226 RuntimeError : 

227 Raised if an unknown weightType string is passed. 

228 """ 

229 if not weightsY: 

230 weightsY = np.ones_like(dataX) 

231 

232 polyFit, polyFitErr, chiSq = fitLeastSq(initialParams, dataX, dataY, function, weightsY=weightsY) 

233 for iteration in range(10): 

234 resid = np.abs(dataY - function(polyFit, dataX)) 

235 if scaleResidual: 

236 resid = resid / np.sqrt(dataY) 

237 if weightType == 'Cauchy': 

238 # Use Cauchy weighting. This is a soft weight. 

239 # At [2, 3, 5, 10] sigma, weights are [.59, .39, .19, .05]. 

240 Z = resid / 2.385 

241 weightsY = 1.0 / (1.0 + np.square(Z)) 

242 elif weightType == 'Anderson': 

243 # Anderson+1972 weighting. This is a hard weight. 

244 # At [2, 3, 5, 10] sigma, weights are [.67, .35, 0.0, 0.0]. 

245 Z = resid / (1.339 * np.pi) 

246 weightsY = np.where(Z < 1.0, np.sinc(Z), 0.0) 

247 elif weightType == 'bisquare': 

248 # Beaton and Tukey (1974) biweight. This is a hard weight. 

249 # At [2, 3, 5, 10] sigma, weights are [.81, .59, 0.0, 0.0]. 

250 Z = resid / 4.685 

251 weightsY = np.where(Z < 1.0, 1.0 - np.square(Z), 0.0) 

252 elif weightType == 'box': 

253 # Hinich and Talwar (1975). This is a hard weight. 

254 # At [2, 3, 5, 10] sigma, weights are [1.0, 0.0, 0.0, 0.0]. 

255 weightsY = np.where(resid < 2.795, 1.0, 0.0) 

256 elif weightType == 'Welsch': 

257 # Dennis and Welsch (1976). This is a hard weight. 

258 # At [2, 3, 5, 10] sigma, weights are [.64, .36, .06, 1e-5]. 

259 Z = resid / 2.985 

260 weightsY = np.exp(-1.0 * np.square(Z)) 

261 elif weightType == 'Huber': 

262 # Huber (1964) weighting. This is a soft weight. 

263 # At [2, 3, 5, 10] sigma, weights are [.67, .45, .27, .13]. 

264 Z = resid / 1.345 

265 weightsY = np.where(Z < 1.0, 1.0, 1 / Z) 

266 elif weightType == 'logistic': 

267 # Logistic weighting. This is a soft weight. 

268 # At [2, 3, 5, 10] sigma, weights are [.56, .40, .24, .12]. 

269 Z = resid / 1.205 

270 weightsY = np.tanh(Z) / Z 

271 elif weightType == 'Fair': 

272 # Fair (1974) weighting. This is a soft weight. 

273 # At [2, 3, 5, 10] sigma, weights are [.41, .32, .22, .12]. 

274 Z = resid / 1.4 

275 weightsY = (1.0 / (1.0 + (Z))) 

276 else: 

277 raise RuntimeError(f"Unknown weighting type: {weightType}") 

278 polyFit, polyFitErr, chiSq = fitLeastSq(initialParams, dataX, dataY, function, weightsY=weightsY) 

279 

280 return polyFit, polyFitErr, chiSq, weightsY 

281 

282 

283def fitLeastSq(initialParams, dataX, dataY, function, weightsY=None): 

284 """Do a fit and estimate the parameter errors using using 

285 scipy.optimize.leastq. 

286 

287 optimize.leastsq returns the fractional covariance matrix. To 

288 estimate the standard deviation of the fit parameters, multiply 

289 the entries of this matrix by the unweighted reduced chi squared 

290 and take the square root of the diagonal elements. 

291 

292 Parameters 

293 ---------- 

294 initialParams : `list` [`float`] 

295 initial values for fit parameters. For ptcFitType=POLYNOMIAL, 

296 its length determines the degree of the polynomial. 

297 dataX : `numpy.array`, (N,) 

298 Data in the abscissa axis. 

299 dataY : `numpy.array`, (N,) 

300 Data in the ordinate axis. 

301 function : callable object (function) 

302 Function to fit the data with. 

303 weightsY : `numpy.array`, (N,) 

304 Weights of the data in the ordinate axis. 

305 

306 Return 

307 ------ 

308 pFitSingleLeastSquares : `list` [`float`] 

309 List with fitted parameters. 

310 pErrSingleLeastSquares : `list` [`float`] 

311 List with errors for fitted parameters. 

312 

313 reducedChiSqSingleLeastSquares : `float` 

314 Reduced chi squared, unweighted if weightsY is not provided. 

315 """ 

316 if weightsY is None: 

317 weightsY = np.ones(len(dataX)) 

318 

319 def errFunc(p, x, y, weightsY=None): 

320 if weightsY is None: 

321 weightsY = np.ones(len(x)) 

322 return (function(p, x) - y)*weightsY 

323 

324 pFit, pCov, infoDict, errMessage, success = leastsq(errFunc, initialParams, 

325 args=(dataX, dataY, weightsY), full_output=1, 

326 epsfcn=0.0001) 

327 

328 if (len(dataY) > len(initialParams)) and pCov is not None: 

329 reducedChiSq = calculateWeightedReducedChi2(dataY, function(pFit, dataX), weightsY, len(dataY), 

330 len(initialParams)) 

331 pCov *= reducedChiSq 

332 else: 

333 pCov = np.zeros((len(initialParams), len(initialParams))) 

334 pCov[:, :] = np.nan 

335 reducedChiSq = np.nan 

336 

337 errorVec = [] 

338 for i in range(len(pFit)): 

339 errorVec.append(np.fabs(pCov[i][i])**0.5) 

340 

341 pFitSingleLeastSquares = pFit 

342 pErrSingleLeastSquares = np.array(errorVec) 

343 

344 return pFitSingleLeastSquares, pErrSingleLeastSquares, reducedChiSq 

345 

346 

347def fitBootstrap(initialParams, dataX, dataY, function, weightsY=None, confidenceSigma=1.): 

348 """Do a fit using least squares and bootstrap to estimate parameter errors. 

349 

350 The bootstrap error bars are calculated by fitting 100 random data sets. 

351 

352 Parameters 

353 ---------- 

354 initialParams : `list` [`float`] 

355 initial values for fit parameters. For ptcFitType=POLYNOMIAL, 

356 its length determines the degree of the polynomial. 

357 dataX : `numpy.array`, (N,) 

358 Data in the abscissa axis. 

359 dataY : `numpy.array`, (N,) 

360 Data in the ordinate axis. 

361 function : callable object (function) 

362 Function to fit the data with. 

363 weightsY : `numpy.array`, (N,), optional. 

364 Weights of the data in the ordinate axis. 

365 confidenceSigma : `float`, optional. 

366 Number of sigmas that determine confidence interval for the 

367 bootstrap errors. 

368 

369 Return 

370 ------ 

371 pFitBootstrap : `list` [`float`] 

372 List with fitted parameters. 

373 pErrBootstrap : `list` [`float`] 

374 List with errors for fitted parameters. 

375 reducedChiSqBootstrap : `float` 

376 Reduced chi squared, unweighted if weightsY is not provided. 

377 """ 

378 if weightsY is None: 

379 weightsY = np.ones(len(dataX)) 

380 

381 def errFunc(p, x, y, weightsY): 

382 if weightsY is None: 

383 weightsY = np.ones(len(x)) 

384 return (function(p, x) - y)*weightsY 

385 

386 # Fit first time 

387 pFit, _ = leastsq(errFunc, initialParams, args=(dataX, dataY, weightsY), full_output=0) 

388 

389 # Get the stdev of the residuals 

390 residuals = errFunc(pFit, dataX, dataY, weightsY) 

391 # 100 random data sets are generated and fitted 

392 pars = [] 

393 for i in range(100): 

394 randomDelta = np.random.normal(0., np.fabs(residuals), len(dataY)) 

395 randomDataY = dataY + randomDelta 

396 randomFit, _ = leastsq(errFunc, initialParams, 

397 args=(dataX, randomDataY, weightsY), full_output=0) 

398 pars.append(randomFit) 

399 pars = np.array(pars) 

400 meanPfit = np.mean(pars, 0) 

401 

402 # confidence interval for parameter estimates 

403 errPfit = confidenceSigma*np.std(pars, 0) 

404 pFitBootstrap = meanPfit 

405 pErrBootstrap = errPfit 

406 

407 reducedChiSq = calculateWeightedReducedChi2(dataY, function(pFitBootstrap, dataX), weightsY, len(dataY), 

408 len(initialParams)) 

409 return pFitBootstrap, pErrBootstrap, reducedChiSq 

410 

411 

412def funcPolynomial(pars, x): 

413 """Polynomial function definition 

414 Parameters 

415 ---------- 

416 params : `list` 

417 Polynomial coefficients. Its length determines the polynomial order. 

418 

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

420 Abscisa array. 

421 

422 Returns 

423 ------- 

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

425 Ordinate array after evaluating polynomial of order 

426 len(pars)-1 at `x`. 

427 """ 

428 return poly.polyval(x, [*pars]) 

429 

430 

431def funcAstier(pars, x): 

432 """Single brighter-fatter parameter model for PTC; Equation 16 of 

433 Astier+19. 

434 

435 Parameters 

436 ---------- 

437 params : `list` 

438 Parameters of the model: a00 (brightter-fatter), gain (e/ADU), 

439 and noise (e^2). 

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

441 Signal mu (ADU). 

442 

443 Returns 

444 ------- 

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

446 C_00 (variance) in ADU^2. 

447 """ 

448 a00, gain, noise = pars 

449 return 0.5/(a00*gain*gain)*(np.exp(2*a00*x*gain)-1) + noise/(gain*gain) # C_00 

450 

451 

452def arrangeFlatsByExpTime(exposureList, exposureIdList, log=None): 

453 """Arrange exposures by exposure time. 

454 

455 Parameters 

456 ---------- 

457 exposureList : `list` [`lsst.pipe.base.connections.DeferredDatasetRef`] 

458 Input list of exposure references. 

459 exposureIdList : `list` [`int`] 

460 List of exposure ids as obtained by dataId[`exposure`]. 

461 log : `lsst.utils.logging.LsstLogAdapter`, optional 

462 Log object. 

463 

464 Returns 

465 ------ 

466 flatsAtExpTime : `dict` [`float`, 

467 `list`[(`lsst.pipe.base.connections.DeferredDatasetRef`, 

468 `int`)]] 

469 Dictionary that groups references to flat-field exposures 

470 (and their IDs) that have the same exposure time (seconds). 

471 """ 

472 flatsAtExpTime = {} 

473 assert len(exposureList) == len(exposureIdList), "Different lengths for exp. list and exp. ID lists" 

474 for expRef, expId in zip(exposureList, exposureIdList): 

475 expTime = expRef.get(component='visitInfo').exposureTime 

476 if not np.isfinite(expTime) and log is not None: 

477 log.warning("Exposure %d has non-finite exposure time.", expId) 

478 listAtExpTime = flatsAtExpTime.setdefault(expTime, []) 

479 listAtExpTime.append((expRef, expId)) 

480 

481 return flatsAtExpTime 

482 

483 

484def arrangeFlatsByExpFlux(exposureList, exposureIdList, fluxKeyword, log=None): 

485 """Arrange exposures by exposure flux. 

486 

487 Parameters 

488 ---------- 

489 exposureList : `list` [`lsst.pipe.base.connections.DeferredDatasetRef`] 

490 Input list of exposure references. 

491 exposureIdList : `list` [`int`] 

492 List of exposure ids as obtained by dataId[`exposure`]. 

493 fluxKeyword : `str` 

494 Header keyword that contains the flux per exposure. 

495 log : `lsst.utils.logging.LsstLogAdapter`, optional 

496 Log object. 

497 

498 Returns 

499 ------- 

500 flatsAtFlux : `dict` [`float`, 

501 `list`[(`lsst.pipe.base.connections.DeferredDatasetRef`, 

502 `int`)]] 

503 Dictionary that groups references to flat-field exposures 

504 (and their IDs) that have the same flux. 

505 """ 

506 flatsAtExpFlux = {} 

507 assert len(exposureList) == len(exposureIdList), "Different lengths for exp. list and exp. ID lists" 

508 for expRef, expId in zip(exposureList, exposureIdList): 

509 # Get flux from header, assuming it is in the metadata. 

510 try: 

511 expFlux = expRef.get().getMetadata()[fluxKeyword] 

512 except KeyError: 

513 # If it's missing from the header, continue; it will 

514 # be caught and rejected when pairing exposures. 

515 expFlux = None 

516 if expFlux is None: 

517 if log is not None: 

518 log.warning("Exposure %d does not have valid header keyword %s.", expId, fluxKeyword) 

519 expFlux = np.nan 

520 listAtExpFlux = flatsAtExpFlux.setdefault(expFlux, []) 

521 listAtExpFlux.append((expRef, expId)) 

522 

523 return flatsAtExpFlux 

524 

525 

526def arrangeFlatsByExpId(exposureList, exposureIdList): 

527 """Arrange exposures by exposure ID. 

528 

529 There is no guarantee that this will properly group exposures, but 

530 allows a sequence of flats that have different illumination 

531 (despite having the same exposure time) to be processed. 

532 

533 Parameters 

534 ---------- 

535 exposureList : `list`[`lsst.pipe.base.connections.DeferredDatasetRef`] 

536 Input list of exposure references. 

537 exposureIdList : `list`[`int`] 

538 List of exposure ids as obtained by dataId[`exposure`]. 

539 

540 Returns 

541 ------ 

542 flatsAtExpId : `dict` [`float`, 

543 `list`[(`lsst.pipe.base.connections.DeferredDatasetRef`, 

544 `int`)]] 

545 Dictionary that groups references to flat-field exposures (and their 

546 IDs) sequentially by their exposure id. 

547 

548 Notes 

549 ----- 

550 

551 This algorithm sorts the input exposure references by their exposure 

552 id, and then assigns each pair of exposure references (exp_j, exp_{j+1}) 

553 to pair k, such that 2*k = j, where j is the python index of one of the 

554 exposure references (starting from zero). By checking for the IndexError 

555 while appending, we can ensure that there will only ever be fully 

556 populated pairs. 

557 """ 

558 flatsAtExpId = {} 

559 assert len(exposureList) == len(exposureIdList), "Different lengths for exp. list and exp. ID lists" 

560 # Sort exposures by expIds, which are in the second list `exposureIdList`. 

561 sortedExposures = sorted(zip(exposureList, exposureIdList), key=lambda pair: pair[1]) 

562 

563 for jPair, expTuple in enumerate(sortedExposures): 

564 if (jPair + 1) % 2: 

565 kPair = jPair // 2 

566 listAtExpId = flatsAtExpId.setdefault(kPair, []) 

567 try: 

568 listAtExpId.append(expTuple) 

569 listAtExpId.append(sortedExposures[jPair + 1]) 

570 except IndexError: 

571 pass 

572 

573 return flatsAtExpId 

574 

575 

576def extractCalibDate(calib): 

577 """Extract common calibration metadata values that will be written to 

578 output header. 

579 

580 Parameters 

581 ---------- 

582 calib : `lsst.afw.image.Exposure` or `lsst.ip.isr.IsrCalib` 

583 Calibration to pull date information from. 

584 

585 Returns 

586 ------- 

587 dateString : `str` 

588 Calibration creation date string to add to header. 

589 """ 

590 if hasattr(calib, "getMetadata"): 

591 if 'CALIB_CREATION_DATE' in calib.getMetadata(): 

592 return " ".join((calib.getMetadata().get("CALIB_CREATION_DATE", "Unknown"), 

593 calib.getMetadata().get("CALIB_CREATION_TIME", "Unknown"))) 

594 else: 

595 return " ".join((calib.getMetadata().get("CALIB_CREATE_DATE", "Unknown"), 

596 calib.getMetadata().get("CALIB_CREATE_TIME", "Unknown"))) 

597 else: 

598 return "Unknown Unknown" 

599 

600 

601class CovFastFourierTransform: 

602 """A class to compute (via FFT) the nearby pixels correlation function. 

603 

604 Implements appendix of Astier+19. 

605 

606 Parameters 

607 ---------- 

608 diff : `numpy.array` 

609 Image where to calculate the covariances (e.g., the difference 

610 image of two flats). 

611 w : `numpy.array` 

612 Weight image (mask): it should consist of 1's (good pixel) and 

613 0's (bad pixels). 

614 fftShape : `tuple` 

615 2d-tuple with the shape of the FFT 

616 maxRangeCov : `int` 

617 Maximum range for the covariances. 

618 """ 

619 

620 def __init__(self, diff, w, fftShape, maxRangeCov): 

621 # check that the zero padding implied by "fft_shape" 

622 # is large enough for the required correlation range 

623 assert fftShape[0] > diff.shape[0]+maxRangeCov+1 

624 assert fftShape[1] > diff.shape[1]+maxRangeCov+1 

625 # for some reason related to numpy.fft.rfftn, 

626 # the second dimension should be even, so 

627 if fftShape[1]%2 == 1: 

628 fftShape = (fftShape[0], fftShape[1]+1) 

629 tIm = np.fft.rfft2(diff*w, fftShape) 

630 tMask = np.fft.rfft2(w, fftShape) 

631 # sum of "squares" 

632 self.pCov = np.fft.irfft2(tIm*tIm.conjugate()) 

633 # sum of values 

634 self.pMean = np.fft.irfft2(tIm*tMask.conjugate()) 

635 # number of w!=0 pixels. 

636 self.pCount = np.fft.irfft2(tMask*tMask.conjugate()) 

637 

638 def cov(self, dx, dy): 

639 """Covariance for dx,dy averaged with dx,-dy if both non zero. 

640 

641 Implements appendix of Astier+19. 

642 

643 Parameters 

644 ---------- 

645 dx : `int` 

646 Lag in x 

647 dy : `int` 

648 Lag in y 

649 

650 Returns 

651 ------- 

652 0.5*(cov1+cov2) : `float` 

653 Covariance at (dx, dy) lag 

654 npix1+npix2 : `int` 

655 Number of pixels used in covariance calculation. 

656 

657 Raises 

658 ------ 

659 ValueError if number of pixels for a given lag is 0. 

660 """ 

661 # compensate rounding errors 

662 nPix1 = int(round(self.pCount[dy, dx])) 

663 if nPix1 == 0: 

664 raise ValueError(f"Could not compute covariance term {dy}, {dx}, as there are no good pixels.") 

665 cov1 = self.pCov[dy, dx]/nPix1-self.pMean[dy, dx]*self.pMean[-dy, -dx]/(nPix1*nPix1) 

666 if (dx == 0 or dy == 0): 

667 return cov1, nPix1 

668 nPix2 = int(round(self.pCount[-dy, dx])) 

669 if nPix2 == 0: 

670 raise ValueError("Could not compute covariance term {dy}, {dx} as there are no good pixels.") 

671 cov2 = self.pCov[-dy, dx]/nPix2-self.pMean[-dy, dx]*self.pMean[dy, -dx]/(nPix2*nPix2) 

672 return 0.5*(cov1+cov2), nPix1+nPix2 

673 

674 def reportCovFastFourierTransform(self, maxRange): 

675 """Produce a list of tuples with covariances. 

676 

677 Implements appendix of Astier+19. 

678 

679 Parameters 

680 ---------- 

681 maxRange : `int` 

682 Maximum range of covariances. 

683 

684 Returns 

685 ------- 

686 tupleVec : `list` 

687 List with covariance tuples. 

688 """ 

689 tupleVec = [] 

690 # (dy,dx) = (0,0) has to be first 

691 for dy in range(maxRange+1): 

692 for dx in range(maxRange+1): 

693 cov, npix = self.cov(dx, dy) 

694 if (dx == 0 and dy == 0): 

695 var = cov 

696 tupleVec.append((dx, dy, var, cov, npix)) 

697 return tupleVec 

698 

699 

700def getFitDataFromCovariances(i, j, mu, fullCov, fullCovModel, fullCovSqrtWeights, gain=1.0, 

701 divideByMu=False, returnMasked=False): 

702 """Get measured signal and covariance, cov model, weigths, and mask at 

703 covariance lag (i, j). 

704 

705 Parameters 

706 ---------- 

707 i : `int` 

708 Lag for covariance matrix. 

709 j : `int` 

710 Lag for covariance matrix. 

711 mu : `list` 

712 Mean signal values. 

713 fullCov : `list` of `numpy.array` 

714 Measured covariance matrices at each mean signal level in mu. 

715 fullCovSqrtWeights : `list` of `numpy.array` 

716 List of square root of measured covariances at each mean 

717 signal level in mu. 

718 fullCovModel : `list` of `numpy.array` 

719 List of modeled covariances at each mean signal level in mu. 

720 gain : `float`, optional 

721 Gain, in e-/ADU. If other than 1.0 (default), the returned 

722 quantities will be in electrons or powers of electrons. 

723 divideByMu : `bool`, optional 

724 Divide returned covariance, model, and weights by the mean 

725 signal mu? 

726 returnMasked : `bool`, optional 

727 Use mask (based on weights) in returned arrays (mu, 

728 covariance, and model)? 

729 

730 Returns 

731 ------- 

732 mu : `numpy.array` 

733 list of signal values at (i, j). 

734 covariance : `numpy.array` 

735 Covariance at (i, j) at each mean signal mu value (fullCov[:, i, j]). 

736 covarianceModel : `numpy.array` 

737 Covariance model at (i, j). 

738 weights : `numpy.array` 

739 Weights at (i, j). 

740 maskFromWeights : `numpy.array`, optional 

741 Boolean mask of the covariance at (i,j), where the weights 

742 differ from 0. 

743 """ 

744 mu = np.array(mu) 

745 fullCov = np.array(fullCov) 

746 fullCovModel = np.array(fullCovModel) 

747 fullCovSqrtWeights = np.array(fullCovSqrtWeights) 

748 covariance = fullCov[:, i, j]*(gain**2) 

749 covarianceModel = fullCovModel[:, i, j]*(gain**2) 

750 weights = fullCovSqrtWeights[:, i, j]/(gain**2) 

751 

752 maskFromWeights = weights != 0 

753 if returnMasked: 

754 weights = weights[maskFromWeights] 

755 covarianceModel = covarianceModel[maskFromWeights] 

756 mu = mu[maskFromWeights] 

757 covariance = covariance[maskFromWeights] 

758 

759 if divideByMu: 

760 covariance /= mu 

761 covarianceModel /= mu 

762 weights *= mu 

763 return mu, covariance, covarianceModel, weights, maskFromWeights 

764 

765 

766def symmetrize(inputArray): 

767 """ Copy array over 4 quadrants prior to convolution. 

768 

769 Parameters 

770 ---------- 

771 inputarray : `numpy.array` 

772 Input array to symmetrize. 

773 

774 Returns 

775 ------- 

776 aSym : `numpy.array` 

777 Symmetrized array. 

778 """ 

779 targetShape = list(inputArray.shape) 

780 r1, r2 = inputArray.shape[-1], inputArray.shape[-2] 

781 targetShape[-1] = 2*r1-1 

782 targetShape[-2] = 2*r2-1 

783 aSym = np.ndarray(tuple(targetShape)) 

784 aSym[..., r2-1:, r1-1:] = inputArray 

785 aSym[..., r2-1:, r1-1::-1] = inputArray 

786 aSym[..., r2-1::-1, r1-1::-1] = inputArray 

787 aSym[..., r2-1::-1, r1-1:] = inputArray 

788 

789 return aSym 

790 

791 

792def ddict2dict(d): 

793 """Convert nested default dictionaries to regular dictionaries. 

794 

795 This is needed to prevent yaml persistence issues. 

796 

797 Parameters 

798 ---------- 

799 d : `defaultdict` 

800 A possibly nested set of `defaultdict`. 

801 

802 Returns 

803 ------- 

804 dict : `dict` 

805 A possibly nested set of `dict`. 

806 """ 

807 for k, v in d.items(): 

808 if isinstance(v, dict): 

809 d[k] = ddict2dict(v) 

810 return dict(d) 

811 

812 

813class Pol2D: 

814 """2D Polynomial Regression. 

815 

816 Parameters 

817 ---------- 

818 x : numpy.ndarray 

819 Input array for the x-coordinate. 

820 y : numpy.ndarray 

821 Input array for the y-coordinate. 

822 z : numpy.ndarray 

823 Input array for the dependent variable. 

824 order : int 

825 Order of the polynomial. 

826 w : numpy.ndarray, optional 

827 Weight array for weighted regression. Default is None. 

828 

829 Notes 

830 ----- 

831 Ported from by https://gitlab.in2p3.fr/astier/bfptc P. Astier. 

832 

833 Example: 

834 >>> x = np.array([1, 2, 3]) 

835 >>> y = np.array([4, 5, 6]) 

836 >>> z = np.array([7, 8, 9]) 

837 >>> order = 2 

838 >>> poly_reg = Pol2D(x, y, z, order) 

839 >>> result = poly_reg.eval(2.5, 5.5) 

840 """ 

841 def __init__(self, x, y, z, order, w=None): 

842 """ 

843 orderx : `int` 

844 Effective order in the x-direction. 

845 ordery : `int` 

846 Effective order in the y-direction. 

847 coeff : `numpy.ndarray` 

848 Coefficients of the polynomial regression. 

849 """ 

850 self.orderx = min(order, x.shape[0] - 1) 

851 self.ordery = min(order, x.shape[1] - 1) 

852 G = self.monomials(x.ravel(), y.ravel()) 

853 if w is None: 

854 self.coeff, _, rank, _ = np.linalg.lstsq(G, z.ravel(), rcond=None) 

855 else: 

856 self.coeff, _, rank, _ = np.linalg.lstsq((w.ravel() * G.T).T, z.ravel() * w.ravel(), rcond=None) 

857 

858 def monomials(self, x, y): 

859 """ 

860 Generate the monomials matrix for the given x and y. 

861 

862 Parameters 

863 ---------- 

864 x : numpy.ndarray 

865 Input array for the x-coordinate. 

866 y : numpy.ndarray 

867 Input array for the y-coordinate. 

868 

869 Returns 

870 ------- 

871 G : numpy.ndarray 

872 Monomials matrix. 

873 """ 

874 ncols = (self.orderx + 1) * (self.ordery + 1) 

875 G = np.zeros(x.shape + (ncols,)) 

876 ij = itertools.product(range(self.orderx + 1), range(self.ordery + 1)) 

877 for k, (i, j) in enumerate(ij): 

878 G[..., k] = x**i * y**j 

879 return G 

880 

881 def eval(self, x, y): 

882 """ 

883 Evaluate the polynomial at the given x and y coordinates. 

884 

885 Parameters 

886 ---------- 

887 x : `float` 

888 x-coordinate for evaluation. 

889 y : `float` 

890 y-coordinate for evaluation. 

891 

892 Returns 

893 ------- 

894 result : `float` 

895 Result of the polynomial evaluation. 

896 """ 

897 G = self.monomials(x, y) 

898 return np.dot(G, self.coeff) 

899 

900 

901class AstierSplineLinearityFitter: 

902 """Class to fit the Astier spline linearity model. 

903 

904 This is a spline fit with photodiode data based on a model 

905 from Pierre Astier, referenced in June 2023 from 

906 https://me.lsst.eu/astier/bot/7224D/model_nonlin.py 

907 

908 This model fits a spline with (optional) nuisance parameters 

909 to allow for different linearity coefficients with different 

910 photodiode settings. The minimization is a least-squares 

911 fit with the residual of 

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

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

914 mean; D_j is the photo-diode measurement corresponding to 

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

916 which is over index j as it is allowed to 

917 be different based on different photodiode settings (e.g. 

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

919 leaks (and is only fit if fit_offset=True). In 

920 addition, if config.doSplineFitTemperature is True then 

921 the fit will adjust mu such that 

922 mu = mu_input*(1 + alpha*(temperature_scaled)) 

923 and temperature_scaled = T - T_ref. 

924 

925 The fit has additional constraints to ensure that the spline 

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

927 condition so that the average of the spline over the full 

928 range is 0. The normalization ensures that the spline only 

929 fits deviations from linearity, rather than the linear 

930 function itself which is degenerate with the gain. 

931 

932 Parameters 

933 ---------- 

934 nodes : `np.ndarray` (N,) 

935 Array of spline node locations. 

936 grouping_values : `np.ndarray` (M,) 

937 Array of values to group values for different proportionality 

938 constants (e.g. CCOBCURR). 

939 pd : `np.ndarray` (M,) 

940 Array of photodiode measurements. 

941 mu : `np.ndarray` (M,) 

942 Array of flat mean values. 

943 mask : `np.ndarray` (M,), optional 

944 Input mask (True is good point, False is bad point). 

945 log : `logging.logger`, optional 

946 Logger object to use for logging. 

947 fit_offset : `bool`, optional 

948 Fit a constant offset to allow for light leaks? 

949 fit_weights : `bool`, optional 

950 Fit for the weight parameters? 

951 weight_pars_start : `list` [ `float` ] 

952 Iterable of 2 weight parameters for weighed fit. These will 

953 be used as input if the weight parameters are not fit. 

954 fit_temperature : `bool`, optional 

955 Fit for temperature scaling? 

956 temperature_scaled : `np.ndarray` (M,), optional 

957 Input scaled temperature values (T - T_ref). 

958 """ 

959 def __init__( 

960 self, 

961 nodes, 

962 grouping_values, 

963 pd, 

964 mu, 

965 mask=None, 

966 log=None, 

967 fit_offset=True, 

968 fit_weights=False, 

969 weight_pars_start=[1.0, 0.0], 

970 fit_temperature=False, 

971 temperature_scaled=None, 

972 ): 

973 self._pd = pd 

974 self._mu = mu 

975 self._grouping_values = grouping_values 

976 self.log = log if log else logging.getLogger(__name__) 

977 self._fit_offset = fit_offset 

978 self._fit_weights = fit_weights 

979 self._weight_pars_start = weight_pars_start 

980 self._fit_temperature = fit_temperature 

981 

982 self._nodes = nodes 

983 if nodes[0] != 0.0: 

984 raise ValueError("First node must be 0.0") 

985 if not np.all(np.diff(nodes) > 0): 

986 raise ValueError("Nodes must be sorted with no repeats.") 

987 

988 # Find the group indices. 

989 u_group_values = np.unique(self._grouping_values) 

990 self.ngroup = len(u_group_values) 

991 

992 self.group_indices = [] 

993 for i in range(self.ngroup): 

994 self.group_indices.append(np.where(self._grouping_values == u_group_values[i])[0]) 

995 

996 # Weight values. Outliers will be set to 0. 

997 if mask is None: 

998 _mask = np.ones(len(mu), dtype=np.bool_) 

999 else: 

1000 _mask = mask 

1001 self._w = self.compute_weights(self._weight_pars_start, self._mu, _mask) 

1002 

1003 if temperature_scaled is None: 

1004 temperature_scaled = np.zeros(len(self._mu)) 

1005 else: 

1006 if len(np.atleast_1d(temperature_scaled)) != len(self._mu): 

1007 raise ValueError("temperature_scaled must be the same length as input mu.") 

1008 self._temperature_scaled = temperature_scaled 

1009 

1010 # Values to regularize spline fit. 

1011 self._x_regularize = np.linspace(0.0, self._mu[self.mask].max(), 100) 

1012 

1013 # Set up the indices for the fit parameters. 

1014 self.par_indices = { 

1015 "values": np.arange(len(self._nodes)), 

1016 "groups": len(self._nodes) + np.arange(self.ngroup), 

1017 "offset": np.zeros(0, dtype=np.int64), 

1018 "weight_pars": np.zeros(0, dtype=np.int64), 

1019 "temperature_coeff": np.zeros(0, dtype=np.int64), 

1020 } 

1021 if self._fit_offset: 

1022 self.par_indices["offset"] = np.arange(1) + ( 

1023 len(self.par_indices["values"]) 

1024 + len(self.par_indices["groups"]) 

1025 ) 

1026 if self._fit_weights: 

1027 self.par_indices["weight_pars"] = np.arange(2) + ( 

1028 len(self.par_indices["values"]) 

1029 + len(self.par_indices["groups"]) 

1030 + len(self.par_indices["offset"]) 

1031 ) 

1032 if self._fit_temperature: 

1033 self.par_indices["temperature_coeff"] = np.arange(1) + ( 

1034 len(self.par_indices["values"]) 

1035 + len(self.par_indices["groups"]) 

1036 + len(self.par_indices["offset"]) 

1037 + len(self.par_indices["weight_pars"]) 

1038 ) 

1039 

1040 @staticmethod 

1041 def compute_weights(weight_pars, mu, mask): 

1042 w = 1./np.sqrt(weight_pars[0]**2. + weight_pars[1]**2./mu) 

1043 w[~mask] = 0.0 

1044 

1045 return w 

1046 

1047 def estimate_p0(self): 

1048 """Estimate initial fit parameters. 

1049 

1050 Returns 

1051 ------- 

1052 p0 : `np.ndarray` 

1053 Parameter array, with spline values (one for each node) followed 

1054 by proportionality constants (one for each group); one extra 

1055 for the offset O (if fit_offset was set to True); two extra 

1056 for the weights (if fit_weights was set to True); and one 

1057 extra for the temperature coefficient (if fit_temperature was 

1058 set to True). 

1059 """ 

1060 npt = (len(self.par_indices["values"]) 

1061 + len(self.par_indices["groups"]) 

1062 + len(self.par_indices["offset"]) 

1063 + len(self.par_indices["weight_pars"]) 

1064 + len(self.par_indices["temperature_coeff"])) 

1065 p0 = np.zeros(npt) 

1066 

1067 # Do a simple linear fit and set all the constants to this. 

1068 linfit = np.polyfit(self._pd[self.mask], self._mu[self.mask], 1) 

1069 p0[self.par_indices["groups"]] = linfit[0] 

1070 

1071 # Look at the residuals... 

1072 ratio_model = self.compute_ratio_model( 

1073 self._nodes, 

1074 self.group_indices, 

1075 self.par_indices, 

1076 p0, 

1077 self._pd, 

1078 self._mu, 

1079 self._temperature_scaled, 

1080 ) 

1081 # ...and adjust the linear parameters accordingly. 

1082 p0[self.par_indices["groups"]] *= np.median(ratio_model[self.mask]) 

1083 

1084 # Re-compute the residuals. 

1085 ratio_model2 = self.compute_ratio_model( 

1086 self._nodes, 

1087 self.group_indices, 

1088 self.par_indices, 

1089 p0, 

1090 self._pd, 

1091 self._mu, 

1092 self._temperature_scaled, 

1093 ) 

1094 

1095 # And compute a first guess of the spline nodes. 

1096 bins = np.searchsorted(self._nodes, self._mu[self.mask]) 

1097 tot_arr = np.zeros(len(self._nodes)) 

1098 n_arr = np.zeros(len(self._nodes), dtype=int) 

1099 np.add.at(tot_arr, bins, ratio_model2[self.mask]) 

1100 np.add.at(n_arr, bins, 1) 

1101 

1102 ratio = np.ones(len(self._nodes)) 

1103 ratio[n_arr > 0] = tot_arr[n_arr > 0]/n_arr[n_arr > 0] 

1104 ratio[0] = 1.0 

1105 p0[self.par_indices["values"]] = (ratio - 1) * self._nodes 

1106 

1107 if self._fit_offset: 

1108 p0[self.par_indices["offset"]] = 0.0 

1109 

1110 if self._fit_weights: 

1111 p0[self.par_indices["weight_pars"]] = self._weight_pars_start 

1112 

1113 return p0 

1114 

1115 @staticmethod 

1116 def compute_ratio_model( 

1117 nodes, 

1118 group_indices, 

1119 par_indices, 

1120 pars, 

1121 pd, 

1122 mu, 

1123 temperature_scaled, 

1124 return_spline=False, 

1125 ): 

1126 """Compute the ratio model values. 

1127 

1128 Parameters 

1129 ---------- 

1130 nodes : `np.ndarray` (M,) 

1131 Array of node positions. 

1132 group_indices : `list` [`np.ndarray`] 

1133 List of group indices, one array for each group. 

1134 par_indices : `dict` 

1135 Dictionary showing which indices in the pars belong to 

1136 each set of fit values. 

1137 pars : `np.ndarray` 

1138 Parameter array, with spline values (one for each node) followed 

1139 by proportionality constants (one for each group); one extra 

1140 for the offset O (if fit_offset was set to True); two extra 

1141 for the weights (if fit_weights was set to True); and one 

1142 extra for the temperature coefficient (if fit_temperature was 

1143 set to True). 

1144 pd : `np.ndarray` (N,) 

1145 Array of photodiode measurements. 

1146 mu : `np.ndarray` (N,) 

1147 Array of flat means. 

1148 temperature_scaled : `np.ndarray` (N,) 

1149 Array of scaled temperature values. 

1150 return_spline : `bool`, optional 

1151 Return the spline interpolation as well as the model ratios? 

1152 

1153 Returns 

1154 ------- 

1155 ratio_models : `np.ndarray` (N,) 

1156 Model ratio, (mu_i - S(mu_i) - O)/(k_j * D_i) 

1157 spl : `lsst.afw.math.thing` 

1158 Spline interpolator (returned if return_spline=True). 

1159 """ 

1160 spl = lsst.afw.math.makeInterpolate( 

1161 nodes, 

1162 pars[par_indices["values"]], 

1163 lsst.afw.math.stringToInterpStyle("AKIMA_SPLINE"), 

1164 ) 

1165 

1166 # Check if we want to do just the left or both with temp scale. 

1167 if len(par_indices["temperature_coeff"]) == 1: 

1168 mu_corr = mu*(1. + pars[par_indices["temperature_coeff"]]*temperature_scaled) 

1169 else: 

1170 mu_corr = mu 

1171 

1172 numerator = mu_corr - spl.interpolate(mu_corr) 

1173 if len(par_indices["offset"]) == 1: 

1174 numerator -= pars[par_indices["offset"][0]] 

1175 denominator = pd.copy() 

1176 ngroup = len(group_indices) 

1177 kj = pars[par_indices["groups"]] 

1178 for j in range(ngroup): 

1179 denominator[group_indices[j]] *= kj[j] 

1180 

1181 if return_spline: 

1182 return numerator / denominator, spl 

1183 else: 

1184 return numerator / denominator 

1185 

1186 def fit(self, p0, min_iter=3, max_iter=20, max_rejection_per_iteration=5, n_sigma_clip=5.0): 

1187 """ 

1188 Perform iterative fit for linear + spline model with offsets. 

1189 

1190 Parameters 

1191 ---------- 

1192 p0 : `np.ndarray` 

1193 Initial parameter array, with spline values (one for each node) 

1194 followed by proportionality constants (one for each group); one 

1195 extra for the offset O (if fit_offset was set to True); two extra 

1196 for the weights (if fit_weights was set to True); and one 

1197 extra for the temperature coefficient (if fit_temperature was 

1198 set to True). 

1199 min_iter : `int`, optional 

1200 Minimum number of fit iterations. 

1201 max_iter : `int`, optional 

1202 Maximum number of fit iterations. 

1203 max_rejection_per_iteration : `int`, optional 

1204 Maximum number of points to reject per iteration. 

1205 n_sigma_clip : `float`, optional 

1206 Number of sigma to do clipping in each iteration. 

1207 """ 

1208 init_params = p0 

1209 for k in range(max_iter): 

1210 params, cov_params, _, msg, ierr = leastsq( 

1211 self, 

1212 init_params, 

1213 full_output=True, 

1214 ftol=1e-5, 

1215 maxfev=12000, 

1216 ) 

1217 init_params = params.copy() 

1218 

1219 # We need to cut off the constraints at the end (there are more 

1220 # residuals than data points.) 

1221 res = self(params)[: len(self._w)] 

1222 std_res = median_abs_deviation(res[self.good_points], scale="normal") 

1223 sample = len(self.good_points) 

1224 

1225 # We don't want to reject too many outliers at once. 

1226 if sample > max_rejection_per_iteration: 

1227 sres = np.sort(np.abs(res)) 

1228 cut = max(sres[-max_rejection_per_iteration], std_res*n_sigma_clip) 

1229 else: 

1230 cut = std_res*n_sigma_clip 

1231 

1232 outliers = np.abs(res) > cut 

1233 self._w[outliers] = 0 

1234 if outliers.sum() != 0: 

1235 self.log.info( 

1236 "After iteration %d there are %d outliers (of %d).", 

1237 k, 

1238 outliers.sum(), 

1239 sample, 

1240 ) 

1241 elif k >= min_iter: 

1242 self.log.info("After iteration %d there are no more outliers.", k) 

1243 break 

1244 

1245 return params 

1246 

1247 @property 

1248 def mask(self): 

1249 return (self._w > 0) 

1250 

1251 @property 

1252 def good_points(self): 

1253 return self.mask.nonzero()[0] 

1254 

1255 def compute_chisq_dof(self, pars): 

1256 """Compute the chi-squared per degree of freedom for a set of pars. 

1257 

1258 Parameters 

1259 ---------- 

1260 pars : `np.ndarray` 

1261 Parameter array. 

1262 

1263 Returns 

1264 ------- 

1265 chisq_dof : `float` 

1266 Chi-squared per degree of freedom. 

1267 """ 

1268 resids = self(pars)[0: len(self.mask)] 

1269 chisq = np.sum(resids[self.mask]**2.) 

1270 dof = self.mask.sum() - self.ngroup 

1271 if self._fit_temperature: 

1272 dof -= 1 

1273 if self._fit_offset: 

1274 dof -= 1 

1275 if self._fit_weights: 

1276 dof -= 2 

1277 

1278 return chisq/dof 

1279 

1280 def __call__(self, pars): 

1281 

1282 ratio_model, spl = self.compute_ratio_model( 

1283 self._nodes, 

1284 self.group_indices, 

1285 self.par_indices, 

1286 pars, 

1287 self._pd, 

1288 self._mu, 

1289 self._temperature_scaled, 

1290 return_spline=True, 

1291 ) 

1292 

1293 _mask = self.mask 

1294 # Update the weights if we are fitting them. 

1295 if self._fit_weights: 

1296 self._w = self.compute_weights(pars[self.par_indices["weight_pars"]], self._mu, _mask) 

1297 resid = self._w*(ratio_model - 1.0) 

1298 

1299 # Ensure masked points have 0 residual. 

1300 resid[~_mask] = 0.0 

1301 

1302 constraint = [1e3 * np.mean(spl.interpolate(self._x_regularize))] 

1303 # 0 should transform to 0 

1304 constraint.append(spl.interpolate(0)*1e10) 

1305 # Use a Jeffreys prior on the weight if we are fitting it. 

1306 if self._fit_weights: 

1307 # This factor ensures that log(fact * w) is negative. 

1308 fact = 1e-3 / self._w.max() 

1309 # We only add non-zero weights to the constraint array. 

1310 log_w = np.sqrt(-2.*np.log(fact*self._w[self._w > 0])) 

1311 constraint = np.hstack([constraint, log_w]) 

1312 

1313 return np.hstack([resid, constraint])