24 from scipy.signal
import fftconvolve
25 from scipy.optimize
import leastsq
26 from .astierCovFitParameters
import FitParameters
34 """Make covariances array from tuple.
38 inputTuple: `numpy.ndarray`
39 Structured array with rows with at least
40 (mu, afwVar, cov, var, i, j, npix), where:
42 mu : 0.5*(m1 + m2), where:
43 mu1: mean value of flat1
44 mu2: mean value of flat2
45 afwVar: variance of difference flat, calculated with afw
46 cov: covariance value at lag(i, j)
47 var: variance(covariance value at lag(0, 0))
50 npix: number of pixels used for covariance calculation.
52 maxRangeFromTuple: `int`
53 Maximum range to select from tuple.
58 Covariance arrays, indexed by mean signal mu.
61 Variance arrays, indexed by mean signal mu.
64 List of mean signal values.
69 The input tuple should contain the following rows:
70 (mu, cov, var, i, j, npix), with one entry per lag, and image pair.
71 Different lags(i.e. different i and j) from the same
72 image pair have the same values of mu1 and mu2. When i==j==0, cov
75 If the input tuple contains several video channels, one should
76 select the data of a given channel *before* entering this
77 routine, as well as apply(e.g.) saturation cuts.
79 The routine returns cov[k_mu, j, i], vcov[(same indices)], and mu[k]
80 where the first index of cov matches the one in mu.
82 This routine implements the loss of variance due to
83 clipping cuts when measuring variances and covariance, but this should happen inside
84 the measurement code, where the cuts are readily available.
87 if maxRangeFromTuple
is not None:
88 cut = (inputTuple[
'i'] < maxRangeFromTuple) & (inputTuple[
'j'] < maxRangeFromTuple)
89 cutTuple = inputTuple[cut]
93 muTemp = cutTuple[
'mu']
94 ind = np.argsort(muTemp)
96 cutTuple = cutTuple[ind]
99 xx = np.hstack(([mu[0]], mu))
100 delta = xx[1:] - xx[:-1]
101 steps, = np.where(delta > 0)
102 ind = np.zeros_like(mu, dtype=int)
106 muVals = np.array(np.unique(mu))
107 i = cutTuple[
'i'].astype(int)
108 j = cutTuple[
'j'].astype(int)
109 c = 0.5*cutTuple[
'cov']
111 v = 0.5*cutTuple[
'var']
113 cov = np.ndarray((len(muVals), np.max(i)+1, np.max(j)+1))
114 var = np.zeros_like(cov)
116 var[ind, i, j] = v**2/n
119 return cov, var, muVals
123 """ Copy array over 4 quadrants prior to convolution.
127 inputarray: `numpy.array`
128 Input array to symmetrize.
136 targetShape =
list(inputArray.shape)
137 r1, r2 = inputArray.shape[-1], inputArray.shape[-2]
138 targetShape[-1] = 2*r1-1
139 targetShape[-2] = 2*r2-1
140 aSym = np.ndarray(tuple(targetShape))
141 aSym[..., r2-1:, r1-1:] = inputArray
142 aSym[..., r2-1:, r1-1::-1] = inputArray
143 aSym[..., r2-1::-1, r1-1::-1] = inputArray
144 aSym[..., r2-1::-1, r1-1:] = inputArray
150 """A class to fit the models in Astier+19 to flat covariances.
152 This code implements the model(and the fit thereof) described in
153 Astier+19: https://arxiv.org/pdf/1905.08677.pdf
157 meanSignals : `list`[`float`]
158 List with means of the difference image of two flats,
159 for a particular amplifier in the detector.
161 covariances : `list`[`numpy.array`]
162 List with 2D covariance arrays at a given mean signal.
164 covsSqrtWeights : `list`[`numpy.array`]
165 List with 2D arrays with weights from `vcov as defined in
166 `makeCovArray`: weight = 1/sqrt(vcov).
168 maxRangeFromTuple: `int`, optional
169 Maximum range to select from tuple.
171 meanSignalMask: `list`[`bool`], optional
172 Mask of mean signal 1D array. Use all entries if empty.
175 def __init__(self, meanSignals, covariances, covsSqrtWeights, maxRangeFromTuple=8, meanSignalsMask=[]):
176 assert (len(meanSignals) == len(covariances))
177 assert (len(covariances) == len(covsSqrtWeights))
178 if len(meanSignalsMask) == 0:
179 meanSignalsMask = np.repeat(
True, len(meanSignals))
180 self.
mumu = meanSignals[meanSignalsMask]
181 self.
covcov = np.nan_to_num(covariances)[meanSignalsMask]
183 self.
sqrtWsqrtW = np.nan_to_num(covsSqrtWeights)[meanSignalsMask]
184 self.
rr = maxRangeFromTuple
185 self.
loggerlogger = lsstLog.Log.getDefaultLogger()
188 """Make a copy of params"""
189 cop = copy.deepcopy(self)
191 if hasattr(self,
'params'):
196 """ Performs a crude parabolic fit of the data in order to start
197 the full fit close to the solution.
200 lenA = self.
rr*self.
rr
203 self.
paramsparams[
'gain'] = 1.
205 self.
paramsparams[
'c'].fix(val=0.)
207 a = self.
paramsparams[
'a'].full.reshape(self.
rr, self.
rr)
208 noise = self.
paramsparams[
'noise'].full.reshape(self.
rr, self.
rr)
209 gain = self.
paramsparams[
'gain'].full[0]
218 for i
in range(self.
rr):
219 for j
in range(self.
rr):
221 parsFit = np.polyfit(self.
mumu, self.
covcov[:, i, j] - model[:, i, j],
222 2, w=self.
sqrtWsqrtW[:, i, j])
224 a[i, j] += parsFit[0]
225 noise[i, j] += parsFit[2]*gain*gain
227 gain = 1./(1/gain+parsFit[1])
228 self.
paramsparams[
'gain'].full[0] = gain
229 chi2 = self.
chi2chi2()
237 """Return an array of free parameter values (it is a copy)."""
238 return self.
paramsparams.free + 0.
241 """Set parameter values."""
242 self.
paramsparams.free = p
246 """Computes full covariances model (Eq. 20 of Astier+19).
250 mu: `numpy.array`, optional
251 List of mean signals.
255 covModel: `numpy.array`
260 By default, computes the covModel for the mu's stored(self.mu).
262 Returns cov[Nmu, self.r, self.r]. The variance for the PTC is cov[:, 0, 0].
263 mu and cov are in ADUs and ADUs squared. To use electrons for both,
264 the gain should be set to 1. This routine implements the model in Astier+19 (1905.08677).
266 The parameters of the full model for C_ij(mu) ("C_ij" and "mu" in ADU^2 and ADU, respectively)
267 in Astier+19 (Eq. 20) are:
269 "a" coefficients (r by r matrix), units: 1/e
270 "b" coefficients (r by r matrix), units: 1/e
271 noise matrix (r by r matrix), units: e^2
274 "b" appears in Eq. 20 only through the "ab" combination, which is defined in this code as "c=ab".
276 sa = (self.
rr, self.
rr)
277 a = self.
paramsparams[
'a'].full.reshape(sa)
278 c = self.
paramsparams[
'c'].full.reshape(sa)
279 gain = self.
paramsparams[
'gain'].full[0]
280 noise = self.
paramsparams[
'noise'].full.reshape(sa)
282 aEnlarged = np.zeros((int(sa[0]*1.5)+1, int(sa[1]*1.5)+1))
283 aEnlarged[0:sa[0], 0:sa[1]] = a
286 cEnlarged = np.zeros((int(sa[0]*1.5)+1, int(sa[1]*1.5)+1))
287 cEnlarged[0:sa[0], 0:sa[1]] = c
289 a2 = fftconvolve(aSym, aSym, mode=
'same')
290 a3 = fftconvolve(a2, aSym, mode=
'same')
291 ac = fftconvolve(aSym, cSym, mode=
'same')
292 (xc, yc) = np.unravel_index(np.abs(aSym).argmax(), a2.shape)
294 a1 = a[np.newaxis, :, :]
295 a2 = a2[np.newaxis, xc:xc + range, yc:yc + range]
296 a3 = a3[np.newaxis, xc:xc + range, yc:yc + range]
297 ac = ac[np.newaxis, xc:xc + range, yc:yc + range]
298 c1 = c[np.newaxis, ::]
302 bigMu = mu[:, np.newaxis, np.newaxis]*gain
304 covModel = (bigMu/(gain*gain)*(a1*bigMu+2./3.*(bigMu*bigMu)*(a2 + c1)
305 + (1./3.*a3 + 5./6.*ac)*(bigMu*bigMu*bigMu)) + noise[np.newaxis, :, :]/gain**2)
307 covModel[:, 0, 0] += mu/gain
312 """'a' matrix from Astier+19(e.g., Eq. 20)"""
313 return self.
paramsparams[
'a'].full.reshape(self.
rr, self.
rr)
316 """'b' matrix from Astier+19(e.g., Eq. 20)"""
317 return self.
paramsparams[
'c'].full.reshape(self.
rr, self.
rr)/self.
getAgetA()
320 """'c'='ab' matrix from Astier+19(e.g., Eq. 20)"""
321 return np.array(self.
paramsparams[
'c'].full.reshape(self.
rr, self.
rr))
323 def _getCovParams(self, what):
324 """Get covariance matrix of parameters from fit"""
325 indices = self.
paramsparams[what].indexof()
326 i1 = indices[:, np.newaxis]
327 i2 = indices[np.newaxis, :]
335 """Get covariance matrix of "a" coefficients from fit"""
337 cova = self.
_getCovParams_getCovParams(
'a').reshape((self.
rr, self.
rr, self.
rr, self.
rr))
343 """Square root of diagonal of the parameter covariance of the fitted "a" matrix"""
345 sigA = np.sqrt(self.
_getCovParams_getCovParams(
'a').diagonal()).reshape((self.
rr, self.
rr))
351 """Get covariance matrix of "a" coefficients from fit
355 aval = self.
getAgetA().flatten()
356 factor = np.outer(aval, aval)
358 return covb.reshape((self.
rr, self.
rr, self.
rr, self.
rr))
361 """Get covariance matrix of "c" coefficients from fit"""
363 return cova.reshape((self.
rr, self.
rr, self.
rr, self.
rr))
366 """Get error on fitted gain parameter"""
368 gainErr = np.sqrt(self.
_getCovParams_getCovParams(
'gain')[0][0])
374 """Get covariances of noise matrix from fit"""
376 return covNoise.reshape((self.
rr, self.
rr, self.
rr, self.
rr))
379 """Square root of diagonal of the parameter covariance of the fitted "noise" matrix"""
382 noise = np.sqrt(covNoise.diagonal()).reshape((self.
rr, self.
rr))
388 """Get gain (e/ADU)"""
389 return self.
paramsparams[
'gain'].full[0]
392 """Get readout noise (e^2)"""
393 return self.
paramsparams[
'noise'].full[0]
396 """Get error on readout noise parameter"""
397 ronSqrt = np.sqrt(np.fabs(self.
getRongetRon()))
400 ronErr = 0.5*(noiseSigma/np.fabs(self.
getRongetRon()))*ronSqrt
406 """Get noise matrix"""
407 return self.
paramsparams[
'noise'].full.reshape(self.
rr, self.
rr)
410 """Get mask of Cov[i,j]"""
411 weights = self.
sqrtWsqrtW[:, i, j]
416 """Set "a" and "b" coeffcients forfull Astier+19 model (Eq. 20). "c=a*b"."""
417 self.
paramsparams[
'a'].full = a.flatten()
418 self.
paramsparams[
'c'].full = a.flatten()*b.flatten()
422 """Calculate weighted chi2 of full-model fit."""
426 """Weighted residuals.
431 c = CovFit(meanSignals, covariances, covsSqrtWeights)
433 coeffs, cov, _, mesg, ierr = leastsq(c.weightedRes, c.getParamValues(), full_output=True)
435 if params
is not None:
437 covModel = np.nan_to_num(self.
evalCovModelevalCovModel())
438 weightedRes = (covModel-self.
covcov)*self.
sqrtWsqrtW
440 return weightedRes.flatten()
443 """Fit measured covariances to full model in Astier+19 (Eq. 20)
448 Initial parameters of the fit.
449 len(pInit) = #entries(a) + #entries(c) + #entries(noise) + 1
450 len(pInit) = r^2 + r^2 + r^2 + 1, where "r" is the maximum lag considered for the
451 covariances calculation, and the extra "1" is the gain.
452 If "b" is 0, then "c" is 0, and len(pInit) will have r^2 fewer entries.
457 Fit parameters (see "Notes" below).
461 The parameters of the full model for C_ij(mu) ("C_ij" and "mu" in ADU^2 and ADU, respectively)
462 in Astier+19 (Eq. 20) are:
464 "a" coefficients (r by r matrix), units: 1/e
465 "b" coefficients (r by r matrix), units: 1/e
466 noise matrix (r by r matrix), units: e^2
469 "b" appears in Eq. 20 only through the "ab" combination, which is defined in this code as "c=ab".
474 params, paramsCov, _, mesg, ierr = leastsq(self.
weightedResweightedRes, pInit, full_output=
True)
479 def getFitData(self, i, j, divideByMu=False, unitsElectrons=False, returnMasked=False):
480 """Get measured signal and covariance, cov model, weigths, and mask at covariance lag (i, j).
485 Lag for covariance matrix.
488 Lag for covariance matrix.
490 divideByMu: `bool`, optional
491 Divide covariance, model, and weights by signal mu?
493 unitsElectrons : `bool`, optional
494 mu, covariance, and model are in ADU (or powers of ADU) If tthis parameter is true, these are
495 multiplied by the adequte factors of the gain to return quantities in electrons
496 (or powers of electrons).
498 returnMasked : `bool`, optional
499 Use mask (based on weights) in returned arrays (mu, covariance, and model)?
504 list of signal values (mu).
506 covariance: `numpy.array`
507 Covariance arrays, indexed by mean signal mu (self.cov[:, i, j]).
509 covarianceModel: `numpy.array`
510 Covariance model (model).
512 weights: `numpy.array`
515 mask : `numpy.array`, optional
516 Boolean mask of the covariance at (i,j).
520 Using a CovFit object, selects from (i, j) and returns
521 mu*gain, self.cov[:, i, j]*gain**2 model*gain**2, and self.sqrtW/gain**2
522 in electrons or ADU if unitsElectrons=False.
530 covariance = self.
covcov[:, i, j]*(gain**2)
531 covarianceModel = self.
evalCovModelevalCovModel()[:, i, j]*(gain**2)
532 weights = self.
sqrtWsqrtW[:, i, j]/(gain**2)
537 weights = weights[mask]
538 covarianceModel = covarianceModel[mask]
540 covariance = covariance[mask]
544 covarianceModel /= mu
547 return mu, covariance, covarianceModel, weights, mask
def _getCovParams(self, what)
def weightedRes(self, params=None)
def evalCovModel(self, mu=None)
def getMaskCov(self, i, j)
def setParamValues(self, p)
def getFitData(self, i, j, divideByMu=False, unitsElectrons=False, returnMasked=False)
def fitFullModel(self, pInit=None)
def __init__(self, meanSignals, covariances, covsSqrtWeights, maxRangeFromTuple=8, meanSignalsMask=[])
daf::base::PropertyList * list
def symmetrize(inputArray)
def makeCovArray(inputTuple, maxRangeFromTuple=8)