Coverage for python/lsst/cp/pipe/astierCovPtcFit.py : 15%

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/>.
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
30import lsst.log as lsstLog
32__all__ = ["CovFit"]
35def computeApproximateAcoeffs(covModel, muEl, gain):
36 """Compute the "a" coefficients of the Antilogus+14 (1402.0725) model as in
37 Guyonnet+15 (1501.01577, eq. 16, the slope of cov/var at a given flux mu in electrons).
39 Eq. 16 of 1501.01577 is an approximation to the more complete model in Astier+19 (1905.08677).
41 Parameters
42 ---------
43 covModel : `list`
44 Covariance model from Eq. 20 in Astier+19.
46 muEl : `np.array`
47 Mean signal in electrons
49 gain : `float`
50 Gain in e-/ADU.
52 Returns
53 -------
54 aCoeffsOld: `numpy.array`
55 Slope of cov/var at a given flux mu in electrons.
57 Notes
58 -----
59 Returns the "a" array, computed this way, to be compared to the actual a_array from the full model
60 (fit.geA()).
61 """
62 covModel = np.array(covModel)
63 var = covModel[0, 0, 0] # ADU^2
64 # For a result in electrons^-1, we have to use mu in electrons.
65 return covModel[0, :, :]/(var*muEl)
68def makeCovArray(inputTuple, maxRangeFromTuple=8):
69 """Make covariances array from tuple.
71 Parameters
72 ----------
73 inputTuple: `numpy.recarray`
74 Recarray with rows with at least (mu, cov, var, i, j, npix), where:
75 mu : 0.5*(m1 + m2), where:
76 mu1: mean value of flat1
77 mu2: mean value of flat2
78 cov: covariance value at lag(i, j)
79 var: variance(covariance value at lag(0, 0))
80 i: lag dimension
81 j: lag dimension
82 npix: number of pixels used for covariance calculation.
84 maxRangeFromTuple: `int`
85 Maximum range to select from tuple.
87 Returns
88 -------
89 cov: `numpy.array`
90 Covariance arrays, indexed by mean signal mu.
92 vCov: `numpy.array`
93 Variance arrays, indexed by mean signal mu.
95 muVals: `numpy.array`
96 List of mean signal values.
98 Notes
99 -----
101 The input tuple should contain the following rows:
102 (mu, cov, var, i, j, npix), with one entry per lag, and image pair.
103 Different lags(i.e. different i and j) from the same
104 image pair have the same values of mu1 and mu2. When i==j==0, cov
105 = var.
107 If the input tuple contains several video channels, one should
108 select the data of a given channel *before* entering this
109 routine, as well as apply(e.g.) saturation cuts.
111 The routine returns cov[k_mu, j, i], vcov[(same indices)], and mu[k]
112 where the first index of cov matches the one in mu.
114 This routine implements the loss of variance due to
115 clipping cuts when measuring variances and covariance, but this should happen inside
116 the measurement code, where the cuts are readily available.
118 """
119 if maxRangeFromTuple is not None:
120 cut = (inputTuple['i'] < maxRangeFromTuple) & (inputTuple['j'] < maxRangeFromTuple)
121 cutTuple = inputTuple[cut]
122 else:
123 cutTuple = inputTuple
124 # increasing mu order, so that we can group measurements with the same mu
125 muTemp = cutTuple['mu']
126 ind = np.argsort(muTemp)
128 cutTuple = cutTuple[ind]
129 # should group measurements on the same image pairs(same average)
130 mu = cutTuple['mu']
131 xx = np.hstack(([mu[0]], mu))
132 delta = xx[1:] - xx[:-1]
133 steps, = np.where(delta > 0)
134 ind = np.zeros_like(mu, dtype=int)
135 ind[steps] = 1
136 ind = np.cumsum(ind) # this acts as an image pair index.
137 # now fill the 3-d cov array(and variance)
138 muVals = np.array(np.unique(mu))
139 i = cutTuple['i'].astype(int)
140 j = cutTuple['j'].astype(int)
141 c = 0.5*cutTuple['cov']
142 n = cutTuple['npix']
143 v = 0.5*cutTuple['var']
144 # book and fill
145 cov = np.ndarray((len(muVals), np.max(i)+1, np.max(j)+1))
146 var = np.zeros_like(cov)
147 cov[ind, i, j] = c
148 var[ind, i, j] = v**2/n
149 var[:, 0, 0] *= 2 # var(v) = 2*v**2/N
150 # compensate for loss of variance and covariance due to outlier elimination(sigma clipping)
151 # when computing variances(cut to 4 sigma): 1 per mill for variances and twice as
152 # much for covariances:
153 fact = 1.0 # 1.00107
154 cov *= fact*fact
155 cov[:, 0, 0] /= fact
157 return cov, var, muVals
160def symmetrize(inputArray):
161 """ Copy array over 4 quadrants prior to convolution.
163 Parameters
164 ----------
165 inputarray: `numpy.array`
166 Input array to symmetrize.
168 Returns
169 -------
170 aSym: `numpy.array`
171 Symmetrized array.
172 """
174 targetShape = list(inputArray.shape)
175 r1, r2 = inputArray.shape[-1], inputArray.shape[-2]
176 targetShape[-1] = 2*r1-1
177 targetShape[-2] = 2*r2-1
178 aSym = np.ndarray(tuple(targetShape))
179 aSym[..., r2-1:, r1-1:] = inputArray
180 aSym[..., r2-1:, r1-1::-1] = inputArray
181 aSym[..., r2-1::-1, r1-1::-1] = inputArray
182 aSym[..., r2-1::-1, r1-1:] = inputArray
184 return aSym
187class Pol2d:
188 """A class to calculate 2D polynomials"""
190 def __init__(self, x, y, z, order, w=None):
191 self.orderx = min(order, x.shape[0]-1)
192 self.ordery = min(order, x.shape[1]-1)
193 G = self.monomials(x.ravel(), y.ravel())
194 if w is None:
195 self.coeff, _, rank, _ = np.linalg.lstsq(G, z.ravel())
196 else:
197 self.coeff, _, rank, _ = np.linalg.lstsq((w.ravel()*G.T).T, z.ravel()*w.ravel())
199 def monomials(self, x, y):
200 ncols = (self.orderx+1)*(self.ordery+1)
201 G = np.zeros(x.shape + (ncols,))
202 ij = itertools.product(range(self.orderx+1), range(self.ordery+1))
203 for k, (i, j) in enumerate(ij):
204 G[..., k] = x**i * y**j
206 return G
208 def eval(self, x, y):
209 G = self.monomials(x, y)
211 return np.dot(G, self.coeff)
214class CovFit:
215 """A class to fit the models in Astier+19 to flat covariances.
217 This code implements the model(and the fit thereof) described in
218 Astier+19: https://arxiv.org/pdf/1905.08677.pdf
219 For the time being it uses as input a numpy recarray (tuple with named tags) which
220 contains one row per covariance and per pair: see the routine makeCovArray.
222 Parameters
223 ----------
224 inputTuple: `numpy.recarray`
225 Tuple with at least (mu, cov, var, i, j, npix), where:
226 mu : 0.5*(m1 + m2), where:
227 mu1: mean value of flat1
228 mu2: mean value of flat2
229 cov: covariance value at lag(i, j)
230 var: variance(covariance value at lag(0, 0))
231 i: lag dimension
232 j: lag dimension
233 npix: number of pixels used for covariance calculation.
235 maxRangeFromTuple: `int`
236 Maximum range to select from tuple.
237 """
239 def __init__(self, inputTuple, maxRangeFromTuple=8):
240 self.cov, self.vcov, self.mu = makeCovArray(inputTuple, maxRangeFromTuple)
241 # make it nan safe, replacing nan's with 0 in weights
242 self.sqrtW = np.nan_to_num(1./np.sqrt(self.vcov))
243 self.r = self.cov.shape[1]
244 self.logger = lsstLog.Log.getDefaultLogger()
246 def subtractDistantOffset(self, maxLag=8, startLag=5, polDegree=1):
247 """Subtract a background/offset to the measured covariances.
249 Parameters
250 ---------
251 maxLag: `int`
252 Maximum lag considered
254 startLag: `int`
255 First lag from where to start the offset subtraction.
257 polDegree: `int`
258 Degree of 2D polynomial to fit to covariance to define offse to be subtracted.
259 """
260 assert(startLag < self.r)
261 for k in range(len(self.mu)):
262 # Make a copy because it is going to be altered
263 w = self.sqrtW[k, ...] + 0.
264 sh = w.shape
265 i, j = np.meshgrid(range(sh[0]), range(sh[1]), indexing='ij')
266 # kill the core for the fit
267 w[:startLag, :startLag] = 0
268 poly = Pol2d(i, j, self.cov[k, ...], polDegree+1, w=w)
269 back = poly.eval(i, j)
270 self.cov[k, ...] -= back
271 self.r = maxLag
272 self.cov = self.cov[:, :maxLag, :maxLag]
273 self.vcov = self.vcov[:, :maxLag, :maxLag]
274 self.sqrtW = self.sqrtW[:, :maxLag, :maxLag]
276 return
278 def copy(self):
279 """Make a copy of params"""
280 cop = copy.deepcopy(self)
281 # deepcopy does not work for FitParameters.
282 if hasattr(self, 'params'):
283 cop.params = self.params.copy()
284 return cop
286 def initFit(self):
287 """ Performs a crude parabolic fit of the data in order to start
288 the full fit close to the solution.
289 """
290 # number of parameters for 'a' array.
291 lenA = self.r*self.r
292 # define parameters: c corresponds to a*b in Astier+19 (Eq. 20).
293 self.params = FitParameters([('a', lenA), ('c', lenA), ('noise', lenA), ('gain', 1)])
294 self.params['gain'] = 1.
295 # c=0 in a first go.
296 self.params['c'].fix(val=0.)
297 # plumbing: extract stuff from the parameter structure
298 a = self.params['a'].full.reshape(self.r, self.r)
299 noise = self.params['noise'].full.reshape(self.r, self.r)
300 gain = self.params['gain'].full[0]
302 # iterate the fit to account for higher orders
303 # the chi2 does not necessarily go down, so one could
304 # stop when it increases
305 oldChi2 = 1e30
306 for _ in range(5):
307 model = self.evalCovModel() # this computes the full model.
308 # loop on lags
309 for i in range(self.r):
310 for j in range(self.r):
311 # fit a parabola for a given lag
312 parsFit = np.polyfit(self.mu, self.cov[:, i, j] - model[:, i, j],
313 2, w=self.sqrtW[:, i, j])
314 # model equation(Eq. 20) in Astier+19:
315 a[i, j] += parsFit[0]
316 noise[i, j] += parsFit[2]*gain*gain
317 if(i + j == 0):
318 gain = 1./(1/gain+parsFit[1])
319 self.params['gain'].full[0] = gain
320 chi2 = self.chi2()
321 if chi2 > oldChi2:
322 break
323 oldChi2 = chi2
325 return
327 def getParamValues(self):
328 """Return an array of free parameter values (it is a copy)."""
329 return self.params.free + 0.
331 def setParamValues(self, p):
332 """Set parameter values."""
333 self.params.free = p
334 return
336 def evalCovModel(self, mu=None):
337 """Computes full covariances model (Eq. 20 of Astier+19).
339 Parameters
340 ----------
341 mu: `numpy.array`, optional
342 List of mean signals.
344 Returns
345 -------
346 covModel: `numpy.array`
347 Covariances model.
349 Notes
350 -----
351 By default, computes the covModel for the mu's stored(self.mu).
353 Returns cov[Nmu, self.r, self.r]. The variance for the PTC is cov[:, 0, 0].
354 mu and cov are in ADUs and ADUs squared. To use electrons for both,
355 the gain should be set to 1. This routine implements the model in Astier+19 (1905.08677).
357 The parameters of the full model for C_ij(mu) ("C_ij" and "mu" in ADU^2 and ADU, respectively)
358 in Astier+19 (Eq. 20) are:
360 "a" coefficients (r by r matrix), units: 1/e
361 "b" coefficients (r by r matrix), units: 1/e
362 noise matrix (r by r matrix), units: e^2
363 gain, units: e/ADU
365 "b" appears in Eq. 20 only through the "ab" combination, which is defined in this code as "c=ab".
366 """
367 sa = (self.r, self.r)
368 a = self.params['a'].full.reshape(sa)
369 c = self.params['c'].full.reshape(sa)
370 gain = self.params['gain'].full[0]
371 noise = self.params['noise'].full.reshape(sa)
372 # pad a with zeros and symmetrize
373 aEnlarged = np.zeros((int(sa[0]*1.5)+1, int(sa[1]*1.5)+1))
374 aEnlarged[0:sa[0], 0:sa[1]] = a
375 aSym = symmetrize(aEnlarged)
376 # pad c with zeros and symmetrize
377 cEnlarged = np.zeros((int(sa[0]*1.5)+1, int(sa[1]*1.5)+1))
378 cEnlarged[0:sa[0], 0:sa[1]] = c
379 cSym = symmetrize(cEnlarged)
380 a2 = fftconvolve(aSym, aSym, mode='same')
381 a3 = fftconvolve(a2, aSym, mode='same')
382 ac = fftconvolve(aSym, cSym, mode='same')
383 (xc, yc) = np.unravel_index(np.abs(aSym).argmax(), a2.shape)
384 range = self.r
385 a1 = a[np.newaxis, :, :]
386 a2 = a2[np.newaxis, xc:xc + range, yc:yc + range]
387 a3 = a3[np.newaxis, xc:xc + range, yc:yc + range]
388 ac = ac[np.newaxis, xc:xc + range, yc:yc + range]
389 c1 = c[np.newaxis, ::]
390 if mu is None:
391 mu = self.mu
392 # assumes that mu is 1d
393 bigMu = mu[:, np.newaxis, np.newaxis]*gain
394 # c(=a*b in Astier+19) also has a contribution to the last term, that is absent for now.
395 covModel = (bigMu/(gain*gain)*(a1*bigMu+2./3.*(bigMu*bigMu)*(a2 + c1) +
396 (1./3.*a3 + 5./6.*ac)*(bigMu*bigMu*bigMu)) + noise[np.newaxis, :, :]/gain**2)
397 # add the Poisson term, and the read out noise (variance)
398 covModel[:, 0, 0] += mu/gain
400 return covModel
402 def getA(self):
403 """'a' matrix from Astier+19(e.g., Eq. 20)"""
404 return self.params['a'].full.reshape(self.r, self.r)
406 def getB(self):
407 """'b' matrix from Astier+19(e.g., Eq. 20)"""
408 return self.params['c'].full.reshape(self.r, self.r)/self.getA()
410 def getC(self):
411 """'c'='ab' matrix from Astier+19(e.g., Eq. 20)"""
412 return np.array(self.params['c'].full.reshape(self.r, self.r))
414 def _getCovParams(self, what):
415 """Get covariance matrix of parameters from fit"""
416 indices = self.params[what].indexof()
417 i1 = indices[:, np.newaxis]
418 i2 = indices[np.newaxis, :]
419 if self.covParams is not None:
420 covp = self.covParams[i1, i2]
421 else:
422 covp = None
423 return covp
425 def getACov(self):
426 """Get covariance matrix of "a" coefficients from fit"""
427 if self._getCovParams('a') is not None:
428 cova = self._getCovParams('a').reshape((self.r, self.r, self.r, self.r))
429 else:
430 cova = None
431 return cova
433 def getASig(self):
434 """Square root of diagonal of the parameter covariance of the fitted "a" matrix"""
435 if self._getCovParams('a') is not None:
436 sigA = np.sqrt(self._getCovParams('a').diagonal()).reshape((self.r, self.r))
437 else:
438 sigA = None
439 return sigA
441 def getBCov(self):
442 """Get covariance matrix of "a" coefficients from fit
443 b = c /a
444 """
445 covb = self._getCovParams('c')
446 aval = self.getA().flatten()
447 factor = np.outer(aval, aval)
448 covb /= factor
449 return covb.reshape((self.r, self.r, self.r, self.r))
451 def getCCov(self):
452 """Get covariance matrix of "c" coefficients from fit"""
453 cova = self._getCovParams('c')
454 return cova.reshape((self.r, self.r, self.r, self.r))
456 def getGainErr(self):
457 """Get error on fitted gain parameter"""
458 if self._getCovParams('gain') is not None:
459 gainErr = np.sqrt(self._getCovParams('gain')[0][0])
460 else:
461 gainErr = 0.0
462 return gainErr
464 def getNoiseCov(self):
465 """Get covariances of noise matrix from fit"""
466 covNoise = self._getCovParams('noise')
467 return covNoise.reshape((self.r, self.r, self.r, self.r))
469 def getNoiseSig(self):
470 """Square root of diagonal of the parameter covariance of the fitted "noise" matrix"""
471 if self._getCovParams('noise') is not None:
472 covNoise = self._getCovParams('noise')
473 noise = np.sqrt(covNoise.diagonal()).reshape((self.r, self.r))
474 else:
475 noise = None
476 return noise
478 def getGain(self):
479 """Get gain (e/ADU)"""
480 return self.params['gain'].full[0]
482 def getRon(self):
483 """Get readout noise (e^2)"""
484 return self.params['noise'].full[0]
486 def getRonErr(self):
487 """Get error on readout noise parameter"""
488 ronSqrt = np.sqrt(np.fabs(self.getRon()))
489 if self.getNoiseSig() is not None:
490 noiseSigma = self.getNoiseSig()[0][0]
491 ronErr = 0.5*(noiseSigma/np.fabs(self.getRon()))*ronSqrt
492 else:
493 ronErr = np.nan
494 return ronErr
496 def getNoise(self):
497 """Get noise matrix"""
498 return self.params['noise'].full.reshape(self.r, self.r)
500 def getMaskCov(self, i, j):
501 """Get mask of Cov[i,j]"""
502 weights = self.sqrtW[:, i, j]
503 mask = weights != 0
504 return mask
506 def setAandB(self, a, b):
507 """Set "a" and "b" coeffcients forfull Astier+19 model (Eq. 20). "c=a*b"."""
508 self.params['a'].full = a.flatten()
509 self.params['c'].full = a.flatten()*b.flatten()
510 return
512 def chi2(self):
513 """Calculate weighted chi2 of full-model fit."""
514 return(self.weightedRes()**2).sum()
516 def wres(self, params=None):
517 """To be used in weightedRes"""
518 if params is not None:
519 self.setParamValues(params)
520 covModel = self.evalCovModel()
521 return((covModel-self.cov)*self.sqrtW)
523 def weightedRes(self, params=None):
524 """Weighted residuas.
526 Notes
527 -----
528 To be used via:
529 c = CovFit(nt)
530 c.initFit()
531 coeffs, cov, _, mesg, ierr = leastsq(c.weightedRes, c.getParamValues(), full_output=True)
532 """
533 return self.wres(params).flatten()
535 def fitFullModel(self, pInit=None, nSigma=5.0, maxFitIter=3):
536 """Fit measured covariances to full model in Astier+19 (Eq. 20)
538 Parameters
539 ----------
540 pInit : `list`
541 Initial parameters of the fit.
542 len(pInit) = #entries(a) + #entries(c) + #entries(noise) + 1
543 len(pInit) = r^2 + r^2 + r^2 + 1, where "r" is the maximum lag considered for the
544 covariances calculation, and the extra "1" is the gain.
545 If "b" is 0, then "c" is 0, and len(pInit) will have r^2 fewer entries.
547 nSigma : `float`, optional
548 Sigma cut to get rid of outliers.
550 maxFitIter : `int`, optional
551 Number of iterations for full model fit.
553 Returns
554 -------
555 params : `np.array`
556 Fit parameters (see "Notes" below).
558 Notes
559 -----
560 The parameters of the full model for C_ij(mu) ("C_ij" and "mu" in ADU^2 and ADU, respectively)
561 in Astier+19 (Eq. 20) are:
563 "a" coefficients (r by r matrix), units: 1/e
564 "b" coefficients (r by r matrix), units: 1/e
565 noise matrix (r by r matrix), units: e^2
566 gain, units: e/ADU
568 "b" appears in Eq. 20 only through the "ab" combination, which is defined in this code as "c=ab".
569 """
571 if pInit is None:
572 pInit = self.getParamValues()
573 nOutliers = 1
574 counter = 0
575 while nOutliers != 0:
576 # If fit fails, None values are returned and caught in getOutputPtcDataCovAstier of ptc.py
577 params, paramsCov, _, mesg, ierr = leastsq(self.weightedRes, pInit, full_output=True)
578 wres = self.weightedRes(params)
579 # Do not count the outliers as significant
580 sig = mad(wres[wres != 0], scale='normal')
581 mask = (np.abs(wres) > (nSigma*sig))
582 self.sqrtW.flat[mask] = 0 # flatten makes a copy
583 nOutliers = mask.sum()
584 counter += 1
585 if counter == maxFitIter:
586 self.logger.warn(f"Max number of iterations,{maxFitIter}, reached in fitFullModel")
587 break
589 self.covParams = paramsCov
590 return params
592 def ndof(self):
593 """Number of degrees of freedom
595 Returns
596 -------
597 mask.sum() - len(self.params.free): `int`
598 Number of usable pixels - number of parameters of fit.
599 """
600 mask = self.sqrtW != 0
602 return mask.sum() - len(self.params.free)
604 def getFitData(self, i, j, divideByMu=False, unitsElectrons=False, returnMasked=False):
605 """Get measured signal and covariance, cov model, weigths, and mask at covariance lag (i, j).
607 Parameters
608 ---------
609 i: `int`
610 Lag for covariance matrix.
612 j: `int`
613 Lag for covariance matrix.
615 divideByMu: `bool`, optional
616 Divide covariance, model, and weights by signal mu?
618 unitsElectrons : `bool`, optional
619 mu, covariance, and model are in ADU (or powers of ADU) If tthis parameter is true, these are
620 multiplied by the adequte factors of the gain to return quantities in electrons
621 (or powers of electrons).
623 returnMasked : `bool`, optional
624 Use mask (based on weights) in returned arrays (mu, covariance, and model)?
626 Returns
627 -------
628 mu: `numpy.array`
629 list of signal values (mu).
631 covariance: `numpy.array`
632 Covariance arrays, indexed by mean signal mu (self.cov[:, i, j]).
634 covarianceModel: `numpy.array`
635 Covariance model (model).
637 weights: `numpy.array`
638 Weights (self.sqrtW)
640 mask : `numpy.array`, optional
641 Boolean mask of the covariance at (i,j).
643 Notes
644 -----
645 Using a CovFit object, selects from (i, j) and returns
646 mu*gain, self.cov[:, i, j]*gain**2 model*gain**2, and self.sqrtW/gain**2
647 in electrons or ADU if unitsElectrons=False.
648 """
649 if unitsElectrons:
650 gain = self.getGain()
651 else:
652 gain = 1.0
654 mu = self.mu*gain
655 covariance = self.cov[:, i, j]*(gain**2)
656 covarianceModel = self.evalCovModel()[:, i, j]*(gain**2)
657 weights = self.sqrtW[:, i, j]/(gain**2)
659 # select data used for the fit:
660 mask = self.getMaskCov(i, j)
661 if returnMasked:
662 weights = weights[mask]
663 covarianceModel = covarianceModel[mask]
664 mu = mu[mask]
665 covariance = covariance[mask]
667 if divideByMu:
668 covariance /= mu
669 covarianceModel /= mu
670 weights *= mu
672 return mu, covariance, covarianceModel, weights, mask
674 def __call__(self, params):
675 self.setParamValues(params)
676 chi2 = self.chi2()
678 return chi2