25 from scipy.stats
import median_abs_deviation
as mad
26 from scipy.signal
import fftconvolve
27 from scipy.optimize
import leastsq
28 from .astierCovFitParameters
import FitParameters
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).
43 Covariance model from Eq. 20 in Astier+19.
46 Mean signal in electrons
53 aCoeffsOld: `numpy.array`
54 Slope of cov/var at a given flux mu in electrons.
58 Returns the "a" array, computed this way, to be compared to the actual a_array from the full model
61 covModel = np.array(covModel)
62 var = covModel[0, 0, 0]
64 return covModel[0, :, :]/(var*muEl)
68 """Make covariances array from tuple.
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))
81 npix: number of pixels used for covariance calculation.
83 maxRangeFromTuple: `int`
84 Maximum range to select from tuple.
89 Covariance arrays, indexed by mean signal mu.
92 Variance arrays, indexed by mean signal mu.
95 List of mean signal values.
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
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.
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.
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.
118 if maxRangeFromTuple
is not None:
119 cut = (inputTuple[
'i'] < maxRangeFromTuple) & (inputTuple[
'j'] < maxRangeFromTuple)
120 cutTuple = inputTuple[cut]
122 cutTuple = inputTuple
124 muTemp = cutTuple[
'mu']
125 ind = np.argsort(muTemp)
127 cutTuple = cutTuple[ind]
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)
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']
142 v = 0.5*cutTuple[
'var']
144 cov = np.ndarray((len(muVals), np.max(i)+1, np.max(j)+1))
145 var = np.zeros_like(cov)
147 var[ind, i, j] = v**2/n
156 return cov, var, muVals
160 """ Copy array over 4 quadrants prior to convolution.
164 inputarray: `numpy.array`
165 Input array to symmetrize.
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
187 """A class to calculate 2D polynomials"""
194 self.coeff, _, rank, _ = np.linalg.lstsq(G, z.ravel())
196 self.coeff, _, rank, _ = np.linalg.lstsq((w.ravel()*G.T).T, z.ravel()*w.ravel())
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
210 return np.dot(G, self.coeff)
214 """A class to fit the models in Astier+19 to flat covariances.
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.
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))
232 npix: number of pixels used for covariance calculation.
234 maxRangeFromTuple: `int`
235 Maximum range to select from tuple.
238 def __init__(self, inputTuple, maxRangeFromTuple=8):
244 """Subtract a background/offset to the measured covariances.
249 Maximum lag considered
252 First lag from where to start the offset subtraction.
255 Degree of 2D polynomial to fit to covariance to define offse to be subtracted.
257 assert(startLag < self.
r)
258 for k
in range(len(self.
mu)):
260 w = self.
sqrtW[k, ...] + 0.
262 i, j = np.meshgrid(range(sh[0]), range(sh[1]), indexing=
'ij')
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
276 """Select signal level based on max average signal in ADU"""
278 index = self.
mu < maxMu
280 self.
mu = self.
mu[:k]
281 self.
cov = self.
cov[:k, ...]
288 """Select signal level based on max average signal in electrons"""
290 kill = (self.
mu*g > maxMuEl)
291 self.
sqrtW[kill, :, :] = 0
296 """Make a copy of params"""
297 cop = copy.deepcopy(self)
299 if hasattr(self,
'params'):
304 """ Performs a crude parabolic fit of the data in order to start
305 the full fit close to the solution.
313 self.
params[
'c'].fix(val=0.)
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]
326 for i
in range(self.
r):
327 for j
in range(self.
r):
329 parsFit = np.polyfit(self.
mu, self.
cov[:, i, j] - model[:, i, j],
330 2, w=self.
sqrtW[:, i, j])
332 a[i, j] += parsFit[0]
333 noise[i, j] += parsFit[2]*gain*gain
335 gain = 1./(1/gain+parsFit[1])
336 self.
params[
'gain'].full[0] = gain
345 """Return an array of free parameter values (it is a copy)."""
346 return self.
params.free + 0.
349 """Set parameter values."""
354 """Computes full covariances model (Eq. 20 of Astier+19).
358 mu: `numpy.array`, optional
359 List of mean signals.
363 covModel: `numpy.array`
368 By default, computes the covModel for the mu's stored(self.mu).
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).
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:
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
382 "b" appears in Eq. 20 only through the "ab" combination, which is defined in this code as "c=ab".
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)
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
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
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)
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, ::]
410 bigMu = mu[:, np.newaxis, np.newaxis]*gain
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)
415 covModel[:, 0, 0] += mu/gain
420 """'a' matrix from Astier+19(e.g., Eq. 20)"""
421 return self.
params[
'a'].full.reshape(self.
r, self.
r)
424 """'b' matrix from Astier+19(e.g., Eq. 20)"""
425 return self.
params[
'c'].full.reshape(self.
r, self.
r)/self.
getA()
428 """'c'='ab' matrix from Astier+19(e.g., Eq. 20)"""
429 return np.array(self.
params[
'c'].full.reshape(self.
r, self.
r))
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, :]
443 """Get covariance matrix of "a" coefficients from fit"""
451 """Square root of diagonal of the parameter covariance of the fitted "a" matrix"""
453 sigA = np.sqrt(self.
_getCovParams(
'a').diagonal()).reshape((self.
r, self.
r))
459 """Get covariance matrix of "a" coefficients from fit
463 aval = self.
getA().flatten()
464 factor = np.outer(aval, aval)
466 return covb.reshape((self.
r, self.
r, self.
r, self.
r))
469 """Get covariance matrix of "c" coefficients from fit"""
471 return cova.reshape((self.
r, self.
r, self.
r, self.
r))
474 """Get error on fitted gain parameter"""
482 """Get covariances of noise matrix from fit"""
484 return covNoise.reshape((self.
r, self.
r, self.
r, self.
r))
487 """Square root of diagonal of the parameter covariance of the fitted "noise" matrix"""
490 noise = np.sqrt(covNoise.diagonal()).reshape((self.
r, self.
r))
496 """Get gain (e/ADU)"""
497 return self.
params[
'gain'].full[0]
500 """Get readout noise (e^2)"""
501 return self.
params[
'noise'].full[0]
504 """Get error on readout noise parameter"""
505 ronSqrt = np.sqrt(np.fabs(self.
getRon()))
508 ronErr = 0.5*(noiseSigma/np.fabs(self.
getRon()))*ronSqrt
514 """Get noise matrix"""
515 return self.
params[
'noise'].full.reshape(self.
r, self.
r)
518 """Get mask of Cov[i,j]"""
519 weights = self.
sqrtW[:, i, j]
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()
530 """Calculate weighte chi2 of full-model fit."""
534 """To be used in weightedRes"""
535 if params
is not None:
538 return((covModel-self.
cov)*self.
sqrtW)
541 """Weighted residuas.
548 coeffs, cov, _, mesg, ierr = leastsq(c.weightedRes, c.getParamValues(), full_output=True)
550 return self.
wres(params).flatten()
553 """Fit measured covariances to full model in Astier+19 (Eq. 20)
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.
564 nSigma : `float`, optional
565 Sigma cut to get rid of outliers.
567 maxFitIter : `int`, optional
568 Number of iterations for full model fit.
573 Fit parameters (see "Notes" below).
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:
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
585 "b" appears in Eq. 20 only through the "ab" combination, which is defined in this code as "c=ab".
592 while nOutliers != 0:
593 params, paramsCov, _, mesg, ierr = leastsq(self.
weightedRes, pInit, full_output=
True)
596 sig = mad(wres[wres != 0], scale=
'normal')
597 mask = (np.abs(wres) > (nSigma*sig))
598 self.
sqrtW.flat[mask] = 0
599 nOutliers = mask.sum()
601 if counter == maxFitIter:
604 if ierr
not in [1, 2, 3, 4]:
605 raise RuntimeError(
"Minimization failed: " + mesg)
611 """Number of degrees of freedom
615 mask.sum() - len(self.params.free): `int`
616 Number of usable pixels - number of parameters of fit.
618 mask = self.
sqrtW != 0
620 return mask.sum() - len(self.
params.free)
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).
628 Lag for covariance matrix.
631 Lag for covariance matrix.
633 divideByMu: `bool`, optional
634 Divide covariance, model, and weights by signal mu?
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).
641 returnMasked : `bool`, optional
642 Use mask (based on weights) in returned arrays (mu, covariance, and model)?
647 list of signal values (mu).
649 covariance: `numpy.array`
650 Covariance arrays, indexed by mean signal mu (self.cov[:, i, j]).
652 covarianceModel: `numpy.array`
653 Covariance model (model).
655 weights: `numpy.array`
658 mask : `numpy.array`, optional
659 Boolean mask of the covariance at (i,j).
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.
673 covariance = self.
cov[:, i, j]*(gain**2)
674 covarianceModel = self.
evalCovModel()[:, i, j]*(gain**2)
675 weights = self.
sqrtW[:, i, j]/(gain**2)
680 weights = weights[mask]
681 covarianceModel = covarianceModel[mask]
683 covariance = covariance[mask]
687 covarianceModel /= mu
690 return mu, covariance, covarianceModel, weights, mask