Coverage for python/lsst/cp/pipe/ptc/astierCovPtcFit.py: 15%
Shortcuts 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
Shortcuts 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
24from scipy.signal import fftconvolve
25from scipy.optimize import leastsq
26from .astierCovFitParameters import FitParameters
28import lsst.log as lsstLog
30__all__ = ["CovFit"]
33def makeCovArray(inputTuple, maxRangeFromTuple=8):
34 """Make covariances array from tuple.
36 Parameters
37 ----------
38 inputTuple : `numpy.ndarray`
39 Structured array with rows with at least
40 (mu, afwVar, cov, var, i, j, npix), where:
42 mu : `float`
43 0.5*(m1 + m2), where mu1 is the mean value of flat1
44 and mu2 is the mean value of flat2.
45 afwVar : `float`
46 Variance of difference flat, calculated with afw.
47 cov : `float`
48 Covariance value at lag(i, j)
49 var : `float`
50 Variance(covariance value at lag(0, 0))
51 i : `int`
52 Lag in dimension "x".
53 j : `int`
54 Lag in dimension "y".
55 npix : `int`
56 Number of pixels used for covariance calculation.
58 maxRangeFromTuple : `int`
59 Maximum range to select from tuple.
61 Returns
62 -------
63 cov : `numpy.array`
64 Covariance arrays, indexed by mean signal mu.
66 vCov : `numpy.array`
67 Variance arrays, indexed by mean signal mu.
69 muVals : `numpy.array`
70 List of mean signal values.
72 Notes
73 -----
75 The input tuple should contain the following rows:
76 (mu, cov, var, i, j, npix), with one entry per lag, and image pair.
77 Different lags(i.e. different i and j) from the same
78 image pair have the same values of mu1 and mu2. When i==j==0, cov
79 = var.
81 If the input tuple contains several video channels, one should
82 select the data of a given channel *before* entering this
83 routine, as well as apply(e.g.) saturation cuts.
85 The routine returns cov[k_mu, j, i], vcov[(same indices)], and mu[k]
86 where the first index of cov matches the one in mu.
88 This routine implements the loss of variance due to clipping cuts
89 when measuring variances and covariance, but this should happen
90 inside the measurement code, where the cuts are readily available.
91 """
92 if maxRangeFromTuple is not None:
93 cut = (inputTuple['i'] < maxRangeFromTuple) & (inputTuple['j'] < maxRangeFromTuple)
94 cutTuple = inputTuple[cut]
95 else:
96 cutTuple = inputTuple
97 # increasing mu order, so that we can group measurements with the same mu
98 muTemp = cutTuple['mu']
99 ind = np.argsort(muTemp)
101 cutTuple = cutTuple[ind]
102 # should group measurements on the same image pairs(same average)
103 mu = cutTuple['mu']
104 xx = np.hstack(([mu[0]], mu))
105 delta = xx[1:] - xx[:-1]
106 steps, = np.where(delta > 0)
107 ind = np.zeros_like(mu, dtype=int)
108 ind[steps] = 1
109 ind = np.cumsum(ind) # this acts as an image pair index.
110 # now fill the 3-d cov array(and variance)
111 muVals = np.array(np.unique(mu))
112 i = cutTuple['i'].astype(int)
113 j = cutTuple['j'].astype(int)
114 c = 0.5*cutTuple['cov']
115 n = cutTuple['npix']
116 v = 0.5*cutTuple['var']
117 # book and fill
118 cov = np.ndarray((len(muVals), np.max(i)+1, np.max(j)+1))
119 var = np.zeros_like(cov)
120 cov[ind, i, j] = c
121 var[ind, i, j] = v**2/n
122 var[:, 0, 0] *= 2 # var(v) = 2*v**2/N
124 return cov, var, muVals
127def symmetrize(inputArray):
128 """ Copy array over 4 quadrants prior to convolution.
130 Parameters
131 ----------
132 inputarray : `numpy.array`
133 Input array to symmetrize.
135 Returns
136 -------
137 aSym : `numpy.array`
138 Symmetrized array.
139 """
140 targetShape = list(inputArray.shape)
141 r1, r2 = inputArray.shape[-1], inputArray.shape[-2]
142 targetShape[-1] = 2*r1-1
143 targetShape[-2] = 2*r2-1
144 aSym = np.ndarray(tuple(targetShape))
145 aSym[..., r2-1:, r1-1:] = inputArray
146 aSym[..., r2-1:, r1-1::-1] = inputArray
147 aSym[..., r2-1::-1, r1-1::-1] = inputArray
148 aSym[..., r2-1::-1, r1-1:] = inputArray
150 return aSym
153class CovFit:
154 """A class to fit the models in Astier+19 to flat covariances.
156 This code implements the model(and the fit thereof) described in
157 Astier+19: https://arxiv.org/pdf/1905.08677.pdf
159 Parameters
160 ----------
161 meanSignals : `list` [`float`]
162 List with means of the difference image of two flats,
163 for a particular amplifier in the detector.
165 covariances : `list` [`numpy.array`, (N, N)]
166 List with 2D covariance arrays at a given mean signal.
168 covsSqrtWeights : `list` [`numpy.array`, (N, N)]
169 List with 2D arrays with weights from `vcov` as defined in
170 `makeCovArray`: weight = 1/sqrt(vcov).
172 maxRangeFromTuple : `int`, optional
173 Maximum range to select from tuple.
175 meanSignalMask : `list` [`bool`], optional
176 Mask of mean signal 1D array. Use all entries if empty.
177 """
179 def __init__(self, meanSignals, covariances, covsSqrtWeights, maxRangeFromTuple=8, meanSignalsMask=[]):
180 assert (len(meanSignals) == len(covariances))
181 assert (len(covariances) == len(covsSqrtWeights))
182 if len(meanSignalsMask) == 0:
183 meanSignalsMask = np.repeat(True, len(meanSignals))
184 self.mu = meanSignals[meanSignalsMask]
185 self.cov = np.nan_to_num(covariances)[meanSignalsMask]
186 # make it nan safe, replacing nan's with 0 in weights
187 self.sqrtW = np.nan_to_num(covsSqrtWeights)[meanSignalsMask]
188 self.r = maxRangeFromTuple
189 self.logger = lsstLog.Log.getDefaultLogger()
191 def copy(self):
192 """Make a copy of params"""
193 cop = copy.deepcopy(self)
194 # deepcopy does not work for FitParameters.
195 if hasattr(self, 'params'):
196 cop.params = self.params.copy()
197 return cop
199 def initFit(self):
200 """ Performs a crude parabolic fit of the data in order to start
201 the full fit close to the solution.
202 """
203 # number of parameters for 'a' array.
204 lenA = self.r*self.r
205 # define parameters: c corresponds to a*b in Astier+19 (Eq. 20).
206 self.params = FitParameters([('a', lenA), ('c', lenA), ('noise', lenA), ('gain', 1)])
207 self.params['gain'] = 1.
208 # c=0 in a first go.
209 self.params['c'].fix(val=0.)
210 # plumbing: extract stuff from the parameter structure
211 a = self.params['a'].full.reshape(self.r, self.r)
212 noise = self.params['noise'].full.reshape(self.r, self.r)
213 gain = self.params['gain'].full[0]
215 # iterate the fit to account for higher orders
216 # the chi2 does not necessarily go down, so one could
217 # stop when it increases
218 oldChi2 = 1e30
219 for _ in range(5):
220 model = np.nan_to_num(self.evalCovModel()) # this computes the full model.
221 # loop on lags
222 for i in range(self.r):
223 for j in range(self.r):
224 # fit a parabola for a given lag
225 parsFit = np.polyfit(self.mu, self.cov[:, i, j] - model[:, i, j],
226 2, w=self.sqrtW[:, i, j])
227 # model equation(Eq. 20) in Astier+19:
228 a[i, j] += parsFit[0]
229 noise[i, j] += parsFit[2]*gain*gain
230 if(i + j == 0):
231 gain = 1./(1/gain+parsFit[1])
232 self.params['gain'].full[0] = gain
233 chi2 = self.chi2()
234 if chi2 > oldChi2:
235 break
236 oldChi2 = chi2
238 return
240 def getParamValues(self):
241 """Return an array of free parameter values (it is a copy)."""
242 return self.params.free + 0.
244 def setParamValues(self, p):
245 """Set parameter values."""
246 self.params.free = p
247 return
249 def evalCovModel(self, mu=None):
250 """Computes full covariances model (Eq. 20 of Astier+19).
252 Parameters
253 ----------
254 mu : `numpy.array`, (N,), optional
255 List of mean signals.
257 Returns
258 -------
259 covModel : `numpy.array`, (N, M, M)
260 Covariances model.
262 Notes
263 -----
264 By default, computes the covModel for the mu's stored(self.mu).
266 Returns cov[Nmu, self.r, self.r]. The variance for the PTC is
267 cov[:, 0, 0]. mu and cov are in ADUs and ADUs squared. To use
268 electrons for both, the gain should be set to 1. This routine
269 implements the model in Astier+19 (1905.08677).
271 The parameters of the full model for C_ij(mu) ("C_ij" and "mu"
272 in ADU^2 and ADU, respectively) in Astier+19 (Eq. 20) are:
274 "a" coefficients (r by r matrix), units: 1/e
275 "b" coefficients (r by r matrix), units: 1/e
276 noise matrix (r by r matrix), units: e^2
277 gain, units: e/ADU
279 "b" appears in Eq. 20 only through the "ab" combination, which
280 is defined in this code as "c=ab".
281 """
282 sa = (self.r, self.r)
283 a = self.params['a'].full.reshape(sa)
284 c = self.params['c'].full.reshape(sa)
285 gain = self.params['gain'].full[0]
286 noise = self.params['noise'].full.reshape(sa)
287 # pad a with zeros and symmetrize
288 aEnlarged = np.zeros((int(sa[0]*1.5)+1, int(sa[1]*1.5)+1))
289 aEnlarged[0:sa[0], 0:sa[1]] = a
290 aSym = symmetrize(aEnlarged)
291 # pad c with zeros and symmetrize
292 cEnlarged = np.zeros((int(sa[0]*1.5)+1, int(sa[1]*1.5)+1))
293 cEnlarged[0:sa[0], 0:sa[1]] = c
294 cSym = symmetrize(cEnlarged)
295 a2 = fftconvolve(aSym, aSym, mode='same')
296 a3 = fftconvolve(a2, aSym, mode='same')
297 ac = fftconvolve(aSym, cSym, mode='same')
298 (xc, yc) = np.unravel_index(np.abs(aSym).argmax(), a2.shape)
299 range = self.r
300 a1 = a[np.newaxis, :, :]
301 a2 = a2[np.newaxis, xc:xc + range, yc:yc + range]
302 a3 = a3[np.newaxis, xc:xc + range, yc:yc + range]
303 ac = ac[np.newaxis, xc:xc + range, yc:yc + range]
304 c1 = c[np.newaxis, ::]
305 if mu is None:
306 mu = self.mu
307 # assumes that mu is 1d
308 bigMu = mu[:, np.newaxis, np.newaxis]*gain
309 # c(=a*b in Astier+19) also has a contribution to the last
310 # term, that is absent for now.
311 covModel = (bigMu/(gain*gain)*(a1*bigMu+2./3.*(bigMu*bigMu)*(a2 + c1)
312 + (1./3.*a3 + 5./6.*ac)*(bigMu*bigMu*bigMu)) + noise[np.newaxis, :, :]/gain**2)
313 # add the Poisson term, and the read out noise (variance)
314 covModel[:, 0, 0] += mu/gain
316 return covModel
318 def getA(self):
319 """'a' matrix from Astier+19(e.g., Eq. 20)"""
320 return self.params['a'].full.reshape(self.r, self.r)
322 def getB(self):
323 """'b' matrix from Astier+19(e.g., Eq. 20)"""
324 return self.params['c'].full.reshape(self.r, self.r)/self.getA()
326 def getC(self):
327 """'c'='ab' matrix from Astier+19(e.g., Eq. 20)"""
328 return np.array(self.params['c'].full.reshape(self.r, self.r))
330 def _getCovParams(self, what):
331 """Get covariance matrix of parameters from fit"""
332 indices = self.params[what].indexof()
333 i1 = indices[:, np.newaxis]
334 i2 = indices[np.newaxis, :]
335 if self.covParams is not None:
336 covp = self.covParams[i1, i2]
337 else:
338 covp = None
339 return covp
341 def getACov(self):
342 """Get covariance matrix of "a" coefficients from fit"""
343 if self._getCovParams('a') is not None:
344 cova = self._getCovParams('a').reshape((self.r, self.r, self.r, self.r))
345 else:
346 cova = None
347 return cova
349 def getASig(self):
350 """Square root of diagonal of the parameter covariance of the fitted
351 "a" matrix
352 """
353 if self._getCovParams('a') is not None:
354 sigA = np.sqrt(self._getCovParams('a').diagonal()).reshape((self.r, self.r))
355 else:
356 sigA = None
357 return sigA
359 def getBCov(self):
360 """Get covariance matrix of "a" coefficients from fit
361 b = c /a
362 """
363 covb = self._getCovParams('c')
364 aval = self.getA().flatten()
365 factor = np.outer(aval, aval)
366 covb /= factor
367 return covb.reshape((self.r, self.r, self.r, self.r))
369 def getCCov(self):
370 """Get covariance matrix of "c" coefficients from fit"""
371 cova = self._getCovParams('c')
372 return cova.reshape((self.r, self.r, self.r, self.r))
374 def getGainErr(self):
375 """Get error on fitted gain parameter"""
376 if self._getCovParams('gain') is not None:
377 gainErr = np.sqrt(self._getCovParams('gain')[0][0])
378 else:
379 gainErr = 0.0
380 return gainErr
382 def getNoiseCov(self):
383 """Get covariances of noise matrix from fit"""
384 covNoise = self._getCovParams('noise')
385 return covNoise.reshape((self.r, self.r, self.r, self.r))
387 def getNoiseSig(self):
388 """Square root of diagonal of the parameter covariance of the fitted
389 "noise" matrix
390 """
391 if self._getCovParams('noise') is not None:
392 covNoise = self._getCovParams('noise')
393 noise = np.sqrt(covNoise.diagonal()).reshape((self.r, self.r))
394 else:
395 noise = None
396 return noise
398 def getGain(self):
399 """Get gain (e/ADU)"""
400 return self.params['gain'].full[0]
402 def getRon(self):
403 """Get readout noise (e^2)"""
404 return self.params['noise'].full[0]
406 def getRonErr(self):
407 """Get error on readout noise parameter"""
408 ronSqrt = np.sqrt(np.fabs(self.getRon()))
409 if self.getNoiseSig() is not None:
410 noiseSigma = self.getNoiseSig()[0][0]
411 ronErr = 0.5*(noiseSigma/np.fabs(self.getRon()))*ronSqrt
412 else:
413 ronErr = np.nan
414 return ronErr
416 def getNoise(self):
417 """Get noise matrix"""
418 return self.params['noise'].full.reshape(self.r, self.r)
420 def getMaskCov(self, i, j):
421 """Get mask of Cov[i,j]"""
422 weights = self.sqrtW[:, i, j]
423 mask = weights != 0
424 return mask
426 def setAandB(self, a, b):
427 """Set "a" and "b" coeffcients forfull Astier+19 model
428 (Eq. 20). "c=a*b".
429 """
430 self.params['a'].full = a.flatten()
431 self.params['c'].full = a.flatten()*b.flatten()
432 return
434 def chi2(self):
435 """Calculate weighted chi2 of full-model fit."""
436 return(self.weightedRes()**2).sum()
438 def weightedRes(self, params=None):
439 """Weighted residuals.
441 Notes
442 -----
443 To be used via:
444 c = CovFit(meanSignals, covariances, covsSqrtWeights)
445 c.initFit()
446 coeffs, cov, _, mesg, ierr = leastsq(c.weightedRes,
447 c.getParamValues(),
448 full_output=True)
449 """
450 if params is not None:
451 self.setParamValues(params)
452 covModel = np.nan_to_num(self.evalCovModel())
453 weightedRes = (covModel-self.cov)*self.sqrtW
455 return weightedRes.flatten()
457 def fitFullModel(self, pInit=None):
458 """Fit measured covariances to full model in Astier+19 (Eq. 20)
460 Parameters
461 ----------
462 pInit : `list`
463 Initial parameters of the fit.
464 len(pInit) = #entries(a) + #entries(c) + #entries(noise) + 1
465 len(pInit) = r^2 + r^2 + r^2 + 1, where "r" is the maximum lag
466 considered for the covariances calculation, and the extra "1"
467 is the gain.
468 If "b" is 0, then "c" is 0, and len(pInit) will have r^2 fewer
469 entries.
471 Returns
472 -------
473 params : `np.array`, (N,)
474 Fit parameters (see "Notes" below).
476 Notes
477 -----
478 The parameters of the full model for C_ij(mu) ("C_ij" and "mu"
479 in ADU^2 and ADU, respectively) in Astier+19 (Eq. 20) are:
481 "a" coefficients (r by r matrix), units: 1/e
482 "b" coefficients (r by r matrix), units: 1/e
483 noise matrix (r by r matrix), units: e^2
484 gain, units: e/ADU
486 "b" appears in Eq. 20 only through the "ab" combination, which
487 is defined in this code as "c=ab".
488 """
489 if pInit is None:
490 pInit = self.getParamValues()
491 params, paramsCov, _, mesg, ierr = leastsq(self.weightedRes, pInit, full_output=True)
492 self.covParams = paramsCov
494 return params
496 def getFitData(self, i, j, divideByMu=False, unitsElectrons=False, returnMasked=False):
497 """Get measured signal and covariance, cov model, weigths, and mask at
498 covariance lag (i, j).
500 Parameters
501 ---------
502 i : `int`
503 Lag for covariance matrix.
505 j : `int`
506 Lag for covariance matrix.
508 divideByMu : `bool`, optional
509 Divide covariance, model, and weights by signal mu?
511 unitsElectrons : `bool`, optional
512 mu, covariance, and model are in ADU (or powers of ADU) If this
513 parameter is true, these are multiplied by the adequate
514 factors of the gain to return quantities in electrons (or
515 powers of electrons).
517 returnMasked : `bool`, optional
518 Use mask (based on weights) in returned arrays (mu,
519 covariance, and model)?
521 Returns
522 -------
523 mu : `numpy.array`, (N,)
524 list of signal values (mu).
526 covariance : `numpy.array`, (N,)
527 Covariance arrays, indexed by mean signal mu (self.cov[:, i, j]).
529 covarianceModel : `numpy.array`, (N,)
530 Covariance model (model).
532 weights : `numpy.array`, (N,)
533 Weights (self.sqrtW)
535 mask : `numpy.array`, (N,), optional
536 Boolean mask of the covariance at (i,j).
538 Notes
539 -----
540 Using a CovFit object, selects from (i, j) and returns
541 mu*gain, self.cov[:, i, j]*gain**2 model*gain**2, and
542 self.sqrtW/gain**2 in electrons or ADU if
543 unitsElectrons=False.
544 """
545 if unitsElectrons:
546 gain = self.getGain()
547 else:
548 gain = 1.0
550 mu = self.mu*gain
551 covariance = self.cov[:, i, j]*(gain**2)
552 covarianceModel = self.evalCovModel()[:, i, j]*(gain**2)
553 weights = self.sqrtW[:, i, j]/(gain**2)
555 # select data used for the fit:
556 mask = self.getMaskCov(i, j)
557 if returnMasked:
558 weights = weights[mask]
559 covarianceModel = covarianceModel[mask]
560 mu = mu[mask]
561 covariance = covariance[mask]
563 if divideByMu:
564 covariance /= mu
565 covarianceModel /= mu
566 weights *= mu
568 return mu, covariance, covarianceModel, weights, mask