LTF
Classes | Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
LTF Class Reference

Linear template fit. More...

#include <LTF.h>

Classes

class  LiTeFit
 Small helper class to collect all results, All input parameters are set by class LTF The user can access all input and output quantities. More...
 

Public Types

enum  Uncertainty { Constrained =1, Unconstrained =0, External =-1 }
 

Public Member Functions

void SetData (const vector< double > &data)
 < More...
 
void SetData (size_t n, const double *data)
 
void AddTemplate (double mval, const vector< double > &dist)
 Add new template for a given reference point. More...
 
void AddTemplate (double mval, size_t n, const double *dist)
 Add new template for a given reference point. More...
 
void AddTemplate (const vector< double > mvals, size_t n, const double *dist)
 Add new template for a given reference point.
 
void AddTemplate (const vector< double > &mvals, const vector< double > &dist)
 
void AddError (std::string name, size_t n, const double *error, double corr=0., LTF::Uncertainty constr=LTF::Uncertainty::Constrained)
 Add a new error source to LTF. More...
 
void AddError (const std::string &name, const std::vector< double > &error, double corr=0., LTF::Uncertainty constr=LTF::Uncertainty::Constrained)
 Add a new error source to LTF.
 
void AddError (const std::string &name, const std::vector< std::vector< double > > &cov, LTF::Uncertainty constr=LTF::Uncertainty::Constrained)
 Set covariance matrix as error source. More...
 
void AddErrorRelative (const std::string &name, const std::vector< std::vector< double > > &cov_rel, LTF::Uncertainty constr=LTF::Uncertainty::Constrained)
 Set covariance matrix as error source. More...
 
void AddErrorRelative (std::string name, size_t n, const double *error, double corr=0., LTF::Uncertainty constr=LTF::Uncertainty::Constrained)
 Add a new error source to LTF with relative values. More...
 
void AddErrorRelative (const std::string &name, std::vector< double > error, double corr=0., LTF::Uncertainty constr=LTF::Uncertainty::Constrained)
 Add a new error source to LTF with relative values.
 
void AddErrorPercent (std::string name, size_t n, const double *error, double corr=0., LTF::Uncertainty constr=LTF::Uncertainty::Constrained)
 Add a new error source to LTF. More...
 
void AddErrorPercent (const std::string &name, std::vector< double > error, double corr=0., LTF::Uncertainty constr=LTF::Uncertainty::Constrained)
 Add a new error source to LTF with relative values.
 
void AddUncorrelatedErrorSquared (const std::string &name, size_t n, const double *error2, LTF::Uncertainty constr=LTF::Uncertainty::Constrained)
 Set uncorrelated error. More...
 
void AddCorrelatedError (const std::string &name, vector< double > error, LTF::Uncertainty constr=LTF::Uncertainty::Constrained)
 Set correlated error. More...
 
void AddTemplateError (const std::string &name, const std::vector< double > &mvals, size_t n, const double *error, double corr=0.)
 Add a new error source of the templates to LTF Note: a repeated call with an already existent name will overwrite a previous uncertainty, or alternatively sets/adds a value if the uncertianty was unspecified or zero before. More...
 
void AddTemplateError (const std::string &name, double mval, size_t n, const double *error, double corr=0.)
 Add a new error source of the templates to LTF.
 
void AddTemplateError (const std::string &name, const std::vector< double > &mvals, const std::vector< double > &error, double corr=0.)
 Add a new error source of the templates to LTF.
 
void AddTemplateError (const std::string &name, double mval, const std::vector< double > &error, double corr=0.)
 Add a new error source of the templates to LTF.
 
void AddTemplateErrorSquared (const std::string &name, const std::vector< double > &mvals, size_t n, const double *error2, double corr=0.)
 Add a new error source of the templates to LTF.
 
void AddTemplateErrorSquared (const std::string &name, double mval, size_t n, const double *error2, double corr=0.)
 Add a new error source of the templates to LTF.
 
void AddTemplateErrorSquared (const std::string &name, const std::vector< double > &mvals, const std::vector< double > &error2, double corr=0.)
 Add a new error source of the templates to LTF.
 
void SetResponseMatrix (const std::vector< std::vector< double > > &A)
 set a (detector) response matrix (optional) A detector response matrix can represent resolution and acceptance effects of the experimental device. The response matrix A 'maps' the predictions to the detector level: y=Ax A response matrix furthermore may account for different binnings of the templates and the data, or it may represent background contributions, etc...
 
void SetResponseMatrix (const Eigen::MatrixXd &A)
 set a (detector) response matrix (optional)
 
void AddResponseMatrixError (const std::string &name, const std::vector< std::vector< double > > &errorA, double corr=0.)
 Add a new error source of the templates to LTF Note: a repeated call with an already existent name will overwrite a previous uncertainty, or alternatively sets/adds a value if the uncertianty was unspecified or zero before. More...
 
void AddResponseMatrixError (const std::string &name, const Eigen::MatrixXd &errorA, double corr=0.)
 Add a new error source of the templates to LTF.
 
void AddResponseMatrixErrorSquared (const std::string &name, std::vector< std::vector< double > > error2A, double corr=0.)
 Add a new error source of the templates to LTF.
 
void AddResponseMatrixErrorSquared (const std::string &name, const Eigen::MatrixXd &error2A, double corr=0.)
 Add a new error source of the templates to LTF.
 
void SetGamma (double gamma=1)
 set exponential factor(s) gamma More...
 
void SetGamma (const vector< double > &gamma)
 set exponential factor(s) gamma
 
void SetFitRange (int binmin, int binmax)
 set fit range, i.e. exclude bins More...
 
const LiTeFitDoLiTeFit ()
 Perform linear template fit (LiTeFit,LTF) More...
 
const LiTeFitDoIterativeFitNewton (int nIter=3, int nPol=2, int nInfrc=1)
 Run the template fit using the iterative newton algorithm and a 2nd order model approximation.
 
void UseLogNormalUncertainties (bool LogNormal=true)
 Select normal distributed or log-normal distributed uncertainties note: log-normal distributed uncertainties are equivalent to normal-distributed relative uncertainties and can be considered as a good approximation to a poisson distribution.
 
bool GetLogNormalUncertainties () const
 get boolean for relative uncertainties (LogNormal)
 
void UseNuisanceParameters (bool UseNuisance=true)
 Consider correlated systematic uncertainties as 'shifts' with nuisance. More...
 
void RescaleInputErrors (std::vector< double > values)
 Rescale uncertainties All uncertainties are stored interally as absolute uncertainties. We assume that these are proportional to the data. In oreder to avoid (or reduce) a bias, one can re-scale all uncertainties, e.g. with the predictions (from a previous fit, or from a selected template.) The same function is availalbe for LTF::LiTeFit and a previous LiTeFit can be repeated with rescaled uncertainties. More...
 

Static Public Member Functions

static Eigen::VectorXd GetDelta_V (const Eigen::MatrixXd &V)
 calculate uncertainty from a covariance matrix
 
static Eigen::MatrixXd GetV_Delta (const Eigen::VectorXd &d, double corr)
 calculate covariance matrix from an uncertainty
 
template<class T >
static std::map< std::string, Eigen::MatrixXd > GetV_DeltaSys (const T &DeltaSys)
 calculate covariance matrices from (systematic) shifts
 
template<class T >
static std::map< std::string, Eigen::VectorXd > GetDelta_Vsource (const T &Vsource)
 calculate size of uncertainty from covariance matrix
 
static Eigen::MatrixXd VSum (const std::map< std::string, Eigen::MatrixXd > &Vs, const std::map< std::string, Eigen::VectorXd > &DeltaSys)
 Calculate covariance matrix as sum of individual matrices and systematic shifts. More...
 
static Eigen::MatrixXd VSum (const std::vector< std::pair< std::string, Eigen::MatrixXd > > &Vs)
 Calculate covariance matrix as sum of individual matrices and systematic shifts. More...
 
static Eigen::MatrixXd Cov_to_Cor (const Eigen::MatrixXd &V)
 Calculate correlation matrix from given covariacne matrix V.
 
static Eigen::MatrixXd Std_to_Mat (const std::vector< std::vector< double > > &A)
 Calculate an Eigen-matrix from std::vector<std::vector<>>. T=double.
 

Protected Member Functions

size_t N () const
 return number of points. This value is not valid after applying apply Fit Range
 
Eigen::MatrixXd Calc_M () const
 Calculate matrix M from template reference points.
 
Eigen::MatrixXd Calc_Y () const
 Calculate matrix Y, representing the templates.
 
Eigen::MatrixXd Calc_X () const
 Calculate matrix X, representing the templates.
 
std::map< std::string, Eigen::MatrixXd > Calc_dY (std::map< std::string, std::map< vector< double >, Eigen::VectorXd > > &fSysY) const
 Calculate template uncertainties dY.
 
void SetLiTeFitInput ()
 Set the input to LiTeFit. More...
 
void SetCorrSys (const std::string &name, double c)
 set correlation coefficient for syst. uncertainties
 

Protected Attributes

Eigen::VectorXd fDt
 data distribution
 
std::vector< std::pair< vector< double >, Eigen::VectorXd > > fTmplDistN
 templates
 
Eigen::MatrixXd fA
 the response matrix (opional)
 
std::vector< std::pair< std::string, Eigen::MatrixXd > > fVs
 Covariance matrices.
 
std::vector< std::pair< std::string, Eigen::VectorXd > > fSys
 correlated (syst.) uncertainty
 
std::map< std::string, LTF::Uncertainty > fSysIsConstr
 uncertainty specification for fSys
 
std::map< std::string, Eigen::MatrixXd > fVsExt
 Covariance matrices, treated as external.
 
std::map< std::string, Eigen::VectorXd > fSysExt
 correlated uncertainties, treated as external
 
std::map< std::string, std::map< vector< double >, Eigen::VectorXd > > fSysY
 uncertainty of the templates (uncorrelated)
 
std::map< std::string, Eigen::MatrixXd > fSysA
 uncertainty of the response matrix (uncorrelated)
 
std::map< std::string, double > fcorrSys
 Correlation coefficient of the systematic uncertainties (dA,dY)
 
Eigen::VectorXd fErrorScale
 Reference values for error rescaling (default: data)
 
std::vector< double > fGamma
 exponential parameter (a+b*x^\gamma)
 
bool fUseLog = false
 Use relative uncertainties.
 
bool fUseNuisance = true
 Nuisance parameter or a covariance matrices.
 
LiTeFit fLTF
 result
 
int fBinMin = 0
 Fit range.
 
int fBinNN = 0
 Fit range.
 

Detailed Description

Linear template fit.

implemented with the linear algebra package eigen

Member Function Documentation

◆ AddCorrelatedError()

void LTF::AddCorrelatedError ( const std::string &  name,
vector< double >  error,
LTF::Uncertainty  constr = LTF::Uncertainty::Constrained 
)

Set correlated error.

Parameters
nameName of the error source
errorsize of the error (or systematic shift), taking care for the sign of the shifts
constrtreat uncertainty constrained or unconstrained in fit, or only propagate it after the fit ('External')

Note: at least one uncorrelated (or covariance) error must be specified. Function used internally for initializing member variables

◆ AddError() [1/2]

void LTF::AddError ( std::string  name,
size_t  n,
const double *  error,
double  corr = 0.,
LTF::Uncertainty  constr = LTF::Uncertainty::Constrained 
)

Add a new error source to LTF.

Parameters
nameName of the error source
nnumber of bins
errorsize of the error
corrsize of correlation (0=uncorrelated, 1=correlated, (0-1) split-up)
constrtreat uncertainty constrained or unconstrained in fit, or only propagate it after the fit ('External')

◆ AddError() [2/2]

void LTF::AddError ( const std::string &  name,
const std::vector< std::vector< double > > &  cov,
LTF::Uncertainty  constr = LTF::Uncertainty::Constrained 
)

Set covariance matrix as error source.

Parameters
nameName of the error source
covcovariance matrix (full-matrix notation, symmetric)
constrtreat uncertainty constrained or unconstrained in fit, or only propagate it after the fit ('External')

note: covariance matrix is always 'constrained' in fit. note: Function used interally for initializing member variables

◆ AddErrorPercent()

void LTF::AddErrorPercent ( std::string  name,
size_t  n,
const double *  error,
double  corr = 0.,
LTF::Uncertainty  constr = LTF::Uncertainty::Constrained 
)
inline

Add a new error source to LTF.

Parameters
errorsize of the relative error in percent (%)

◆ AddErrorRelative() [1/2]

void LTF::AddErrorRelative ( const std::string &  name,
const std::vector< std::vector< double > > &  cov_rel,
LTF::Uncertainty  constr = LTF::Uncertainty::Constrained 
)

Set covariance matrix as error source.

Parameters
nameName of the error source
covcovariance matrix (full-matrix notation, symmetric, with relative uncertainties)
constrtreat uncertainty constrained or unconstrained in fit, or only propagate it after the fit ('External')

note: covariance matrix is always 'constrained' in fit. note: internally, the relative uncertainties are first multiplied by data, and later on again divided by data.

◆ AddErrorRelative() [2/2]

void LTF::AddErrorRelative ( std::string  name,
size_t  n,
const double *  error,
double  corr = 0.,
LTF::Uncertainty  constr = LTF::Uncertainty::Constrained 
)
inline

Add a new error source to LTF with relative values.

Parameters
errorsize of the relative error Pass a relative uncertainty to LTF Note: This setter function only passes the uncertainty values as relative uncertainties to LTF, whereas the probability density function in the linear template fit can be chosen with LTF::UseLogNormalUncertainties()

◆ AddResponseMatrixError()

void LTF::AddResponseMatrixError ( const std::string &  name,
const std::vector< std::vector< double > > &  errorA,
double  corr = 0. 
)
inline

Add a new error source of the templates to LTF Note: a repeated call with an already existent name will overwrite a previous uncertainty, or alternatively sets/adds a value if the uncertianty was unspecified or zero before.

Parameters
nameName of the error source
errorsize of the error of the elements of A
corrsize of correlation (0=uncorrelated, 1=correlated, (0-1) split-up [not implemented])

◆ AddTemplate() [1/2]

void LTF::AddTemplate ( double  mval,
const vector< double > &  dist 
)
inline

Add new template for a given reference point.

Parameters
mvalreference point(s) of the template distribution
nnumber of bins
disttemplate distribution

◆ AddTemplate() [2/2]

void LTF::AddTemplate ( double  mval,
size_t  n,
const double *  dist 
)
inline

Add new template for a given reference point.

Add new template for a given reference point

Add a new template distribution for a given reference values mval

◆ AddTemplateError()

void LTF::AddTemplateError ( const std::string &  name,
const std::vector< double > &  mvals,
size_t  n,
const double *  error,
double  corr = 0. 
)

Add a new error source of the templates to LTF Note: a repeated call with an already existent name will overwrite a previous uncertainty, or alternatively sets/adds a value if the uncertianty was unspecified or zero before.

Parameters
nameName of the error source
mvalsreference point(s) of the template distribution
nnumber of bins
errorsize of the error
corrsize of correlation (0=uncorrelated, 1=correlated, (0-1) split-up [not implemented])

◆ AddUncorrelatedErrorSquared()

void LTF::AddUncorrelatedErrorSquared ( const std::string &  name,
size_t  n,
const double *  error2,
LTF::Uncertainty  constr = LTF::Uncertainty::Constrained 
)
inline

Set uncorrelated error.

Parameters
nameName of the error source
nnumber of bins
error2errors squared
constrtreat uncertainty constrained or unconstrained in fit, or only propagate it after the fit ('External')

Note: at least one uncorrelated (or covariance) error must be specified.

◆ DoLiTeFit()

const LTF::LiTeFit & LTF::DoLiTeFit ( )

Perform linear template fit (LiTeFit,LTF)

Run the linear template fit and return the LiTeFit object with all results.

Returns
LiTeFit-object with all relevant quantities

◆ RescaleInputErrors()

void LTF::RescaleInputErrors ( std::vector< double >  values)
inline

Rescale uncertainties All uncertainties are stored interally as absolute uncertainties. We assume that these are proportional to the data. In oreder to avoid (or reduce) a bias, one can re-scale all uncertainties, e.g. with the predictions (from a previous fit, or from a selected template.) The same function is availalbe for LTF::LiTeFit and a previous LiTeFit can be repeated with rescaled uncertainties.

Note: this option cannot be applied when using relative uncertainties in the chi2.

Note: all errors passed to LTF are expected to scale with data, or as relative uncertainties. This method then 're-scales' the size of the absolute uncertainties.

Note: template and response-matrix uncertainties are omitted. Note: also statistical and uncorrelated uncertainties are rescaled!

◆ SetData()

void LTF::SetData ( const vector< double > &  data)
inline

<

Set data distribution Set data distribution

◆ SetFitRange()

void LTF::SetFitRange ( int  binmin,
int  binmax 
)
inline

set fit range, i.e. exclude bins

Parameters
binminfirst bin to be included in the fit
binmaxlast bin that is included in the fit

◆ SetGamma()

void LTF::SetGamma ( double  gamma = 1)
inline

set exponential factor(s) gamma

Parameters
gamma
binmaxlast bin that is included in the fit Note: better use: SetGamma({gamma1,gamma2,...,});

◆ SetLiTeFitInput()

void LTF::SetLiTeFitInput ( )
protected

Set the input to LiTeFit.

<Apply fit range, according to input parameters specified by SetFitRange() Initialize 'input' parameters of output object fLTF.

◆ UseNuisanceParameters()

void LTF::UseNuisanceParameters ( bool  UseNuisance = true)
inline

Consider correlated systematic uncertainties as 'shifts' with nuisance.

parameter or as a covariance matrix. default: true (use nuisance parameters). Note: this setting does not apply to already present uncertainties, but only to correlated uncertainties set in the following.

◆ VSum() [1/2]

Eigen::MatrixXd LTF::VSum ( const std::map< std::string, Eigen::MatrixXd > &  Vs,
const std::map< std::string, Eigen::VectorXd > &  DeltaSys 
)
static

Calculate covariance matrix as sum of individual matrices and systematic shifts.

Parameters
VsCovariance matrices, i.e. uncorr. and/or cov. uncertainties
DeltaSyscorrelated systematic uncertainty

◆ VSum() [2/2]

Eigen::MatrixXd LTF::VSum ( const std::vector< std::pair< std::string, Eigen::MatrixXd > > &  Vs)
static

Calculate covariance matrix as sum of individual matrices and systematic shifts.

Parameters
VsCovariance matrices, i.e. uncorr. and/or cov. uncertainties

The documentation for this class was generated from the following files: