Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1# This file is part of cp_pipe. 

2# 

3# Developed for the LSST Data Management System. 

4# This product includes software developed by the LSST Project 

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

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

7# for details of code ownership. 

8# 

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

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

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

12# (at your option) any later version. 

13# 

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

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

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

17# GNU General Public License for more details. 

18# 

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

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

21 

22import numpy as np 

23import copy 

24from scipy.signal import fftconvolve 

25from scipy.optimize import leastsq 

26from .astierCovFitParameters import FitParameters 

27 

28import lsst.log as lsstLog 

29 

30__all__ = ["CovFit"] 

31 

32 

33def makeCovArray(inputTuple, maxRangeFromTuple=8): 

34 """Make covariances array from tuple. 

35 

36 Parameters 

37 ---------- 

38 inputTuple: `numpy.ndarray` 

39 Structured array with rows with at least 

40 (mu, afwVar, cov, var, i, j, npix), where: 

41 

42 mu : 0.5*(m1 + m2), where: 

43 mu1: mean value of flat1 

44 mu2: mean value of flat2 

45 afwVar: variance of difference flat, calculated with afw 

46 cov: covariance value at lag(i, j) 

47 var: variance(covariance value at lag(0, 0)) 

48 i: lag dimension 

49 j: lag dimension 

50 npix: number of pixels used for covariance calculation. 

51 

52 maxRangeFromTuple: `int` 

53 Maximum range to select from tuple. 

54 

55 Returns 

56 ------- 

57 cov: `numpy.array` 

58 Covariance arrays, indexed by mean signal mu. 

59 

60 vCov: `numpy.array` 

61 Variance arrays, indexed by mean signal mu. 

62 

63 muVals: `numpy.array` 

64 List of mean signal values. 

65 

66 Notes 

67 ----- 

68 

69 The input tuple should contain the following rows: 

70 (mu, cov, var, i, j, npix), with one entry per lag, and image pair. 

71 Different lags(i.e. different i and j) from the same 

72 image pair have the same values of mu1 and mu2. When i==j==0, cov 

73 = var. 

74 

75 If the input tuple contains several video channels, one should 

76 select the data of a given channel *before* entering this 

77 routine, as well as apply(e.g.) saturation cuts. 

78 

79 The routine returns cov[k_mu, j, i], vcov[(same indices)], and mu[k] 

80 where the first index of cov matches the one in mu. 

81 

82 This routine implements the loss of variance due to 

83 clipping cuts when measuring variances and covariance, but this should happen inside 

84 the measurement code, where the cuts are readily available. 

85 

86 """ 

87 if maxRangeFromTuple is not None: 

88 cut = (inputTuple['i'] < maxRangeFromTuple) & (inputTuple['j'] < maxRangeFromTuple) 

89 cutTuple = inputTuple[cut] 

90 else: 

91 cutTuple = inputTuple 

92 # increasing mu order, so that we can group measurements with the same mu 

93 muTemp = cutTuple['mu'] 

94 ind = np.argsort(muTemp) 

95 

96 cutTuple = cutTuple[ind] 

97 # should group measurements on the same image pairs(same average) 

98 mu = cutTuple['mu'] 

99 xx = np.hstack(([mu[0]], mu)) 

100 delta = xx[1:] - xx[:-1] 

101 steps, = np.where(delta > 0) 

102 ind = np.zeros_like(mu, dtype=int) 

103 ind[steps] = 1 

104 ind = np.cumsum(ind) # this acts as an image pair index. 

105 # now fill the 3-d cov array(and variance) 

106 muVals = np.array(np.unique(mu)) 

107 i = cutTuple['i'].astype(int) 

108 j = cutTuple['j'].astype(int) 

109 c = 0.5*cutTuple['cov'] 

110 n = cutTuple['npix'] 

111 v = 0.5*cutTuple['var'] 

112 # book and fill 

113 cov = np.ndarray((len(muVals), np.max(i)+1, np.max(j)+1)) 

114 var = np.zeros_like(cov) 

115 cov[ind, i, j] = c 

116 var[ind, i, j] = v**2/n 

117 var[:, 0, 0] *= 2 # var(v) = 2*v**2/N 

118 # compensate for loss of variance and covariance due to outlier elimination(sigma clipping) 

119 # when computing variances(cut to 4 sigma): 1 per mill for variances and twice as 

120 # much for covariances: 

121 fact = 1.0 # 1.00107 

122 cov *= fact*fact 

123 cov[:, 0, 0] /= fact 

124 

125 return cov, var, muVals 

126 

127 

128def symmetrize(inputArray): 

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

130 

131 Parameters 

132 ---------- 

133 inputarray: `numpy.array` 

134 Input array to symmetrize. 

135 

136 Returns 

137 ------- 

138 aSym: `numpy.array` 

139 Symmetrized array. 

140 """ 

141 

142 targetShape = list(inputArray.shape) 

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

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

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

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

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

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

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

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

151 

152 return aSym 

153 

154 

155class CovFit: 

156 """A class to fit the models in Astier+19 to flat covariances. 

157 

158 This code implements the model(and the fit thereof) described in 

159 Astier+19: https://arxiv.org/pdf/1905.08677.pdf 

160 

161 Parameters 

162 ---------- 

163 meanSignals : `list`[`float`] 

164 List with means of the difference image of two flats, 

165 for a particular amplifier in the detector. 

166 

167 covariances : `list`[`numpy.array`] 

168 List with 2D covariance arrays at a given mean signal. 

169 

170 covsSqrtWeights : `list`[`numpy.array`] 

171 List with 2D arrays with weights from `vcov as defined in 

172 `makeCovArray`: weight = 1/sqrt(vcov). 

173 

174 maxRangeFromTuple: `int`, optional 

175 Maximum range to select from tuple. 

176 

177 meanSignalMask: `list`[`bool`], optional 

178 Mask of mean signal 1D array. Use all entries if empty. 

179 """ 

180 

181 def __init__(self, meanSignals, covariances, covsSqrtWeights, maxRangeFromTuple=8, meanSignalsMask=[]): 

182 assert (len(meanSignals) == len(covariances)) 

183 assert (len(covariances) == len(covsSqrtWeights)) 

184 if len(meanSignalsMask) == 0: 

185 meanSignalsMask = np.repeat(True, len(meanSignals)) 

186 self.mu = meanSignals[meanSignalsMask] 

187 self.cov = np.nan_to_num(covariances)[meanSignalsMask] 

188 # make it nan safe, replacing nan's with 0 in weights 

189 self.sqrtW = np.nan_to_num(covsSqrtWeights)[meanSignalsMask] 

190 self.r = maxRangeFromTuple 

191 self.logger = lsstLog.Log.getDefaultLogger() 

192 

193 def copy(self): 

194 """Make a copy of params""" 

195 cop = copy.deepcopy(self) 

196 # deepcopy does not work for FitParameters. 

197 if hasattr(self, 'params'): 

198 cop.params = self.params.copy() 

199 return cop 

200 

201 def initFit(self): 

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

203 the full fit close to the solution. 

204 """ 

205 # number of parameters for 'a' array. 

206 lenA = self.r*self.r 

207 # define parameters: c corresponds to a*b in Astier+19 (Eq. 20). 

208 self.params = FitParameters([('a', lenA), ('c', lenA), ('noise', lenA), ('gain', 1)]) 

209 self.params['gain'] = 1. 

210 # c=0 in a first go. 

211 self.params['c'].fix(val=0.) 

212 # plumbing: extract stuff from the parameter structure 

213 a = self.params['a'].full.reshape(self.r, self.r) 

214 noise = self.params['noise'].full.reshape(self.r, self.r) 

215 gain = self.params['gain'].full[0] 

216 

217 # iterate the fit to account for higher orders 

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

219 # stop when it increases 

220 oldChi2 = 1e30 

221 for _ in range(5): 

222 model = np.nan_to_num(self.evalCovModel()) # this computes the full model. 

223 # loop on lags 

224 for i in range(self.r): 

225 for j in range(self.r): 

226 # fit a parabola for a given lag 

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

228 2, w=self.sqrtW[:, i, j]) 

229 # model equation(Eq. 20) in Astier+19: 

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

231 noise[i, j] += parsFit[2]*gain*gain 

232 if(i + j == 0): 

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

234 self.params['gain'].full[0] = gain 

235 chi2 = self.chi2() 

236 if chi2 > oldChi2: 

237 break 

238 oldChi2 = chi2 

239 

240 return 

241 

242 def getParamValues(self): 

243 """Return an array of free parameter values (it is a copy).""" 

244 return self.params.free + 0. 

245 

246 def setParamValues(self, p): 

247 """Set parameter values.""" 

248 self.params.free = p 

249 return 

250 

251 def evalCovModel(self, mu=None): 

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

253 

254 Parameters 

255 ---------- 

256 mu: `numpy.array`, optional 

257 List of mean signals. 

258 

259 Returns 

260 ------- 

261 covModel: `numpy.array` 

262 Covariances model. 

263 

264 Notes 

265 ----- 

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

267 

268 Returns cov[Nmu, self.r, self.r]. The variance for the PTC is cov[:, 0, 0]. 

269 mu and cov are in ADUs and ADUs squared. To use electrons for both, 

270 the gain should be set to 1. This routine implements the model in Astier+19 (1905.08677). 

271 

272 The parameters of the full model for C_ij(mu) ("C_ij" and "mu" in ADU^2 and ADU, respectively) 

273 in Astier+19 (Eq. 20) are: 

274 

275 "a" coefficients (r by r matrix), units: 1/e 

276 "b" coefficients (r by r matrix), units: 1/e 

277 noise matrix (r by r matrix), units: e^2 

278 gain, units: e/ADU 

279 

280 "b" appears in Eq. 20 only through the "ab" combination, which is defined in this code as "c=ab". 

281 """ 

282 sa = (self.r, self.r) 

283 a = self.params['a'].full.reshape(sa) 

284 c = self.params['c'].full.reshape(sa) 

285 gain = self.params['gain'].full[0] 

286 noise = self.params['noise'].full.reshape(sa) 

287 # pad a with zeros and symmetrize 

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

289 aEnlarged[0:sa[0], 0:sa[1]] = a 

290 aSym = symmetrize(aEnlarged) 

291 # pad c with zeros and symmetrize 

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

293 cEnlarged[0:sa[0], 0:sa[1]] = c 

294 cSym = symmetrize(cEnlarged) 

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

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

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

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

299 range = self.r 

300 a1 = a[np.newaxis, :, :] 

301 a2 = a2[np.newaxis, xc:xc + range, yc:yc + range] 

302 a3 = a3[np.newaxis, xc:xc + range, yc:yc + range] 

303 ac = ac[np.newaxis, xc:xc + range, yc:yc + range] 

304 c1 = c[np.newaxis, ::] 

305 if mu is None: 

306 mu = self.mu 

307 # assumes that mu is 1d 

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

309 # c(=a*b in Astier+19) also has a contribution to the last term, that is absent for now. 

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

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

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

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

314 

315 return covModel 

316 

317 def getA(self): 

318 """'a' matrix from Astier+19(e.g., Eq. 20)""" 

319 return self.params['a'].full.reshape(self.r, self.r) 

320 

321 def getB(self): 

322 """'b' matrix from Astier+19(e.g., Eq. 20)""" 

323 return self.params['c'].full.reshape(self.r, self.r)/self.getA() 

324 

325 def getC(self): 

326 """'c'='ab' matrix from Astier+19(e.g., Eq. 20)""" 

327 return np.array(self.params['c'].full.reshape(self.r, self.r)) 

328 

329 def _getCovParams(self, what): 

330 """Get covariance matrix of parameters from fit""" 

331 indices = self.params[what].indexof() 

332 i1 = indices[:, np.newaxis] 

333 i2 = indices[np.newaxis, :] 

334 if self.covParams is not None: 

335 covp = self.covParams[i1, i2] 

336 else: 

337 covp = None 

338 return covp 

339 

340 def getACov(self): 

341 """Get covariance matrix of "a" coefficients from fit""" 

342 if self._getCovParams('a') is not None: 

343 cova = self._getCovParams('a').reshape((self.r, self.r, self.r, self.r)) 

344 else: 

345 cova = None 

346 return cova 

347 

348 def getASig(self): 

349 """Square root of diagonal of the parameter covariance of the fitted "a" matrix""" 

350 if self._getCovParams('a') is not None: 

351 sigA = np.sqrt(self._getCovParams('a').diagonal()).reshape((self.r, self.r)) 

352 else: 

353 sigA = None 

354 return sigA 

355 

356 def getBCov(self): 

357 """Get covariance matrix of "a" coefficients from fit 

358 b = c /a 

359 """ 

360 covb = self._getCovParams('c') 

361 aval = self.getA().flatten() 

362 factor = np.outer(aval, aval) 

363 covb /= factor 

364 return covb.reshape((self.r, self.r, self.r, self.r)) 

365 

366 def getCCov(self): 

367 """Get covariance matrix of "c" coefficients from fit""" 

368 cova = self._getCovParams('c') 

369 return cova.reshape((self.r, self.r, self.r, self.r)) 

370 

371 def getGainErr(self): 

372 """Get error on fitted gain parameter""" 

373 if self._getCovParams('gain') is not None: 

374 gainErr = np.sqrt(self._getCovParams('gain')[0][0]) 

375 else: 

376 gainErr = 0.0 

377 return gainErr 

378 

379 def getNoiseCov(self): 

380 """Get covariances of noise matrix from fit""" 

381 covNoise = self._getCovParams('noise') 

382 return covNoise.reshape((self.r, self.r, self.r, self.r)) 

383 

384 def getNoiseSig(self): 

385 """Square root of diagonal of the parameter covariance of the fitted "noise" matrix""" 

386 if self._getCovParams('noise') is not None: 

387 covNoise = self._getCovParams('noise') 

388 noise = np.sqrt(covNoise.diagonal()).reshape((self.r, self.r)) 

389 else: 

390 noise = None 

391 return noise 

392 

393 def getGain(self): 

394 """Get gain (e/ADU)""" 

395 return self.params['gain'].full[0] 

396 

397 def getRon(self): 

398 """Get readout noise (e^2)""" 

399 return self.params['noise'].full[0] 

400 

401 def getRonErr(self): 

402 """Get error on readout noise parameter""" 

403 ronSqrt = np.sqrt(np.fabs(self.getRon())) 

404 if self.getNoiseSig() is not None: 

405 noiseSigma = self.getNoiseSig()[0][0] 

406 ronErr = 0.5*(noiseSigma/np.fabs(self.getRon()))*ronSqrt 

407 else: 

408 ronErr = np.nan 

409 return ronErr 

410 

411 def getNoise(self): 

412 """Get noise matrix""" 

413 return self.params['noise'].full.reshape(self.r, self.r) 

414 

415 def getMaskCov(self, i, j): 

416 """Get mask of Cov[i,j]""" 

417 weights = self.sqrtW[:, i, j] 

418 mask = weights != 0 

419 return mask 

420 

421 def setAandB(self, a, b): 

422 """Set "a" and "b" coeffcients forfull Astier+19 model (Eq. 20). "c=a*b".""" 

423 self.params['a'].full = a.flatten() 

424 self.params['c'].full = a.flatten()*b.flatten() 

425 return 

426 

427 def chi2(self): 

428 """Calculate weighted chi2 of full-model fit.""" 

429 return(self.weightedRes()**2).sum() 

430 

431 def weightedRes(self, params=None): 

432 """Weighted residuals. 

433 

434 Notes 

435 ----- 

436 To be used via: 

437 c = CovFit(meanSignals, covariances, covsSqrtWeights) 

438 c.initFit() 

439 coeffs, cov, _, mesg, ierr = leastsq(c.weightedRes, c.getParamValues(), full_output=True) 

440 """ 

441 if params is not None: 

442 self.setParamValues(params) 

443 covModel = np.nan_to_num(self.evalCovModel()) 

444 weightedRes = (covModel-self.cov)*self.sqrtW 

445 

446 return weightedRes.flatten() 

447 

448 def fitFullModel(self, pInit=None): 

449 """Fit measured covariances to full model in Astier+19 (Eq. 20) 

450 

451 Parameters 

452 ---------- 

453 pInit : `list` 

454 Initial parameters of the fit. 

455 len(pInit) = #entries(a) + #entries(c) + #entries(noise) + 1 

456 len(pInit) = r^2 + r^2 + r^2 + 1, where "r" is the maximum lag considered for the 

457 covariances calculation, and the extra "1" is the gain. 

458 If "b" is 0, then "c" is 0, and len(pInit) will have r^2 fewer entries. 

459 

460 Returns 

461 ------- 

462 params : `np.array` 

463 Fit parameters (see "Notes" below). 

464 

465 Notes 

466 ----- 

467 The parameters of the full model for C_ij(mu) ("C_ij" and "mu" in ADU^2 and ADU, respectively) 

468 in Astier+19 (Eq. 20) are: 

469 

470 "a" coefficients (r by r matrix), units: 1/e 

471 "b" coefficients (r by r matrix), units: 1/e 

472 noise matrix (r by r matrix), units: e^2 

473 gain, units: e/ADU 

474 

475 "b" appears in Eq. 20 only through the "ab" combination, which is defined in this code as "c=ab". 

476 """ 

477 

478 if pInit is None: 

479 pInit = self.getParamValues() 

480 params, paramsCov, _, mesg, ierr = leastsq(self.weightedRes, pInit, full_output=True) 

481 self.covParams = paramsCov 

482 

483 return params 

484 

485 def getFitData(self, i, j, divideByMu=False, unitsElectrons=False, returnMasked=False): 

486 """Get measured signal and covariance, cov model, weigths, and mask at covariance lag (i, j). 

487 

488 Parameters 

489 --------- 

490 i: `int` 

491 Lag for covariance matrix. 

492 

493 j: `int` 

494 Lag for covariance matrix. 

495 

496 divideByMu: `bool`, optional 

497 Divide covariance, model, and weights by signal mu? 

498 

499 unitsElectrons : `bool`, optional 

500 mu, covariance, and model are in ADU (or powers of ADU) If tthis parameter is true, these are 

501 multiplied by the adequte factors of the gain to return quantities in electrons 

502 (or powers of electrons). 

503 

504 returnMasked : `bool`, optional 

505 Use mask (based on weights) in returned arrays (mu, covariance, and model)? 

506 

507 Returns 

508 ------- 

509 mu: `numpy.array` 

510 list of signal values (mu). 

511 

512 covariance: `numpy.array` 

513 Covariance arrays, indexed by mean signal mu (self.cov[:, i, j]). 

514 

515 covarianceModel: `numpy.array` 

516 Covariance model (model). 

517 

518 weights: `numpy.array` 

519 Weights (self.sqrtW) 

520 

521 mask : `numpy.array`, optional 

522 Boolean mask of the covariance at (i,j). 

523 

524 Notes 

525 ----- 

526 Using a CovFit object, selects from (i, j) and returns 

527 mu*gain, self.cov[:, i, j]*gain**2 model*gain**2, and self.sqrtW/gain**2 

528 in electrons or ADU if unitsElectrons=False. 

529 """ 

530 if unitsElectrons: 

531 gain = self.getGain() 

532 else: 

533 gain = 1.0 

534 

535 mu = self.mu*gain 

536 covariance = self.cov[:, i, j]*(gain**2) 

537 covarianceModel = self.evalCovModel()[:, i, j]*(gain**2) 

538 weights = self.sqrtW[:, i, j]/(gain**2) 

539 

540 # select data used for the fit: 

541 mask = self.getMaskCov(i, j) 

542 if returnMasked: 

543 weights = weights[mask] 

544 covarianceModel = covarianceModel[mask] 

545 mu = mu[mask] 

546 covariance = covariance[mask] 

547 

548 if divideByMu: 

549 covariance /= mu 

550 covarianceModel /= mu 

551 weights *= mu 

552 

553 return mu, covariance, covarianceModel, weights, mask