LSSTApplications  21.0.0+1b62c9342b,21.0.0+45a059f35e,21.0.0-1-ga51b5d4+ceb9cf20a3,21.0.0-19-g7c7630f+a88ebbf2d9,21.0.0-2-g103fe59+3522cf3bc7,21.0.0-2-g1367e85+571a348718,21.0.0-2-g2909d54+45a059f35e,21.0.0-2-g45278ab+1b62c9342b,21.0.0-2-g4bc9b9f+35a70d5868,21.0.0-2-g5242d73+571a348718,21.0.0-2-g54e2caa+aa129c4686,21.0.0-2-g66bcc37+3caef57c29,21.0.0-2-g7f82c8f+6f9059e2fe,21.0.0-2-g8dde007+5d1b9cb3f5,21.0.0-2-g8f08a60+73884b2cf5,21.0.0-2-g973f35b+1d054a08b9,21.0.0-2-ga326454+6f9059e2fe,21.0.0-2-ga63a54e+3d2c655db6,21.0.0-2-gc738bc1+a567cb0f17,21.0.0-2-gde069b7+5a8f2956b8,21.0.0-2-ge17e5af+571a348718,21.0.0-2-ge712728+834f2a3ece,21.0.0-2-gecfae73+dfe6e80958,21.0.0-2-gfc62afb+571a348718,21.0.0-21-g006371a9+88174a2081,21.0.0-3-g4c5b185+7fd31a6834,21.0.0-3-g6d51c4a+3caef57c29,21.0.0-3-gaa929c8+55f5a6a5c9,21.0.0-3-gd222c45+afc8332dbe,21.0.0-3-gd5de2f2+3caef57c29,21.0.0-4-g3300ddd+1b62c9342b,21.0.0-4-g5873dc9+9a92674037,21.0.0-4-g8a80011+5955f0fd15,21.0.0-5-gb7080ec+8658c79ec4,21.0.0-5-gcff38f6+89f2a0074d,21.0.0-6-gd3283ba+55f5a6a5c9,21.0.0-8-g19111d86+2c4b0a9f47,21.0.0-9-g7bed000b9+c7d3cce47e,w.2021.03
LSSTDataManagementBasePackage
Public Member Functions | Public Attributes | List of all members
lsst.cp.pipe.astierCovPtcFit.CovFit Class Reference

Public Member Functions

def __init__ (self, inputTuple, maxRangeFromTuple=8, meanSignalMask=[])
 
def subtractDistantOffset (self, maxLag=8, startLag=5, polDegree=1)
 
def copy (self)
 
def initFit (self)
 
def getParamValues (self)
 
def setParamValues (self, p)
 
def evalCovModel (self, mu=None)
 
def getA (self)
 
def getB (self)
 
def getC (self)
 
def getACov (self)
 
def getASig (self)
 
def getBCov (self)
 
def getCCov (self)
 
def getGainErr (self)
 
def getNoiseCov (self)
 
def getNoiseSig (self)
 
def getGain (self)
 
def getRon (self)
 
def getRonErr (self)
 
def getNoise (self)
 
def getMaskCov (self, i, j)
 
def setAandB (self, a, b)
 
def chi2 (self)
 
def wres (self, params=None)
 
def weightedRes (self, params=None)
 
def fitFullModel (self, pInit=None)
 
def ndof (self)
 
def getFitData (self, i, j, divideByMu=False, unitsElectrons=False, returnMasked=False)
 
def __call__ (self, params)
 

Public Attributes

 mu
 
 sqrtW
 
 r
 
 logger
 
 maskMu
 
 cov
 
 vcov
 
 params
 
 covParams
 

Detailed Description

A class to fit the models in Astier+19 to flat covariances.

This code implements the model(and the fit thereof) described in
Astier+19: https://arxiv.org/pdf/1905.08677.pdf
For the time being it uses as input a numpy recarray (tuple with named tags) which
contains one row per covariance and per pair: see the routine makeCovArray.

Parameters
----------
inputTuple: `numpy.recarray`
    Tuple with at least (mu, cov, var, i, j, npix), where:
    mu : 0.5*(m1 + m2), where:
        mu1: mean value of flat1
        mu2: mean value of flat2
    cov: covariance value at lag(i, j)
    var: variance(covariance value at lag(0, 0))
    i: lag dimension
    j: lag dimension
    npix: number of pixels used for covariance calculation.

maxRangeFromTuple: `int`, optional
    Maximum range to select from tuple.

meanSignalMask: `list`[`bool`], optional
    Mask of mean signal 1D array. Use all entries if empty.

Definition at line 213 of file astierCovPtcFit.py.

Constructor & Destructor Documentation

◆ __init__()

def lsst.cp.pipe.astierCovPtcFit.CovFit.__init__ (   self,
  inputTuple,
  maxRangeFromTuple = 8,
  meanSignalMask = [] 
)

Definition at line 241 of file astierCovPtcFit.py.

241  def __init__(self, inputTuple, maxRangeFromTuple=8, meanSignalMask=[]):
242  self.cov, self.vcov, self.mu = makeCovArray(inputTuple, maxRangeFromTuple)
243  # make it nan safe, replacing nan's with 0 in weights
244  self.sqrtW = np.nan_to_num(1./np.sqrt(self.vcov))
245  self.r = self.cov.shape[1]
246  self.logger = lsstLog.Log.getDefaultLogger()
247  if len(meanSignalMask):
248  self.maskMu = meanSignalMask
249  else:
250  self.maskMu = np.repeat(True, len(self.mu))
251 

Member Function Documentation

◆ __call__()

def lsst.cp.pipe.astierCovPtcFit.CovFit.__call__ (   self,
  params 
)

Definition at line 664 of file astierCovPtcFit.py.

664  def __call__(self, params):
665  self.setParamValues(params)
666  chi2 = self.chi2()
667 
668  return chi2

◆ chi2()

def lsst.cp.pipe.astierCovPtcFit.CovFit.chi2 (   self)
Calculate weighted chi2 of full-model fit.

Definition at line 518 of file astierCovPtcFit.py.

518  def chi2(self):
519  """Calculate weighted chi2 of full-model fit."""
520  return(self.weightedRes()**2).sum()
521 

◆ copy()

def lsst.cp.pipe.astierCovPtcFit.CovFit.copy (   self)
Make a copy of params

Definition at line 284 of file astierCovPtcFit.py.

284  def copy(self):
285  """Make a copy of params"""
286  cop = copy.deepcopy(self)
287  # deepcopy does not work for FitParameters.
288  if hasattr(self, 'params'):
289  cop.params = self.params.copy()
290  return cop
291 

◆ evalCovModel()

def lsst.cp.pipe.astierCovPtcFit.CovFit.evalCovModel (   self,
  mu = None 
)
Computes full covariances model (Eq. 20 of Astier+19).

Parameters
----------
mu: `numpy.array`, optional
    List of mean signals.

Returns
-------
covModel: `numpy.array`
    Covariances model.

Notes
-----
By default, computes the covModel for the mu's stored(self.mu).

Returns cov[Nmu, self.r, self.r]. The variance for the PTC is cov[:, 0, 0].
mu and cov are in ADUs and ADUs squared. To use electrons for both,
the gain should be set to 1. This routine implements the model in Astier+19 (1905.08677).

The parameters of the full model for C_ij(mu) ("C_ij" and "mu" in ADU^2 and ADU, respectively)
in Astier+19 (Eq. 20) are:

"a" coefficients (r by r matrix), units: 1/e
"b" coefficients (r by r matrix), units: 1/e
noise matrix (r by r matrix), units: e^2
gain, units: e/ADU

"b" appears in Eq. 20 only through the "ab" combination, which is defined in this code as "c=ab".

Definition at line 342 of file astierCovPtcFit.py.

342  def evalCovModel(self, mu=None):
343  """Computes full covariances model (Eq. 20 of Astier+19).
344 
345  Parameters
346  ----------
347  mu: `numpy.array`, optional
348  List of mean signals.
349 
350  Returns
351  -------
352  covModel: `numpy.array`
353  Covariances model.
354 
355  Notes
356  -----
357  By default, computes the covModel for the mu's stored(self.mu).
358 
359  Returns cov[Nmu, self.r, self.r]. The variance for the PTC is cov[:, 0, 0].
360  mu and cov are in ADUs and ADUs squared. To use electrons for both,
361  the gain should be set to 1. This routine implements the model in Astier+19 (1905.08677).
362 
363  The parameters of the full model for C_ij(mu) ("C_ij" and "mu" in ADU^2 and ADU, respectively)
364  in Astier+19 (Eq. 20) are:
365 
366  "a" coefficients (r by r matrix), units: 1/e
367  "b" coefficients (r by r matrix), units: 1/e
368  noise matrix (r by r matrix), units: e^2
369  gain, units: e/ADU
370 
371  "b" appears in Eq. 20 only through the "ab" combination, which is defined in this code as "c=ab".
372  """
373  sa = (self.r, self.r)
374  a = self.params['a'].full.reshape(sa)
375  c = self.params['c'].full.reshape(sa)
376  gain = self.params['gain'].full[0]
377  noise = self.params['noise'].full.reshape(sa)
378  # pad a with zeros and symmetrize
379  aEnlarged = np.zeros((int(sa[0]*1.5)+1, int(sa[1]*1.5)+1))
380  aEnlarged[0:sa[0], 0:sa[1]] = a
381  aSym = symmetrize(aEnlarged)
382  # pad c with zeros and symmetrize
383  cEnlarged = np.zeros((int(sa[0]*1.5)+1, int(sa[1]*1.5)+1))
384  cEnlarged[0:sa[0], 0:sa[1]] = c
385  cSym = symmetrize(cEnlarged)
386  a2 = fftconvolve(aSym, aSym, mode='same')
387  a3 = fftconvolve(a2, aSym, mode='same')
388  ac = fftconvolve(aSym, cSym, mode='same')
389  (xc, yc) = np.unravel_index(np.abs(aSym).argmax(), a2.shape)
390  range = self.r
391  a1 = a[np.newaxis, :, :]
392  a2 = a2[np.newaxis, xc:xc + range, yc:yc + range]
393  a3 = a3[np.newaxis, xc:xc + range, yc:yc + range]
394  ac = ac[np.newaxis, xc:xc + range, yc:yc + range]
395  c1 = c[np.newaxis, ::]
396  if mu is None:
397  mu = self.mu
398  # assumes that mu is 1d
399  bigMu = mu[:, np.newaxis, np.newaxis]*gain
400  # c(=a*b in Astier+19) also has a contribution to the last term, that is absent for now.
401  covModel = (bigMu/(gain*gain)*(a1*bigMu+2./3.*(bigMu*bigMu)*(a2 + c1) +
402  (1./3.*a3 + 5./6.*ac)*(bigMu*bigMu*bigMu)) + noise[np.newaxis, :, :]/gain**2)
403  # add the Poisson term, and the read out noise (variance)
404  covModel[:, 0, 0] += mu/gain
405 
406  return covModel
407 

◆ fitFullModel()

def lsst.cp.pipe.astierCovPtcFit.CovFit.fitFullModel (   self,
  pInit = None 
)
Fit measured covariances to full model in Astier+19 (Eq. 20)

Parameters
----------
pInit : `list`
    Initial parameters of the fit.
    len(pInit) = #entries(a) + #entries(c) + #entries(noise) + 1
    len(pInit) = r^2 + r^2 + r^2 + 1, where "r" is the maximum lag considered for the
    covariances calculation, and the extra "1" is the gain.
    If "b" is 0, then "c" is 0, and len(pInit) will have r^2 fewer entries.

Returns
-------
params : `np.array`
    Fit parameters (see "Notes" below).

Notes
-----
The parameters of the full model for C_ij(mu) ("C_ij" and "mu" in ADU^2 and ADU, respectively)
in Astier+19 (Eq. 20) are:

    "a" coefficients (r by r matrix), units: 1/e
    "b" coefficients (r by r matrix), units: 1/e
    noise matrix (r by r matrix), units: e^2
    gain, units: e/ADU

"b" appears in Eq. 20 only through the "ab" combination, which is defined in this code as "c=ab".

Definition at line 544 of file astierCovPtcFit.py.

544  def fitFullModel(self, pInit=None):
545  """Fit measured covariances to full model in Astier+19 (Eq. 20)
546 
547  Parameters
548  ----------
549  pInit : `list`
550  Initial parameters of the fit.
551  len(pInit) = #entries(a) + #entries(c) + #entries(noise) + 1
552  len(pInit) = r^2 + r^2 + r^2 + 1, where "r" is the maximum lag considered for the
553  covariances calculation, and the extra "1" is the gain.
554  If "b" is 0, then "c" is 0, and len(pInit) will have r^2 fewer entries.
555 
556  Returns
557  -------
558  params : `np.array`
559  Fit parameters (see "Notes" below).
560 
561  Notes
562  -----
563  The parameters of the full model for C_ij(mu) ("C_ij" and "mu" in ADU^2 and ADU, respectively)
564  in Astier+19 (Eq. 20) are:
565 
566  "a" coefficients (r by r matrix), units: 1/e
567  "b" coefficients (r by r matrix), units: 1/e
568  noise matrix (r by r matrix), units: e^2
569  gain, units: e/ADU
570 
571  "b" appears in Eq. 20 only through the "ab" combination, which is defined in this code as "c=ab".
572  """
573 
574  if pInit is None:
575  pInit = self.getParamValues()
576 
577  params, paramsCov, _, mesg, ierr = leastsq(self.weightedRes, pInit, full_output=True)
578  self.covParams = paramsCov
579 
580  return params
581 

◆ getA()

def lsst.cp.pipe.astierCovPtcFit.CovFit.getA (   self)
'a' matrix from Astier+19(e.g., Eq. 20)

Definition at line 408 of file astierCovPtcFit.py.

408  def getA(self):
409  """'a' matrix from Astier+19(e.g., Eq. 20)"""
410  return self.params['a'].full.reshape(self.r, self.r)
411 

◆ getACov()

def lsst.cp.pipe.astierCovPtcFit.CovFit.getACov (   self)
Get covariance matrix of "a" coefficients from fit

Definition at line 431 of file astierCovPtcFit.py.

431  def getACov(self):
432  """Get covariance matrix of "a" coefficients from fit"""
433  if self._getCovParams('a') is not None:
434  cova = self._getCovParams('a').reshape((self.r, self.r, self.r, self.r))
435  else:
436  cova = None
437  return cova
438 

◆ getASig()

def lsst.cp.pipe.astierCovPtcFit.CovFit.getASig (   self)
Square root of diagonal of the parameter covariance of the fitted "a" matrix

Definition at line 439 of file astierCovPtcFit.py.

439  def getASig(self):
440  """Square root of diagonal of the parameter covariance of the fitted "a" matrix"""
441  if self._getCovParams('a') is not None:
442  sigA = np.sqrt(self._getCovParams('a').diagonal()).reshape((self.r, self.r))
443  else:
444  sigA = None
445  return sigA
446 

◆ getB()

def lsst.cp.pipe.astierCovPtcFit.CovFit.getB (   self)
'b' matrix from Astier+19(e.g., Eq. 20)

Definition at line 412 of file astierCovPtcFit.py.

412  def getB(self):
413  """'b' matrix from Astier+19(e.g., Eq. 20)"""
414  return self.params['c'].full.reshape(self.r, self.r)/self.getA()
415 

◆ getBCov()

def lsst.cp.pipe.astierCovPtcFit.CovFit.getBCov (   self)
Get covariance matrix of "a" coefficients from fit
b = c /a

Definition at line 447 of file astierCovPtcFit.py.

447  def getBCov(self):
448  """Get covariance matrix of "a" coefficients from fit
449  b = c /a
450  """
451  covb = self._getCovParams('c')
452  aval = self.getA().flatten()
453  factor = np.outer(aval, aval)
454  covb /= factor
455  return covb.reshape((self.r, self.r, self.r, self.r))
456 

◆ getC()

def lsst.cp.pipe.astierCovPtcFit.CovFit.getC (   self)
'c'='ab' matrix from Astier+19(e.g., Eq. 20)

Definition at line 416 of file astierCovPtcFit.py.

416  def getC(self):
417  """'c'='ab' matrix from Astier+19(e.g., Eq. 20)"""
418  return np.array(self.params['c'].full.reshape(self.r, self.r))
419 

◆ getCCov()

def lsst.cp.pipe.astierCovPtcFit.CovFit.getCCov (   self)
Get covariance matrix of "c" coefficients from fit

Definition at line 457 of file astierCovPtcFit.py.

457  def getCCov(self):
458  """Get covariance matrix of "c" coefficients from fit"""
459  cova = self._getCovParams('c')
460  return cova.reshape((self.r, self.r, self.r, self.r))
461 

◆ getFitData()

def lsst.cp.pipe.astierCovPtcFit.CovFit.getFitData (   self,
  i,
  j,
  divideByMu = False,
  unitsElectrons = False,
  returnMasked = False 
)
Get measured signal and covariance, cov model, weigths, and mask at covariance lag (i, j).

Parameters
---------
i: `int`
    Lag for covariance matrix.

j: `int`
    Lag for covariance matrix.

divideByMu: `bool`, optional
    Divide covariance, model, and weights by signal mu?

unitsElectrons : `bool`, optional
    mu, covariance, and model are in ADU (or powers of ADU) If tthis parameter is true, these are
    multiplied by the adequte factors of the gain to return quantities in electrons
    (or powers of electrons).

returnMasked : `bool`, optional
    Use mask (based on weights) in returned arrays (mu, covariance, and model)?

Returns
-------
mu: `numpy.array`
    list of signal values (mu).

covariance: `numpy.array`
    Covariance arrays, indexed by mean signal mu (self.cov[:, i, j]).

covarianceModel: `numpy.array`
    Covariance model (model).

weights: `numpy.array`
    Weights (self.sqrtW)

mask : `numpy.array`, optional
    Boolean mask of the covariance at (i,j).

Notes
-----
Using a CovFit object, selects from (i, j) and returns
mu*gain, self.cov[:, i, j]*gain**2 model*gain**2, and self.sqrtW/gain**2
in electrons or ADU if unitsElectrons=False.

Definition at line 594 of file astierCovPtcFit.py.

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).
596 
597  Parameters
598  ---------
599  i: `int`
600  Lag for covariance matrix.
601 
602  j: `int`
603  Lag for covariance matrix.
604 
605  divideByMu: `bool`, optional
606  Divide covariance, model, and weights by signal mu?
607 
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).
612 
613  returnMasked : `bool`, optional
614  Use mask (based on weights) in returned arrays (mu, covariance, and model)?
615 
616  Returns
617  -------
618  mu: `numpy.array`
619  list of signal values (mu).
620 
621  covariance: `numpy.array`
622  Covariance arrays, indexed by mean signal mu (self.cov[:, i, j]).
623 
624  covarianceModel: `numpy.array`
625  Covariance model (model).
626 
627  weights: `numpy.array`
628  Weights (self.sqrtW)
629 
630  mask : `numpy.array`, optional
631  Boolean mask of the covariance at (i,j).
632 
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
643 
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)
648 
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]
656 
657  if divideByMu:
658  covariance /= mu
659  covarianceModel /= mu
660  weights *= mu
661 
662  return mu, covariance, covarianceModel, weights, mask
663 

◆ getGain()

def lsst.cp.pipe.astierCovPtcFit.CovFit.getGain (   self)
Get gain (e/ADU)

Definition at line 484 of file astierCovPtcFit.py.

484  def getGain(self):
485  """Get gain (e/ADU)"""
486  return self.params['gain'].full[0]
487 

◆ getGainErr()

def lsst.cp.pipe.astierCovPtcFit.CovFit.getGainErr (   self)
Get error on fitted gain parameter

Definition at line 462 of file astierCovPtcFit.py.

462  def getGainErr(self):
463  """Get error on fitted gain parameter"""
464  if self._getCovParams('gain') is not None:
465  gainErr = np.sqrt(self._getCovParams('gain')[0][0])
466  else:
467  gainErr = 0.0
468  return gainErr
469 

◆ getMaskCov()

def lsst.cp.pipe.astierCovPtcFit.CovFit.getMaskCov (   self,
  i,
  j 
)
Get mask of Cov[i,j]

Definition at line 506 of file astierCovPtcFit.py.

506  def getMaskCov(self, i, j):
507  """Get mask of Cov[i,j]"""
508  weights = self.sqrtW[:, i, j]
509  mask = weights != 0
510  return mask
511 

◆ getNoise()

def lsst.cp.pipe.astierCovPtcFit.CovFit.getNoise (   self)
Get noise matrix

Definition at line 502 of file astierCovPtcFit.py.

502  def getNoise(self):
503  """Get noise matrix"""
504  return self.params['noise'].full.reshape(self.r, self.r)
505 

◆ getNoiseCov()

def lsst.cp.pipe.astierCovPtcFit.CovFit.getNoiseCov (   self)
Get covariances of noise matrix from fit

Definition at line 470 of file astierCovPtcFit.py.

470  def getNoiseCov(self):
471  """Get covariances of noise matrix from fit"""
472  covNoise = self._getCovParams('noise')
473  return covNoise.reshape((self.r, self.r, self.r, self.r))
474 

◆ getNoiseSig()

def lsst.cp.pipe.astierCovPtcFit.CovFit.getNoiseSig (   self)
Square root of diagonal of the parameter covariance of the fitted "noise" matrix

Definition at line 475 of file astierCovPtcFit.py.

475  def getNoiseSig(self):
476  """Square root of diagonal of the parameter covariance of the fitted "noise" matrix"""
477  if self._getCovParams('noise') is not None:
478  covNoise = self._getCovParams('noise')
479  noise = np.sqrt(covNoise.diagonal()).reshape((self.r, self.r))
480  else:
481  noise = None
482  return noise
483 

◆ getParamValues()

def lsst.cp.pipe.astierCovPtcFit.CovFit.getParamValues (   self)
Return an array of free parameter values (it is a copy).

Definition at line 333 of file astierCovPtcFit.py.

333  def getParamValues(self):
334  """Return an array of free parameter values (it is a copy)."""
335  return self.params.free + 0.
336 

◆ getRon()

def lsst.cp.pipe.astierCovPtcFit.CovFit.getRon (   self)
Get readout noise (e^2)

Definition at line 488 of file astierCovPtcFit.py.

488  def getRon(self):
489  """Get readout noise (e^2)"""
490  return self.params['noise'].full[0]
491 

◆ getRonErr()

def lsst.cp.pipe.astierCovPtcFit.CovFit.getRonErr (   self)
Get error on readout noise parameter

Definition at line 492 of file astierCovPtcFit.py.

492  def getRonErr(self):
493  """Get error on readout noise parameter"""
494  ronSqrt = np.sqrt(np.fabs(self.getRon()))
495  if self.getNoiseSig() is not None:
496  noiseSigma = self.getNoiseSig()[0][0]
497  ronErr = 0.5*(noiseSigma/np.fabs(self.getRon()))*ronSqrt
498  else:
499  ronErr = np.nan
500  return ronErr
501 

◆ initFit()

def lsst.cp.pipe.astierCovPtcFit.CovFit.initFit (   self)
 Performs a crude parabolic fit of the data in order to start
the full fit close to the solution.

Definition at line 292 of file astierCovPtcFit.py.

292  def initFit(self):
293  """ Performs a crude parabolic fit of the data in order to start
294  the full fit close to the solution.
295  """
296  # number of parameters for 'a' array.
297  lenA = self.r*self.r
298  # define parameters: c corresponds to a*b in Astier+19 (Eq. 20).
299  self.params = FitParameters([('a', lenA), ('c', lenA), ('noise', lenA), ('gain', 1)])
300  self.params['gain'] = 1.
301  # c=0 in a first go.
302  self.params['c'].fix(val=0.)
303  # plumbing: extract stuff from the parameter structure
304  a = self.params['a'].full.reshape(self.r, self.r)
305  noise = self.params['noise'].full.reshape(self.r, self.r)
306  gain = self.params['gain'].full[0]
307 
308  # iterate the fit to account for higher orders
309  # the chi2 does not necessarily go down, so one could
310  # stop when it increases
311  oldChi2 = 1e30
312  for _ in range(5):
313  model = self.evalCovModel() # this computes the full model.
314  # loop on lags
315  for i in range(self.r):
316  for j in range(self.r):
317  # fit a parabola for a given lag
318  parsFit = np.polyfit(self.mu, self.cov[:, i, j] - model[:, i, j],
319  2, w=self.sqrtW[:, i, j])
320  # model equation(Eq. 20) in Astier+19:
321  a[i, j] += parsFit[0]
322  noise[i, j] += parsFit[2]*gain*gain
323  if(i + j == 0):
324  gain = 1./(1/gain+parsFit[1])
325  self.params['gain'].full[0] = gain
326  chi2 = self.chi2()
327  if chi2 > oldChi2:
328  break
329  oldChi2 = chi2
330 
331  return
332 

◆ ndof()

def lsst.cp.pipe.astierCovPtcFit.CovFit.ndof (   self)
Number of degrees of freedom

Returns
-------
mask.sum() - len(self.params.free): `int`
    Number of usable pixels - number of parameters of fit.

Definition at line 582 of file astierCovPtcFit.py.

582  def ndof(self):
583  """Number of degrees of freedom
584 
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
591 
592  return mask.sum() - len(self.params.free)
593 

◆ setAandB()

def lsst.cp.pipe.astierCovPtcFit.CovFit.setAandB (   self,
  a,
  b 
)
Set "a" and "b" coeffcients forfull Astier+19 model (Eq. 20). "c=a*b".

Definition at line 512 of file astierCovPtcFit.py.

512  def setAandB(self, a, b):
513  """Set "a" and "b" coeffcients forfull Astier+19 model (Eq. 20). "c=a*b"."""
514  self.params['a'].full = a.flatten()
515  self.params['c'].full = a.flatten()*b.flatten()
516  return
517 

◆ setParamValues()

def lsst.cp.pipe.astierCovPtcFit.CovFit.setParamValues (   self,
  p 
)
Set parameter values.

Definition at line 337 of file astierCovPtcFit.py.

337  def setParamValues(self, p):
338  """Set parameter values."""
339  self.params.free = p
340  return
341 

◆ subtractDistantOffset()

def lsst.cp.pipe.astierCovPtcFit.CovFit.subtractDistantOffset (   self,
  maxLag = 8,
  startLag = 5,
  polDegree = 1 
)
Subtract a background/offset to the measured covariances.

Parameters
---------
maxLag: `int`
    Maximum lag considered

startLag: `int`
    First lag from where to start the offset subtraction.

polDegree: `int`
    Degree of 2D polynomial to fit to covariance to define offse to be subtracted.

Definition at line 252 of file astierCovPtcFit.py.

252  def subtractDistantOffset(self, maxLag=8, startLag=5, polDegree=1):
253  """Subtract a background/offset to the measured covariances.
254 
255  Parameters
256  ---------
257  maxLag: `int`
258  Maximum lag considered
259 
260  startLag: `int`
261  First lag from where to start the offset subtraction.
262 
263  polDegree: `int`
264  Degree of 2D polynomial to fit to covariance to define offse to be subtracted.
265  """
266  assert(startLag < self.r)
267  for k in range(len(self.mu)):
268  # Make a copy because it is going to be altered
269  w = self.sqrtW[k, ...] + 0.
270  sh = w.shape
271  i, j = np.meshgrid(range(sh[0]), range(sh[1]), indexing='ij')
272  # kill the core for the fit
273  w[:startLag, :startLag] = 0
274  poly = Pol2d(i, j, self.cov[k, ...], polDegree+1, w=w)
275  back = poly.eval(i, j)
276  self.cov[k, ...] -= back
277  self.r = maxLag
278  self.cov = self.cov[:, :maxLag, :maxLag]
279  self.vcov = self.vcov[:, :maxLag, :maxLag]
280  self.sqrtW = self.sqrtW[:, :maxLag, :maxLag]
281 
282  return
283 

◆ weightedRes()

def lsst.cp.pipe.astierCovPtcFit.CovFit.weightedRes (   self,
  params = None 
)
Weighted residuas.

Notes
-----
To be used via:
c = CovFit(nt)
c.initFit()
coeffs, cov, _, mesg, ierr = leastsq(c.weightedRes, c.getParamValues(), full_output=True)

Definition at line 532 of file astierCovPtcFit.py.

532  def weightedRes(self, params=None):
533  """Weighted residuas.
534 
535  Notes
536  -----
537  To be used via:
538  c = CovFit(nt)
539  c.initFit()
540  coeffs, cov, _, mesg, ierr = leastsq(c.weightedRes, c.getParamValues(), full_output=True)
541  """
542  return self.wres(params).flatten()
543 

◆ wres()

def lsst.cp.pipe.astierCovPtcFit.CovFit.wres (   self,
  params = None 
)
To be used in weightedRes

Definition at line 522 of file astierCovPtcFit.py.

522  def wres(self, params=None):
523  """To be used in weightedRes"""
524  if params is not None:
525  self.setParamValues(params)
526  covModel = self.evalCovModel()
527  weightedRes = (covModel-self.cov)*self.sqrtW
528  maskedWeightedRes = weightedRes[self.maskMu]
529 
530  return maskedWeightedRes
531 

Member Data Documentation

◆ cov

lsst.cp.pipe.astierCovPtcFit.CovFit.cov

Definition at line 278 of file astierCovPtcFit.py.

◆ covParams

lsst.cp.pipe.astierCovPtcFit.CovFit.covParams

Definition at line 578 of file astierCovPtcFit.py.

◆ logger

lsst.cp.pipe.astierCovPtcFit.CovFit.logger

Definition at line 246 of file astierCovPtcFit.py.

◆ maskMu

lsst.cp.pipe.astierCovPtcFit.CovFit.maskMu

Definition at line 248 of file astierCovPtcFit.py.

◆ mu

lsst.cp.pipe.astierCovPtcFit.CovFit.mu

Definition at line 242 of file astierCovPtcFit.py.

◆ params

lsst.cp.pipe.astierCovPtcFit.CovFit.params

Definition at line 299 of file astierCovPtcFit.py.

◆ r

lsst.cp.pipe.astierCovPtcFit.CovFit.r

Definition at line 245 of file astierCovPtcFit.py.

◆ sqrtW

lsst.cp.pipe.astierCovPtcFit.CovFit.sqrtW

Definition at line 244 of file astierCovPtcFit.py.

◆ vcov

lsst.cp.pipe.astierCovPtcFit.CovFit.vcov

Definition at line 279 of file astierCovPtcFit.py.


The documentation for this class was generated from the following file:
lsst.cp.pipe.astierCovPtcFit.makeCovArray
def makeCovArray(inputTuple, maxRangeFromTuple=8)
Definition: astierCovPtcFit.py:67
chi2
lsst.cp.pipe.astierCovPtcFit.symmetrize
def symmetrize(inputArray)
Definition: astierCovPtcFit.py:159