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 

119 return cov, var, muVals 

120 

121 

122def symmetrize(inputArray): 

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

124 

125 Parameters 

126 ---------- 

127 inputarray: `numpy.array` 

128 Input array to symmetrize. 

129 

130 Returns 

131 ------- 

132 aSym: `numpy.array` 

133 Symmetrized array. 

134 """ 

135 

136 targetShape = list(inputArray.shape) 

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

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

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

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

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

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

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

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

145 

146 return aSym 

147 

148 

149class CovFit: 

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

151 

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

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

154 

155 Parameters 

156 ---------- 

157 meanSignals : `list`[`float`] 

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

159 for a particular amplifier in the detector. 

160 

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

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

163 

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

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

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

167 

168 maxRangeFromTuple: `int`, optional 

169 Maximum range to select from tuple. 

170 

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

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

173 """ 

174 

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

176 assert (len(meanSignals) == len(covariances)) 

177 assert (len(covariances) == len(covsSqrtWeights)) 

178 if len(meanSignalsMask) == 0: 

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

180 self.mu = meanSignals[meanSignalsMask] 

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

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

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

184 self.r = maxRangeFromTuple 

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

186 

187 def copy(self): 

188 """Make a copy of params""" 

189 cop = copy.deepcopy(self) 

190 # deepcopy does not work for FitParameters. 

191 if hasattr(self, 'params'): 

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

193 return cop 

194 

195 def initFit(self): 

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

197 the full fit close to the solution. 

198 """ 

199 # number of parameters for 'a' array. 

200 lenA = self.r*self.r 

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

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

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

204 # c=0 in a first go. 

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

206 # plumbing: extract stuff from the parameter structure 

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

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

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

210 

211 # iterate the fit to account for higher orders 

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

213 # stop when it increases 

214 oldChi2 = 1e30 

215 for _ in range(5): 

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

217 # loop on lags 

218 for i in range(self.r): 

219 for j in range(self.r): 

220 # fit a parabola for a given lag 

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

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

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

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

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

226 if(i + j == 0): 

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

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

229 chi2 = self.chi2() 

230 if chi2 > oldChi2: 

231 break 

232 oldChi2 = chi2 

233 

234 return 

235 

236 def getParamValues(self): 

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

238 return self.params.free + 0. 

239 

240 def setParamValues(self, p): 

241 """Set parameter values.""" 

242 self.params.free = p 

243 return 

244 

245 def evalCovModel(self, mu=None): 

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

247 

248 Parameters 

249 ---------- 

250 mu: `numpy.array`, optional 

251 List of mean signals. 

252 

253 Returns 

254 ------- 

255 covModel: `numpy.array` 

256 Covariances model. 

257 

258 Notes 

259 ----- 

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

261 

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

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

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

265 

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

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

268 

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

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

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

272 gain, units: e/ADU 

273 

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

275 """ 

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

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

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

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

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

281 # pad a with zeros and symmetrize 

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

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

284 aSym = symmetrize(aEnlarged) 

285 # pad c with zeros and symmetrize 

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

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

288 cSym = symmetrize(cEnlarged) 

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

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

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

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

293 range = self.r 

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

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

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

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

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

299 if mu is None: 

300 mu = self.mu 

301 # assumes that mu is 1d 

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

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

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

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

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

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

308 

309 return covModel 

310 

311 def getA(self): 

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

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

314 

315 def getB(self): 

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

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

318 

319 def getC(self): 

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

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

322 

323 def _getCovParams(self, what): 

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

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

326 i1 = indices[:, np.newaxis] 

327 i2 = indices[np.newaxis, :] 

328 if self.covParams is not None: 

329 covp = self.covParams[i1, i2] 

330 else: 

331 covp = None 

332 return covp 

333 

334 def getACov(self): 

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

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

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

338 else: 

339 cova = None 

340 return cova 

341 

342 def getASig(self): 

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

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

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

346 else: 

347 sigA = None 

348 return sigA 

349 

350 def getBCov(self): 

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

352 b = c /a 

353 """ 

354 covb = self._getCovParams('c') 

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

356 factor = np.outer(aval, aval) 

357 covb /= factor 

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

359 

360 def getCCov(self): 

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

362 cova = self._getCovParams('c') 

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

364 

365 def getGainErr(self): 

366 """Get error on fitted gain parameter""" 

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

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

369 else: 

370 gainErr = 0.0 

371 return gainErr 

372 

373 def getNoiseCov(self): 

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

375 covNoise = self._getCovParams('noise') 

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

377 

378 def getNoiseSig(self): 

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

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

381 covNoise = self._getCovParams('noise') 

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

383 else: 

384 noise = None 

385 return noise 

386 

387 def getGain(self): 

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

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

390 

391 def getRon(self): 

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

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

394 

395 def getRonErr(self): 

396 """Get error on readout noise parameter""" 

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

398 if self.getNoiseSig() is not None: 

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

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

401 else: 

402 ronErr = np.nan 

403 return ronErr 

404 

405 def getNoise(self): 

406 """Get noise matrix""" 

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

408 

409 def getMaskCov(self, i, j): 

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

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

412 mask = weights != 0 

413 return mask 

414 

415 def setAandB(self, a, b): 

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

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

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

419 return 

420 

421 def chi2(self): 

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

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

424 

425 def weightedRes(self, params=None): 

426 """Weighted residuals. 

427 

428 Notes 

429 ----- 

430 To be used via: 

431 c = CovFit(meanSignals, covariances, covsSqrtWeights) 

432 c.initFit() 

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

434 """ 

435 if params is not None: 

436 self.setParamValues(params) 

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

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

439 

440 return weightedRes.flatten() 

441 

442 def fitFullModel(self, pInit=None): 

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

444 

445 Parameters 

446 ---------- 

447 pInit : `list` 

448 Initial parameters of the fit. 

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

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

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

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

453 

454 Returns 

455 ------- 

456 params : `np.array` 

457 Fit parameters (see "Notes" below). 

458 

459 Notes 

460 ----- 

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

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

463 

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

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

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

467 gain, units: e/ADU 

468 

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

470 """ 

471 

472 if pInit is None: 

473 pInit = self.getParamValues() 

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

475 self.covParams = paramsCov 

476 

477 return params 

478 

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

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

481 

482 Parameters 

483 --------- 

484 i: `int` 

485 Lag for covariance matrix. 

486 

487 j: `int` 

488 Lag for covariance matrix. 

489 

490 divideByMu: `bool`, optional 

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

492 

493 unitsElectrons : `bool`, optional 

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

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

496 (or powers of electrons). 

497 

498 returnMasked : `bool`, optional 

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

500 

501 Returns 

502 ------- 

503 mu: `numpy.array` 

504 list of signal values (mu). 

505 

506 covariance: `numpy.array` 

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

508 

509 covarianceModel: `numpy.array` 

510 Covariance model (model). 

511 

512 weights: `numpy.array` 

513 Weights (self.sqrtW) 

514 

515 mask : `numpy.array`, optional 

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

517 

518 Notes 

519 ----- 

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

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

522 in electrons or ADU if unitsElectrons=False. 

523 """ 

524 if unitsElectrons: 

525 gain = self.getGain() 

526 else: 

527 gain = 1.0 

528 

529 mu = self.mu*gain 

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

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

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

533 

534 # select data used for the fit: 

535 mask = self.getMaskCov(i, j) 

536 if returnMasked: 

537 weights = weights[mask] 

538 covarianceModel = covarianceModel[mask] 

539 mu = mu[mask] 

540 covariance = covariance[mask] 

541 

542 if divideByMu: 

543 covariance /= mu 

544 covarianceModel /= mu 

545 weights *= mu 

546 

547 return mu, covariance, covarianceModel, weights, mask