Coverage for python/lsst/cp/pipe/ptc/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.signal import fftconvolve
26from scipy.optimize import leastsq
27from .astierCovFitParameters import FitParameters
29import lsst.log as lsstLog
31__all__ = ["CovFit"]
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).
38 Eq. 16 of 1501.01577 is an approximation to the more complete model in Astier+19 (1905.08677).
40 Parameters
41 ---------
42 covModel : `list`
43 Covariance model from Eq. 20 in Astier+19.
45 muEl : `np.array`
46 Mean signal in electrons
48 gain : `float`
49 Gain in e-/ADU.
51 Returns
52 -------
53 aCoeffsOld: `numpy.array`
54 Slope of cov/var at a given flux mu in electrons.
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)
67def makeCovArray(inputTuple, maxRangeFromTuple=8):
68 """Make covariances array from tuple.
70 Parameters
71 ----------
72 inputTuple: `numpy.ndarray`
73 Structured array with rows with at least
74 (mu, afwVar, cov, var, i, j, npix), where:
76 mu : 0.5*(m1 + m2), where:
77 mu1: mean value of flat1
78 mu2: mean value of flat2
79 afwVar: variance of difference flat, calculated with afw
80 cov: covariance value at lag(i, j)
81 var: variance(covariance value at lag(0, 0))
82 i: lag dimension
83 j: lag dimension
84 npix: number of pixels used for covariance calculation.
86 maxRangeFromTuple: `int`
87 Maximum range to select from tuple.
89 Returns
90 -------
91 cov: `numpy.array`
92 Covariance arrays, indexed by mean signal mu.
94 vCov: `numpy.array`
95 Variance arrays, indexed by mean signal mu.
97 muVals: `numpy.array`
98 List of mean signal values.
100 Notes
101 -----
103 The input tuple should contain the following rows:
104 (mu, cov, var, i, j, npix), with one entry per lag, and image pair.
105 Different lags(i.e. different i and j) from the same
106 image pair have the same values of mu1 and mu2. When i==j==0, cov
107 = var.
109 If the input tuple contains several video channels, one should
110 select the data of a given channel *before* entering this
111 routine, as well as apply(e.g.) saturation cuts.
113 The routine returns cov[k_mu, j, i], vcov[(same indices)], and mu[k]
114 where the first index of cov matches the one in mu.
116 This routine implements the loss of variance due to
117 clipping cuts when measuring variances and covariance, but this should happen inside
118 the measurement code, where the cuts are readily available.
120 """
121 if maxRangeFromTuple is not None:
122 cut = (inputTuple['i'] < maxRangeFromTuple) & (inputTuple['j'] < maxRangeFromTuple)
123 cutTuple = inputTuple[cut]
124 else:
125 cutTuple = inputTuple
126 # increasing mu order, so that we can group measurements with the same mu
127 muTemp = cutTuple['mu']
128 ind = np.argsort(muTemp)
130 cutTuple = cutTuple[ind]
131 # should group measurements on the same image pairs(same average)
132 mu = cutTuple['mu']
133 xx = np.hstack(([mu[0]], mu))
134 delta = xx[1:] - xx[:-1]
135 steps, = np.where(delta > 0)
136 ind = np.zeros_like(mu, dtype=int)
137 ind[steps] = 1
138 ind = np.cumsum(ind) # this acts as an image pair index.
139 # now fill the 3-d cov array(and variance)
140 muVals = np.array(np.unique(mu))
141 i = cutTuple['i'].astype(int)
142 j = cutTuple['j'].astype(int)
143 c = 0.5*cutTuple['cov']
144 n = cutTuple['npix']
145 v = 0.5*cutTuple['var']
146 # book and fill
147 cov = np.ndarray((len(muVals), np.max(i)+1, np.max(j)+1))
148 var = np.zeros_like(cov)
149 cov[ind, i, j] = c
150 var[ind, i, j] = v**2/n
151 var[:, 0, 0] *= 2 # var(v) = 2*v**2/N
152 # compensate for loss of variance and covariance due to outlier elimination(sigma clipping)
153 # when computing variances(cut to 4 sigma): 1 per mill for variances and twice as
154 # much for covariances:
155 fact = 1.0 # 1.00107
156 cov *= fact*fact
157 cov[:, 0, 0] /= fact
159 return cov, var, muVals
162def symmetrize(inputArray):
163 """ Copy array over 4 quadrants prior to convolution.
165 Parameters
166 ----------
167 inputarray: `numpy.array`
168 Input array to symmetrize.
170 Returns
171 -------
172 aSym: `numpy.array`
173 Symmetrized array.
174 """
176 targetShape = list(inputArray.shape)
177 r1, r2 = inputArray.shape[-1], inputArray.shape[-2]
178 targetShape[-1] = 2*r1-1
179 targetShape[-2] = 2*r2-1
180 aSym = np.ndarray(tuple(targetShape))
181 aSym[..., r2-1:, r1-1:] = inputArray
182 aSym[..., r2-1:, r1-1::-1] = inputArray
183 aSym[..., r2-1::-1, r1-1::-1] = inputArray
184 aSym[..., r2-1::-1, r1-1:] = inputArray
186 return aSym
189class Pol2d:
190 """A class to calculate 2D polynomials"""
192 def __init__(self, x, y, z, order, w=None):
193 self.orderx = min(order, x.shape[0]-1)
194 self.ordery = min(order, x.shape[1]-1)
195 G = self.monomials(x.ravel(), y.ravel())
196 if w is None:
197 self.coeff, _, rank, _ = np.linalg.lstsq(G, z.ravel())
198 else:
199 self.coeff, _, rank, _ = np.linalg.lstsq((w.ravel()*G.T).T, z.ravel()*w.ravel())
201 def monomials(self, x, y):
202 ncols = (self.orderx+1)*(self.ordery+1)
203 G = np.zeros(x.shape + (ncols,))
204 ij = itertools.product(range(self.orderx+1), range(self.ordery+1))
205 for k, (i, j) in enumerate(ij):
206 G[..., k] = x**i * y**j
208 return G
210 def eval(self, x, y):
211 G = self.monomials(x, y)
213 return np.dot(G, self.coeff)
216class CovFit:
217 """A class to fit the models in Astier+19 to flat covariances.
219 This code implements the model(and the fit thereof) described in
220 Astier+19: https://arxiv.org/pdf/1905.08677.pdf
222 Parameters
223 ----------
224 meanSignals : `list`[`float`]
225 List with means of the difference image of two flats,
226 for a particular amplifier in the detector.
228 covariances : `list`[`numpy.array`]
229 List with 2D covariance arrays at a given mean signal.
231 covsSqrtWeights : `list`[`numpy.array`]
232 List with 2D arrays with weights from `vcov as defined in
233 `makeCovArray`: weight = 1/sqrt(vcov).
235 maxRangeFromTuple: `int`, optional
236 Maximum range to select from tuple.
238 meanSignalMask: `list`[`bool`], optional
239 Mask of mean signal 1D array. Use all entries if empty.
240 """
242 def __init__(self, meanSignals, covariances, covsSqrtWeights, maxRangeFromTuple=8, meanSignalsMask=[]):
243 assert (len(meanSignals) == len(covariances))
244 assert (len(covariances) == len(covsSqrtWeights))
245 if len(meanSignalsMask) == 0:
246 meanSignalsMask = np.repeat(True, len(meanSignals))
247 self.mu = meanSignals[meanSignalsMask]
248 self.cov = np.nan_to_num(covariances)[meanSignalsMask]
249 # make it nan safe, replacing nan's with 0 in weights
250 self.sqrtW = np.nan_to_num(covsSqrtWeights)[meanSignalsMask]
251 self.r = maxRangeFromTuple
252 self.logger = lsstLog.Log.getDefaultLogger()
254 def subtractDistantOffset(self, maxLag=8, startLag=5, polDegree=1):
255 """Subtract a background/offset to the measured covariances.
257 Parameters
258 ---------
259 maxLag: `int`
260 Maximum lag considered
262 startLag: `int`
263 First lag from where to start the offset subtraction.
265 polDegree: `int`
266 Degree of 2D polynomial to fit to covariance to define offse to be subtracted.
267 """
268 assert(startLag < self.r)
269 for k in range(len(self.mu)):
270 # Make a copy because it is going to be altered
271 w = self.sqrtW[k, ...] + 0.
272 sh = w.shape
273 i, j = np.meshgrid(range(sh[0]), range(sh[1]), indexing='ij')
274 # kill the core for the fit
275 w[:startLag, :startLag] = 0
276 poly = Pol2d(i, j, self.cov[k, ...], polDegree+1, w=w)
277 back = poly.eval(i, j)
278 self.cov[k, ...] -= back
279 self.r = maxLag
280 self.cov = self.cov[:, :maxLag, :maxLag]
281 self.sqrtW = self.sqrtW[:, :maxLag, :maxLag]
283 return
285 def copy(self):
286 """Make a copy of params"""
287 cop = copy.deepcopy(self)
288 # deepcopy does not work for FitParameters.
289 if hasattr(self, 'params'):
290 cop.params = self.params.copy()
291 return cop
293 def initFit(self):
294 """ Performs a crude parabolic fit of the data in order to start
295 the full fit close to the solution.
296 """
297 # number of parameters for 'a' array.
298 lenA = self.r*self.r
299 # define parameters: c corresponds to a*b in Astier+19 (Eq. 20).
300 self.params = FitParameters([('a', lenA), ('c', lenA), ('noise', lenA), ('gain', 1)])
301 self.params['gain'] = 1.
302 # c=0 in a first go.
303 self.params['c'].fix(val=0.)
304 # plumbing: extract stuff from the parameter structure
305 a = self.params['a'].full.reshape(self.r, self.r)
306 noise = self.params['noise'].full.reshape(self.r, self.r)
307 gain = self.params['gain'].full[0]
309 # iterate the fit to account for higher orders
310 # the chi2 does not necessarily go down, so one could
311 # stop when it increases
312 oldChi2 = 1e30
313 for _ in range(5):
314 model = np.nan_to_num(self.evalCovModel()) # this computes the full model.
315 # loop on lags
316 for i in range(self.r):
317 for j in range(self.r):
318 # fit a parabola for a given lag
319 parsFit = np.polyfit(self.mu, self.cov[:, i, j] - model[:, i, j],
320 2, w=self.sqrtW[:, i, j])
321 # model equation(Eq. 20) in Astier+19:
322 a[i, j] += parsFit[0]
323 noise[i, j] += parsFit[2]*gain*gain
324 if(i + j == 0):
325 gain = 1./(1/gain+parsFit[1])
326 self.params['gain'].full[0] = gain
327 chi2 = self.chi2()
328 if chi2 > oldChi2:
329 break
330 oldChi2 = chi2
332 return
334 def getParamValues(self):
335 """Return an array of free parameter values (it is a copy)."""
336 return self.params.free + 0.
338 def setParamValues(self, p):
339 """Set parameter values."""
340 self.params.free = p
341 return
343 def evalCovModel(self, mu=None):
344 """Computes full covariances model (Eq. 20 of Astier+19).
346 Parameters
347 ----------
348 mu: `numpy.array`, optional
349 List of mean signals.
351 Returns
352 -------
353 covModel: `numpy.array`
354 Covariances model.
356 Notes
357 -----
358 By default, computes the covModel for the mu's stored(self.mu).
360 Returns cov[Nmu, self.r, self.r]. The variance for the PTC is cov[:, 0, 0].
361 mu and cov are in ADUs and ADUs squared. To use electrons for both,
362 the gain should be set to 1. This routine implements the model in Astier+19 (1905.08677).
364 The parameters of the full model for C_ij(mu) ("C_ij" and "mu" in ADU^2 and ADU, respectively)
365 in Astier+19 (Eq. 20) are:
367 "a" coefficients (r by r matrix), units: 1/e
368 "b" coefficients (r by r matrix), units: 1/e
369 noise matrix (r by r matrix), units: e^2
370 gain, units: e/ADU
372 "b" appears in Eq. 20 only through the "ab" combination, which is defined in this code as "c=ab".
373 """
374 sa = (self.r, self.r)
375 a = self.params['a'].full.reshape(sa)
376 c = self.params['c'].full.reshape(sa)
377 gain = self.params['gain'].full[0]
378 noise = self.params['noise'].full.reshape(sa)
379 # pad a with zeros and symmetrize
380 aEnlarged = np.zeros((int(sa[0]*1.5)+1, int(sa[1]*1.5)+1))
381 aEnlarged[0:sa[0], 0:sa[1]] = a
382 aSym = symmetrize(aEnlarged)
383 # pad c with zeros and symmetrize
384 cEnlarged = np.zeros((int(sa[0]*1.5)+1, int(sa[1]*1.5)+1))
385 cEnlarged[0:sa[0], 0:sa[1]] = c
386 cSym = symmetrize(cEnlarged)
387 a2 = fftconvolve(aSym, aSym, mode='same')
388 a3 = fftconvolve(a2, aSym, mode='same')
389 ac = fftconvolve(aSym, cSym, mode='same')
390 (xc, yc) = np.unravel_index(np.abs(aSym).argmax(), a2.shape)
391 range = self.r
392 a1 = a[np.newaxis, :, :]
393 a2 = a2[np.newaxis, xc:xc + range, yc:yc + range]
394 a3 = a3[np.newaxis, xc:xc + range, yc:yc + range]
395 ac = ac[np.newaxis, xc:xc + range, yc:yc + range]
396 c1 = c[np.newaxis, ::]
397 if mu is None:
398 mu = self.mu
399 # assumes that mu is 1d
400 bigMu = mu[:, np.newaxis, np.newaxis]*gain
401 # c(=a*b in Astier+19) also has a contribution to the last term, that is absent for now.
402 covModel = (bigMu/(gain*gain)*(a1*bigMu+2./3.*(bigMu*bigMu)*(a2 + c1)
403 + (1./3.*a3 + 5./6.*ac)*(bigMu*bigMu*bigMu)) + noise[np.newaxis, :, :]/gain**2)
404 # add the Poisson term, and the read out noise (variance)
405 covModel[:, 0, 0] += mu/gain
407 return covModel
409 def getA(self):
410 """'a' matrix from Astier+19(e.g., Eq. 20)"""
411 return self.params['a'].full.reshape(self.r, self.r)
413 def getB(self):
414 """'b' matrix from Astier+19(e.g., Eq. 20)"""
415 return self.params['c'].full.reshape(self.r, self.r)/self.getA()
417 def getC(self):
418 """'c'='ab' matrix from Astier+19(e.g., Eq. 20)"""
419 return np.array(self.params['c'].full.reshape(self.r, self.r))
421 def _getCovParams(self, what):
422 """Get covariance matrix of parameters from fit"""
423 indices = self.params[what].indexof()
424 i1 = indices[:, np.newaxis]
425 i2 = indices[np.newaxis, :]
426 if self.covParams is not None:
427 covp = self.covParams[i1, i2]
428 else:
429 covp = None
430 return covp
432 def getACov(self):
433 """Get covariance matrix of "a" coefficients from fit"""
434 if self._getCovParams('a') is not None:
435 cova = self._getCovParams('a').reshape((self.r, self.r, self.r, self.r))
436 else:
437 cova = None
438 return cova
440 def getASig(self):
441 """Square root of diagonal of the parameter covariance of the fitted "a" matrix"""
442 if self._getCovParams('a') is not None:
443 sigA = np.sqrt(self._getCovParams('a').diagonal()).reshape((self.r, self.r))
444 else:
445 sigA = None
446 return sigA
448 def getBCov(self):
449 """Get covariance matrix of "a" coefficients from fit
450 b = c /a
451 """
452 covb = self._getCovParams('c')
453 aval = self.getA().flatten()
454 factor = np.outer(aval, aval)
455 covb /= factor
456 return covb.reshape((self.r, self.r, self.r, self.r))
458 def getCCov(self):
459 """Get covariance matrix of "c" coefficients from fit"""
460 cova = self._getCovParams('c')
461 return cova.reshape((self.r, self.r, self.r, self.r))
463 def getGainErr(self):
464 """Get error on fitted gain parameter"""
465 if self._getCovParams('gain') is not None:
466 gainErr = np.sqrt(self._getCovParams('gain')[0][0])
467 else:
468 gainErr = 0.0
469 return gainErr
471 def getNoiseCov(self):
472 """Get covariances of noise matrix from fit"""
473 covNoise = self._getCovParams('noise')
474 return covNoise.reshape((self.r, self.r, self.r, self.r))
476 def getNoiseSig(self):
477 """Square root of diagonal of the parameter covariance of the fitted "noise" matrix"""
478 if self._getCovParams('noise') is not None:
479 covNoise = self._getCovParams('noise')
480 noise = np.sqrt(covNoise.diagonal()).reshape((self.r, self.r))
481 else:
482 noise = None
483 return noise
485 def getGain(self):
486 """Get gain (e/ADU)"""
487 return self.params['gain'].full[0]
489 def getRon(self):
490 """Get readout noise (e^2)"""
491 return self.params['noise'].full[0]
493 def getRonErr(self):
494 """Get error on readout noise parameter"""
495 ronSqrt = np.sqrt(np.fabs(self.getRon()))
496 if self.getNoiseSig() is not None:
497 noiseSigma = self.getNoiseSig()[0][0]
498 ronErr = 0.5*(noiseSigma/np.fabs(self.getRon()))*ronSqrt
499 else:
500 ronErr = np.nan
501 return ronErr
503 def getNoise(self):
504 """Get noise matrix"""
505 return self.params['noise'].full.reshape(self.r, self.r)
507 def getMaskCov(self, i, j):
508 """Get mask of Cov[i,j]"""
509 weights = self.sqrtW[:, i, j]
510 mask = weights != 0
511 return mask
513 def setAandB(self, a, b):
514 """Set "a" and "b" coeffcients forfull Astier+19 model (Eq. 20). "c=a*b"."""
515 self.params['a'].full = a.flatten()
516 self.params['c'].full = a.flatten()*b.flatten()
517 return
519 def chi2(self):
520 """Calculate weighted chi2 of full-model fit."""
521 return(self.weightedRes()**2).sum()
523 def wres(self, params=None):
524 """To be used in weightedRes"""
525 if params is not None:
526 self.setParamValues(params)
527 covModel = np.nan_to_num(self.evalCovModel())
528 weightedRes = (covModel-self.cov)*self.sqrtW
529 maskedWeightedRes = weightedRes
531 return maskedWeightedRes
533 def weightedRes(self, params=None):
534 """Weighted residuals.
536 Notes
537 -----
538 To be used via:
539 c = CovFit(meanSignals, covariances, covsSqrtWeights)
540 c.initFit()
541 coeffs, cov, _, mesg, ierr = leastsq(c.weightedRes, c.getParamValues(), full_output=True)
542 """
543 return self.wres(params).flatten()
545 def fitFullModel(self, pInit=None):
546 """Fit measured covariances to full model in Astier+19 (Eq. 20)
548 Parameters
549 ----------
550 pInit : `list`
551 Initial parameters of the fit.
552 len(pInit) = #entries(a) + #entries(c) + #entries(noise) + 1
553 len(pInit) = r^2 + r^2 + r^2 + 1, where "r" is the maximum lag considered for the
554 covariances calculation, and the extra "1" is the gain.
555 If "b" is 0, then "c" is 0, and len(pInit) will have r^2 fewer entries.
557 Returns
558 -------
559 params : `np.array`
560 Fit parameters (see "Notes" below).
562 Notes
563 -----
564 The parameters of the full model for C_ij(mu) ("C_ij" and "mu" in ADU^2 and ADU, respectively)
565 in Astier+19 (Eq. 20) are:
567 "a" coefficients (r by r matrix), units: 1/e
568 "b" coefficients (r by r matrix), units: 1/e
569 noise matrix (r by r matrix), units: e^2
570 gain, units: e/ADU
572 "b" appears in Eq. 20 only through the "ab" combination, which is defined in this code as "c=ab".
573 """
575 if pInit is None:
576 pInit = self.getParamValues()
577 params, paramsCov, _, mesg, ierr = leastsq(self.weightedRes, pInit, full_output=True)
578 self.covParams = paramsCov
580 return params
582 def ndof(self):
583 """Number of degrees of freedom
585 Returns
586 -------
587 mask.sum() - len(self.params.free): `int`
588 Number of usable pixels - number of parameters of fit.
589 """
590 mask = self.sqrtW != 0
592 return mask.sum() - len(self.params.free)
594 def getFitData(self, i, j, divideByMu=False, unitsElectrons=False, returnMasked=False):
595 """Get measured signal and covariance, cov model, weigths, and mask at covariance lag (i, j).
597 Parameters
598 ---------
599 i: `int`
600 Lag for covariance matrix.
602 j: `int`
603 Lag for covariance matrix.
605 divideByMu: `bool`, optional
606 Divide covariance, model, and weights by signal mu?
608 unitsElectrons : `bool`, optional
609 mu, covariance, and model are in ADU (or powers of ADU) If tthis parameter is true, these are
610 multiplied by the adequte factors of the gain to return quantities in electrons
611 (or powers of electrons).
613 returnMasked : `bool`, optional
614 Use mask (based on weights) in returned arrays (mu, covariance, and model)?
616 Returns
617 -------
618 mu: `numpy.array`
619 list of signal values (mu).
621 covariance: `numpy.array`
622 Covariance arrays, indexed by mean signal mu (self.cov[:, i, j]).
624 covarianceModel: `numpy.array`
625 Covariance model (model).
627 weights: `numpy.array`
628 Weights (self.sqrtW)
630 mask : `numpy.array`, optional
631 Boolean mask of the covariance at (i,j).
633 Notes
634 -----
635 Using a CovFit object, selects from (i, j) and returns
636 mu*gain, self.cov[:, i, j]*gain**2 model*gain**2, and self.sqrtW/gain**2
637 in electrons or ADU if unitsElectrons=False.
638 """
639 if unitsElectrons:
640 gain = self.getGain()
641 else:
642 gain = 1.0
644 mu = self.mu*gain
645 covariance = self.cov[:, i, j]*(gain**2)
646 covarianceModel = self.evalCovModel()[:, i, j]*(gain**2)
647 weights = self.sqrtW[:, i, j]/(gain**2)
649 # select data used for the fit:
650 mask = self.getMaskCov(i, j)
651 if returnMasked:
652 weights = weights[mask]
653 covarianceModel = covarianceModel[mask]
654 mu = mu[mask]
655 covariance = covariance[mask]
657 if divideByMu:
658 covariance /= mu
659 covarianceModel /= mu
660 weights *= mu
662 return mu, covariance, covarianceModel, weights, mask
664 def __call__(self, params):
665 self.setParamValues(params)
666 chi2 = self.chi2()
668 return chi2