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_absolute_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(fit, muEl): 

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 fit: `lsst.cp.pipe.astierCovPtcFit.CovFit` 

43 CovFit object 

44 

45 Returns 

46 ------- 

47 aCoeffsOld: `numpy.array` 

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

49 

50 Notes 

51 ----- 

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

53 (fit.geA()). 

54 """ 

55 

56 gain = fit.getGain() 

57 muAdu = np.array([muEl/gain]) 

58 model = fit.evalCovModel(muAdu) 

59 var = model[0, 0, 0] # ADU^2 

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

61 

62 return model[0, :, :]/(var*muEl) 

63 

64 

65def makeCovArray(inputTuple, maxRangeFromTuple=8): 

66 """Make covariances array from tuple. 

67 

68 Parameters 

69 ---------- 

70 inputTuple: `numpy.recarray` 

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

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

73 mu1: mean value of flat1 

74 mu2: mean value of flat2 

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

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

77 i: lag dimension 

78 j: lag dimension 

79 npix: number of pixels used for covariance calculation. 

80 

81 maxRangeFromTuple: `int` 

82 Maximum range to select from tuple. 

83 

84 Returns 

85 ------- 

86 cov: `numpy.array` 

87 Covariance arrays, indexed by mean signal mu. 

88 

89 vCov: `numpy.array` 

90 Variance arrays, indexed by mean signal mu. 

91 

92 muVals: `numpy.array` 

93 List of mean signal values. 

94 

95 Notes 

96 ----- 

97 

98 The input tuple should contain the following rows: 

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

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

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

102 = var. 

103 

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

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

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

107 

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

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

110 

111 This routine implements the loss of variance due to 

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

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

114 

115 """ 

116 if maxRangeFromTuple is not None: 

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

118 cutTuple = inputTuple[cut] 

119 else: 

120 cutTuple = inputTuple 

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

122 muTemp = cutTuple['mu'] 

123 ind = np.argsort(muTemp) 

124 

125 cutTuple = cutTuple[ind] 

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

127 mu = cutTuple['mu'] 

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

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

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

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

132 ind[steps] = 1 

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

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

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

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

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

138 c = 0.5*cutTuple['cov'] 

139 n = cutTuple['npix'] 

140 v = 0.5*cutTuple['var'] 

141 # book and fill 

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

143 var = np.zeros_like(cov) 

144 cov[ind, i, j] = c 

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

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

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

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

149 # much for covariances: 

150 fact = 1.0 # 1.00107 

151 cov *= fact*fact 

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

153 

154 return cov, var, muVals 

155 

156 

157def symmetrize(inputArray): 

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

159 

160 Parameters 

161 ---------- 

162 inputarray: `numpy.array` 

163 Input array to symmetrize. 

164 

165 Returns 

166 ------- 

167 aSym: `numpy.array` 

168 Symmetrized array. 

169 """ 

170 

171 targetShape = list(inputArray.shape) 

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

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

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

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

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

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

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

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

180 

181 return aSym 

182 

183 

184class Pol2d: 

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

186 

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

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

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

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

191 if w is None: 

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

193 else: 

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

195 

196 def monomials(self, x, y): 

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

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

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

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

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

202 

203 return G 

204 

205 def eval(self, x, y): 

206 G = self.monomials(x, y) 

207 

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

209 

210 

211class CovFit: 

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

213 

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

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

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

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

218 

219 Parameters 

220 ---------- 

221 inputTuple: `numpy.recarray` 

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

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

224 mu1: mean value of flat1 

225 mu2: mean value of flat2 

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

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

228 i: lag dimension 

229 j: lag dimension 

230 npix: number of pixels used for covariance calculation. 

231 

232 maxRangeFromTuple: `int` 

233 Maximum range to select from tuple. 

234 """ 

235 

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

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

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

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

240 

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

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

243 

244 Parameters 

245 --------- 

246 maxLag: `int` 

247 Maximum lag considered 

248 

249 startLag: `int` 

250 First lag from where to start the offset subtraction. 

251 

252 polDegree: `int` 

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

254 """ 

255 assert(startLag < self.r) 

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

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

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

259 sh = w.shape 

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

261 # kill the core for the fit 

262 w[:startLag, :startLag] = 0 

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

264 back = poly.eval(i, j) 

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

266 self.r = maxLag 

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

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

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

270 

271 return 

272 

273 def setMaxMu(self, maxMu): 

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

275 # mus are sorted at construction 

276 index = self.mu < maxMu 

277 k = index.sum() 

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

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

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

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

282 

283 return 

284 

285 def setMaxMuElectrons(self, maxMuEl): 

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

287 g = self.getGain() 

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

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

290 

291 return 

292 

293 def copy(self): 

294 """Make a copy of params""" 

295 cop = copy.deepcopy(self) 

296 # deepcopy does not work for FitParameters. 

297 if hasattr(self, 'params'): 

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

299 return cop 

300 

301 def initFit(self): 

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

303 the full fit close to the solution. 

304 """ 

305 # number of parameters for 'a' array. 

306 lenA = self.r*self.r 

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

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

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

310 # c=0 in a first go. 

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

312 # plumbing: extract stuff from the parameter structure 

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

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

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

316 

317 # iterate the fit to account for higher orders 

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

319 # stop when it increases 

320 oldChi2 = 1e30 

321 for _ in range(5): 

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

323 # loop on lags 

324 for i in range(self.r): 

325 for j in range(self.r): 

326 # fit a parabola for a given lag 

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

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

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

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

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

332 if(i + j == 0): 

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

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

335 chi2 = self.chi2() 

336 if chi2 > oldChi2: 

337 break 

338 oldChi2 = chi2 

339 

340 return 

341 

342 def getParamValues(self): 

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

344 return self.params.free + 0. 

345 

346 def setParamValues(self, p): 

347 """Set parameter values.""" 

348 self.params.free = p 

349 return 

350 

351 def evalCovModel(self, mu=None): 

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

353 

354 Parameters 

355 ---------- 

356 mu: `numpy.array`, optional 

357 List of mean signals. 

358 

359 Returns 

360 ------- 

361 covModel: `numpy.array` 

362 Covariances model. 

363 

364 Notes 

365 ----- 

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

367 

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

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

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

371 

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

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

374 

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

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

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

378 gain, units: e/ADU 

379 

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

381 """ 

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

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

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

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

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

387 # pad a with zeros and symmetrize 

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

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

390 aSym = symmetrize(aEnlarged) 

391 # pad c with zeros and symmetrize 

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

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

394 cSym = symmetrize(cEnlarged) 

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

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

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

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

399 range = self.r 

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

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

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

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

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

405 if mu is None: 

406 mu = self.mu 

407 # assumes that mu is 1d 

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

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

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

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

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

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

414 

415 return covModel 

416 

417 def getA(self): 

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

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

420 

421 def getB(self): 

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

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

424 

425 def getC(self): 

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

427 return self.params['c'].full.reshape(self.r, self.r) 

428 

429 def _getCovParams(self, what): 

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

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

432 i1 = indices[:, np.newaxis] 

433 i2 = indices[np.newaxis, :] 

434 if self.covParams is not None: 

435 covp = self.covParams[i1, i2] 

436 else: 

437 covp = None 

438 return covp 

439 

440 def getACov(self): 

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

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

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

444 else: 

445 cova = None 

446 return cova 

447 

448 def getASig(self): 

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

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

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

452 else: 

453 sigA = None 

454 return sigA 

455 

456 def getBCov(self): 

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

458 b = c /a 

459 """ 

460 covb = self._getCovParams('c') 

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

462 factor = np.outer(aval, aval) 

463 covb /= factor 

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

465 

466 def getCCov(self): 

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

468 cova = self._getCovParams('c') 

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

470 

471 def getGainErr(self): 

472 """Get error on fitted gain parameter""" 

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

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

475 else: 

476 gainErr = 0.0 

477 return gainErr 

478 

479 def getNoiseCov(self): 

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

481 covNoise = self._getCovParams('noise') 

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

483 

484 def getNoiseSig(self): 

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

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

487 covNoise = self._getCovParams('noise') 

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

489 else: 

490 noise = None 

491 return noise 

492 

493 def getGain(self): 

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

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

496 

497 def getRon(self): 

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

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

500 

501 def getRonErr(self): 

502 """Get error on readout noise parameter""" 

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

504 if self.getNoiseSig() is not None: 

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

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

507 else: 

508 ronErr = np.nan 

509 return ronErr 

510 

511 def getNoise(self): 

512 """Get noise matrix""" 

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

514 

515 def getMaskCov(self, i, j): 

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

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

518 mask = weights != 0 

519 return mask 

520 

521 def setAandB(self, a, b): 

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

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

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

525 return 

526 

527 def chi2(self): 

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

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

530 

531 def wres(self, params=None): 

532 """To be used in weightedRes""" 

533 if params is not None: 

534 self.setParamValues(params) 

535 covModel = self.evalCovModel() 

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

537 

538 def weightedRes(self, params=None): 

539 """Weighted residuas. 

540 

541 Notes 

542 ----- 

543 To be used via: 

544 c = CovFit(nt) 

545 c.initFit() 

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

547 """ 

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

549 

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

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

552 

553 Parameters 

554 ---------- 

555 pInit : `list` 

556 Initial parameters of the fit. 

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

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

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

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

561 

562 nSigma : `float`, optional 

563 Sigma cut to get rid of outliers. 

564 

565 maxFitIter : `int`, optional 

566 Number of iterations for full model fit. 

567 

568 Returns 

569 ------- 

570 params : `np.array` 

571 Fit parameters (see "Notes" below). 

572 

573 Notes 

574 ----- 

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

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

577 

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

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

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

581 gain, units: e/ADU 

582 

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

584 """ 

585 

586 if pInit is None: 

587 pInit = self.getParamValues() 

588 nOutliers = 1 

589 counter = 0 

590 while nOutliers != 0: 

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

592 wres = self.weightedRes(params) 

593 # Do not count the outliers as significant 

594 sig = mad(wres[wres != 0]) 

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

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

597 nOutliers = mask.sum() 

598 counter += 1 

599 if counter == maxFitIter: 

600 break 

601 

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

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

604 self.covParams = paramsCov 

605 

606 return params 

607 

608 def ndof(self): 

609 """Number of degrees of freedom 

610 

611 Returns 

612 ------- 

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

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

615 """ 

616 mask = self.sqrtW != 0 

617 

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

619 

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

621 """Get measured signal and covariance, cov model, weigths, and mask. 

622 

623 Parameters 

624 --------- 

625 i: `int` 

626 Lag for covariance 

627 

628 j: `int` 

629 Lag for covariance 

630 

631 divideByMu: `bool`, optional 

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

633 

634 unitsElectrons : `bool`, optional 

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

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

637 (or powers of electrons). 

638 

639 returneMasked : `bool`, optional 

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

641 

642 Returns 

643 ------- 

644 mu: `numpy.array` 

645 list of signal values (mu). 

646 

647 covariance: `numpy.array` 

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

649 

650 model: `numpy.array` 

651 Covariance model (model). 

652 

653 weights: `numpy.array` 

654 Weights (self.sqrtW) 

655 

656 mask : `numpy.array`, optional 

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

658 

659 Notes 

660 ----- 

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

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

663 in electrons or ADU if unitsElectrons=False. 

664 """ 

665 if unitsElectrons: 

666 gain = self.getGain() 

667 else: 

668 gain = 1.0 

669 

670 mu = self.mu*gain 

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

672 model = self.evalCovModel()[:, i, j]*(gain**2) 

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

674 

675 # select data used for the fit: 

676 mask = self.getMaskCov(i, j) 

677 if returnMasked: 

678 weights = weights[mask] 

679 model = model[mask] 

680 mu = mu[mask] 

681 covariance = covariance[mask] 

682 

683 if divideByMu: 

684 covariance /= mu 

685 model /= mu 

686 weights *= mu 

687 

688 return mu, covariance, model, weights, mask 

689 

690 def __call__(self, params): 

691 self.setParamValues(params) 

692 chi2 = self.chi2() 

693 

694 return chi2