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 

24import itertools 

25from scipy.stats import median_abs_deviation as mad 

26from scipy.signal import fftconvolve 

27from scipy.optimize import leastsq 

28from .astierCovFitParameters import FitParameters 

29 

30 

31__all__ = ["CovFit"] 

32 

33 

34def computeApproximateAcoeffs(covModel, muEl, gain): 

35 """Compute the "a" coefficients of the Antilogus+14 (1402.0725) model as in 

36 Guyonnet+15 (1501.01577, eq. 16, the slope of cov/var at a given flux mu in electrons). 

37 

38 Eq. 16 of 1501.01577 is an approximation to the more complete model in Astier+19 (1905.08677). 

39 

40 Parameters 

41 --------- 

42 covModel : `list` 

43 Covariance model from Eq. 20 in Astier+19. 

44 

45 muEl : `np.array` 

46 Mean signal in electrons 

47 

48 gain : `float` 

49 Gain in e-/ADU. 

50 

51 Returns 

52 ------- 

53 aCoeffsOld: `numpy.array` 

54 Slope of cov/var at a given flux mu in electrons. 

55 

56 Notes 

57 ----- 

58 Returns the "a" array, computed this way, to be compared to the actual a_array from the full model 

59 (fit.geA()). 

60 """ 

61 covModel = np.array(covModel) 

62 var = covModel[0, 0, 0] # ADU^2 

63 # For a result in electrons^-1, we have to use mu in electrons. 

64 return covModel[0, :, :]/(var*muEl) 

65 

66 

67def makeCovArray(inputTuple, maxRangeFromTuple=8): 

68 """Make covariances array from tuple. 

69 

70 Parameters 

71 ---------- 

72 inputTuple: `numpy.recarray` 

73 Recarray with rows with at least (mu, cov, var, i, j, npix), where: 

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

75 mu1: mean value of flat1 

76 mu2: mean value of flat2 

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

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

79 i: lag dimension 

80 j: lag dimension 

81 npix: number of pixels used for covariance calculation. 

82 

83 maxRangeFromTuple: `int` 

84 Maximum range to select from tuple. 

85 

86 Returns 

87 ------- 

88 cov: `numpy.array` 

89 Covariance arrays, indexed by mean signal mu. 

90 

91 vCov: `numpy.array` 

92 Variance arrays, indexed by mean signal mu. 

93 

94 muVals: `numpy.array` 

95 List of mean signal values. 

96 

97 Notes 

98 ----- 

99 

100 The input tuple should contain the following rows: 

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

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

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

104 = var. 

105 

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

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

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

109 

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

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

112 

113 This routine implements the loss of variance due to 

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

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

116 

117 """ 

118 if maxRangeFromTuple is not None: 

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

120 cutTuple = inputTuple[cut] 

121 else: 

122 cutTuple = inputTuple 

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

124 muTemp = cutTuple['mu'] 

125 ind = np.argsort(muTemp) 

126 

127 cutTuple = cutTuple[ind] 

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

129 mu = cutTuple['mu'] 

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

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

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

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

134 ind[steps] = 1 

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

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

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

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

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

140 c = 0.5*cutTuple['cov'] 

141 n = cutTuple['npix'] 

142 v = 0.5*cutTuple['var'] 

143 # book and fill 

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

145 var = np.zeros_like(cov) 

146 cov[ind, i, j] = c 

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

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

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

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

151 # much for covariances: 

152 fact = 1.0 # 1.00107 

153 cov *= fact*fact 

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

155 

156 return cov, var, muVals 

157 

158 

159def symmetrize(inputArray): 

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

161 

162 Parameters 

163 ---------- 

164 inputarray: `numpy.array` 

165 Input array to symmetrize. 

166 

167 Returns 

168 ------- 

169 aSym: `numpy.array` 

170 Symmetrized array. 

171 """ 

172 

173 targetShape = list(inputArray.shape) 

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

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

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

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

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

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

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

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

182 

183 return aSym 

184 

185 

186class Pol2d: 

187 """A class to calculate 2D polynomials""" 

188 

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

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

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

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

193 if w is None: 

194 self.coeff, _, rank, _ = np.linalg.lstsq(G, z.ravel()) 

195 else: 

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

197 

198 def monomials(self, x, y): 

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

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

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

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

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

204 

205 return G 

206 

207 def eval(self, x, y): 

208 G = self.monomials(x, y) 

209 

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

211 

212 

213class CovFit: 

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

215 

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

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

218 For the time being it uses as input a numpy recarray (tuple with named tags) which 

219 contains one row per covariance and per pair: see the routine makeCovArray. 

220 

221 Parameters 

222 ---------- 

223 inputTuple: `numpy.recarray` 

224 Tuple with at least (mu, cov, var, i, j, npix), where: 

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

226 mu1: mean value of flat1 

227 mu2: mean value of flat2 

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

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

230 i: lag dimension 

231 j: lag dimension 

232 npix: number of pixels used for covariance calculation. 

233 

234 maxRangeFromTuple: `int` 

235 Maximum range to select from tuple. 

236 """ 

237 

238 def __init__(self, inputTuple, maxRangeFromTuple=8): 

239 self.cov, self.vcov, self.mu = makeCovArray(inputTuple, maxRangeFromTuple) 

240 self.sqrtW = 1./np.sqrt(self.vcov) 

241 self.r = self.cov.shape[1] 

242 

243 def subtractDistantOffset(self, maxLag=8, startLag=5, polDegree=1): 

244 """Subtract a background/offset to the measured covariances. 

245 

246 Parameters 

247 --------- 

248 maxLag: `int` 

249 Maximum lag considered 

250 

251 startLag: `int` 

252 First lag from where to start the offset subtraction. 

253 

254 polDegree: `int` 

255 Degree of 2D polynomial to fit to covariance to define offse to be subtracted. 

256 """ 

257 assert(startLag < self.r) 

258 for k in range(len(self.mu)): 

259 # Make a copy because it is going to be altered 

260 w = self.sqrtW[k, ...] + 0. 

261 sh = w.shape 

262 i, j = np.meshgrid(range(sh[0]), range(sh[1]), indexing='ij') 

263 # kill the core for the fit 

264 w[:startLag, :startLag] = 0 

265 poly = Pol2d(i, j, self.cov[k, ...], polDegree+1, w=w) 

266 back = poly.eval(i, j) 

267 self.cov[k, ...] -= back 

268 self.r = maxLag 

269 self.cov = self.cov[:, :maxLag, :maxLag] 

270 self.vcov = self.vcov[:, :maxLag, :maxLag] 

271 self.sqrtW = self.sqrtW[:, :maxLag, :maxLag] 

272 

273 return 

274 

275 def setMaxMu(self, maxMu): 

276 """Select signal level based on max average signal in ADU""" 

277 # mus are sorted at construction 

278 index = self.mu < maxMu 

279 k = index.sum() 

280 self.mu = self.mu[:k] 

281 self.cov = self.cov[:k, ...] 

282 self.vcov = self.vcov[:k, ...] 

283 self.sqrtW = self.sqrtW[:k, ...] 

284 

285 return 

286 

287 def setMaxMuElectrons(self, maxMuEl): 

288 """Select signal level based on max average signal in electrons""" 

289 g = self.getGain() 

290 kill = (self.mu*g > maxMuEl) 

291 self.sqrtW[kill, :, :] = 0 

292 

293 return 

294 

295 def copy(self): 

296 """Make a copy of params""" 

297 cop = copy.deepcopy(self) 

298 # deepcopy does not work for FitParameters. 

299 if hasattr(self, 'params'): 

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

301 return cop 

302 

303 def initFit(self): 

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

305 the full fit close to the solution. 

306 """ 

307 # number of parameters for 'a' array. 

308 lenA = self.r*self.r 

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

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

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

312 # c=0 in a first go. 

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

314 # plumbing: extract stuff from the parameter structure 

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

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

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

318 

319 # iterate the fit to account for higher orders 

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

321 # stop when it increases 

322 oldChi2 = 1e30 

323 for _ in range(5): 

324 model = self.evalCovModel() # this computes the full model. 

325 # loop on lags 

326 for i in range(self.r): 

327 for j in range(self.r): 

328 # fit a parabola for a given lag 

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

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

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

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

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

334 if(i + j == 0): 

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

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

337 chi2 = self.chi2() 

338 if chi2 > oldChi2: 

339 break 

340 oldChi2 = chi2 

341 

342 return 

343 

344 def getParamValues(self): 

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

346 return self.params.free + 0. 

347 

348 def setParamValues(self, p): 

349 """Set parameter values.""" 

350 self.params.free = p 

351 return 

352 

353 def evalCovModel(self, mu=None): 

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

355 

356 Parameters 

357 ---------- 

358 mu: `numpy.array`, optional 

359 List of mean signals. 

360 

361 Returns 

362 ------- 

363 covModel: `numpy.array` 

364 Covariances model. 

365 

366 Notes 

367 ----- 

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

369 

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

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

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

373 

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

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

376 

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

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

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

380 gain, units: e/ADU 

381 

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

383 """ 

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

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

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

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

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

389 # pad a with zeros and symmetrize 

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

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

392 aSym = symmetrize(aEnlarged) 

393 # pad c with zeros and symmetrize 

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

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

396 cSym = symmetrize(cEnlarged) 

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

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

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

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

401 range = self.r 

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

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

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

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

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

407 if mu is None: 

408 mu = self.mu 

409 # assumes that mu is 1d 

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

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

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

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

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

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

416 

417 return covModel 

418 

419 def getA(self): 

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

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

422 

423 def getB(self): 

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

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

426 

427 def getC(self): 

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

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

430 

431 def _getCovParams(self, what): 

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

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

434 i1 = indices[:, np.newaxis] 

435 i2 = indices[np.newaxis, :] 

436 if self.covParams is not None: 

437 covp = self.covParams[i1, i2] 

438 else: 

439 covp = None 

440 return covp 

441 

442 def getACov(self): 

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

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

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

446 else: 

447 cova = None 

448 return cova 

449 

450 def getASig(self): 

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

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

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

454 else: 

455 sigA = None 

456 return sigA 

457 

458 def getBCov(self): 

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

460 b = c /a 

461 """ 

462 covb = self._getCovParams('c') 

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

464 factor = np.outer(aval, aval) 

465 covb /= factor 

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

467 

468 def getCCov(self): 

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

470 cova = self._getCovParams('c') 

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

472 

473 def getGainErr(self): 

474 """Get error on fitted gain parameter""" 

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

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

477 else: 

478 gainErr = 0.0 

479 return gainErr 

480 

481 def getNoiseCov(self): 

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

483 covNoise = self._getCovParams('noise') 

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

485 

486 def getNoiseSig(self): 

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

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

489 covNoise = self._getCovParams('noise') 

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

491 else: 

492 noise = None 

493 return noise 

494 

495 def getGain(self): 

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

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

498 

499 def getRon(self): 

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

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

502 

503 def getRonErr(self): 

504 """Get error on readout noise parameter""" 

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

506 if self.getNoiseSig() is not None: 

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

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

509 else: 

510 ronErr = np.nan 

511 return ronErr 

512 

513 def getNoise(self): 

514 """Get noise matrix""" 

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

516 

517 def getMaskCov(self, i, j): 

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

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

520 mask = weights != 0 

521 return mask 

522 

523 def setAandB(self, a, b): 

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

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

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

527 return 

528 

529 def chi2(self): 

530 """Calculate weighte chi2 of full-model fit.""" 

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

532 

533 def wres(self, params=None): 

534 """To be used in weightedRes""" 

535 if params is not None: 

536 self.setParamValues(params) 

537 covModel = self.evalCovModel() 

538 return((covModel-self.cov)*self.sqrtW) 

539 

540 def weightedRes(self, params=None): 

541 """Weighted residuas. 

542 

543 Notes 

544 ----- 

545 To be used via: 

546 c = CovFit(nt) 

547 c.initFit() 

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

549 """ 

550 return self.wres(params).flatten() 

551 

552 def fitFullModel(self, pInit=None, nSigma=5.0, maxFitIter=3): 

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

554 

555 Parameters 

556 ---------- 

557 pInit : `list` 

558 Initial parameters of the fit. 

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

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

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

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

563 

564 nSigma : `float`, optional 

565 Sigma cut to get rid of outliers. 

566 

567 maxFitIter : `int`, optional 

568 Number of iterations for full model fit. 

569 

570 Returns 

571 ------- 

572 params : `np.array` 

573 Fit parameters (see "Notes" below). 

574 

575 Notes 

576 ----- 

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

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

579 

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

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

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

583 gain, units: e/ADU 

584 

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

586 """ 

587 

588 if pInit is None: 

589 pInit = self.getParamValues() 

590 nOutliers = 1 

591 counter = 0 

592 while nOutliers != 0: 

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

594 wres = self.weightedRes(params) 

595 # Do not count the outliers as significant 

596 sig = mad(wres[wres != 0], scale='normal') 

597 mask = (np.abs(wres) > (nSigma*sig)) 

598 self.sqrtW.flat[mask] = 0 # flatten makes a copy 

599 nOutliers = mask.sum() 

600 counter += 1 

601 if counter == maxFitIter: 

602 break 

603 

604 if ierr not in [1, 2, 3, 4]: 

605 raise RuntimeError("Minimization failed: " + mesg) 

606 self.covParams = paramsCov 

607 

608 return params 

609 

610 def ndof(self): 

611 """Number of degrees of freedom 

612 

613 Returns 

614 ------- 

615 mask.sum() - len(self.params.free): `int` 

616 Number of usable pixels - number of parameters of fit. 

617 """ 

618 mask = self.sqrtW != 0 

619 

620 return mask.sum() - len(self.params.free) 

621 

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

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

624 

625 Parameters 

626 --------- 

627 i: `int` 

628 Lag for covariance matrix. 

629 

630 j: `int` 

631 Lag for covariance matrix. 

632 

633 divideByMu: `bool`, optional 

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

635 

636 unitsElectrons : `bool`, optional 

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

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

639 (or powers of electrons). 

640 

641 returnMasked : `bool`, optional 

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

643 

644 Returns 

645 ------- 

646 mu: `numpy.array` 

647 list of signal values (mu). 

648 

649 covariance: `numpy.array` 

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

651 

652 covarianceModel: `numpy.array` 

653 Covariance model (model). 

654 

655 weights: `numpy.array` 

656 Weights (self.sqrtW) 

657 

658 mask : `numpy.array`, optional 

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

660 

661 Notes 

662 ----- 

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

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

665 in electrons or ADU if unitsElectrons=False. 

666 """ 

667 if unitsElectrons: 

668 gain = self.getGain() 

669 else: 

670 gain = 1.0 

671 

672 mu = self.mu*gain 

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

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

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

676 

677 # select data used for the fit: 

678 mask = self.getMaskCov(i, j) 

679 if returnMasked: 

680 weights = weights[mask] 

681 covarianceModel = covarianceModel[mask] 

682 mu = mu[mask] 

683 covariance = covariance[mask] 

684 

685 if divideByMu: 

686 covariance /= mu 

687 covarianceModel /= mu 

688 weights *= mu 

689 

690 return mu, covariance, covarianceModel, weights, mask 

691 

692 def __call__(self, params): 

693 self.setParamValues(params) 

694 chi2 = self.chi2() 

695 

696 return chi2