![]() |
AliRoot Core
3dc7879 (3dc7879)
|
The AliTMinuitToolkit serves as an easy to use interface for the TMinuit. More...
#include <AliTMinuitToolkit.h>
Public Types | |
enum | EVerboseFlags { kPrintMinuit =0x00001, kSysInfo =0x00002, kPrintAll =0x00004, kStreamFcn =0x00010, kStreamFcnPoint =0x00020, kFCNInfo =0x000400 } |
enum | EFitFlags { kCovarFailed =0x0001 } |
Public Member Functions | |
AliTMinuitToolkit (const char *streamerName=nullptr, const char *mode="recreate") | |
virtual | ~AliTMinuitToolkit () |
void | SetStreamer (const char *streamerName, const char *mode="recreate") |
Set stream for writing info into TTree. see examples(STEER/STEERBase/TTreeStream.cxx) More... | |
TTreeSRedirector * | GetStreamer () |
void | ClearData () |
void | Fit (Option_t *option="") |
Internal function that calls the fitter. More... | |
void | FitHistogram (const TH1 *his, Option_t *option="") |
Fit a one dimensional histogram. More... | |
void | FitGraph (const TGraph *gr, Option_t *option="") |
Fit a TGraph object. More... | |
Long64_t | FillFitter (TTree *inputTree, TString values, TString variables, TString selection, Int_t firstEntry, Int_t nEntries, Bool_t doReset=kTRUE) |
Fill input matrices of the fitter. More... | |
TString | GetFitFunctionAsAlias (Option_t *option="", TTree *tree=nullptr) |
Format alias string for the for tree alias or for fit title. More... | |
void | Bootstrap (ULong_t nIter, const char *reportName, Option_t *option=nullptr) |
Bootstrap minimization. More... | |
void | TwoFoldCrossValidation (UInt_t nIter, const char *reportName, Option_t *option=nullptr) |
void | MISAC (Int_t nFitPoints, UInt_t nIter, const char *reportName, Option_t *option=nullptr) |
Fitting with outlier removal inspired by RANSAC - see https://en.wikipedia.org/w/index.php?title=Random_sample_consensus&oldid=767657742. More... | |
void | SetFitFunction (TF1 *formula, Bool_t doReset) |
Set formula of the fitted function and reset parameter. See doReset. More... | |
void | SetInitialParam (const TMatrixD *initialParam) |
Set initial parameters matrix. More... | |
TMatrixD * | GetInitialParam () |
void | SetLogLikelihoodFunction (TF1 *logLikelihood) |
void | SetFitAlgorithm (Char_t *const name) |
void | SetMaxCalls (Int_t calls) |
void | SetTolerance (Double_t tol) |
void | SetPoints (TMatrixD *points) |
void | SetValues (TMatrixD *values) |
void | EnableRobust (Bool_t b) |
const TMatrixD * | GetPoints () const |
const TMatrixD * | GetValues () const |
TF1 * | GetFormula () const |
const TF1 * | GetLogLikelihoodFunction () const |
const TVectorD * | GetParameters () const |
const TVectorD * | GetRMSEstimator () const |
const TMatrixD * | GetCovarianceMatrix () const |
const TMatrixD * | GetMISACEstimators () const |
Double_t | GetChisquare () |
Int_t | GetStatus () |
Bool_t | IsHuberCost () const |
TObjArray * | GetListOfVariables () |
void | SetVerbose (Int_t verbose) |
Bool_t | IsVerbose (Int_t flag) |
Static Public Member Functions | |
static void | FitterFCN (int &nParam, double *info, double &chi2, double *gin, int iflag) |
Internal function used in the TMinuit minimization. More... | |
static AliTMinuitToolkit * | Fit (TH1 *his, const char *fitterName, Option_t *fitOption="", Option_t *hisOption="", Option_t *funOption="", Double_t xMin=0, Double_t xMax=0) |
Make histogram fit. More... | |
static AliTMinuitToolkit * | Fit (TGraph *graph, const char *fitterName, Option_t *fitOption="", Option_t *graphOption="", Option_t *funOption="", Double_t xMin=0, Double_t xMax=0) |
Make graph fit. More... | |
static AliTMinuitToolkit * | GetPredefinedFitter (const std::string &fitterName) |
static void | SetPredefinedFitter (const std::string &fitterName, AliTMinuitToolkit *fitter) |
static void | RegisterDefaultFitters () |
static AliTMinuitToolkit * | RegisterPlaneFitter (Int_t nPlanes, Int_t logLikeType) |
Register default fitters. More... | |
static Double_t | GaussCachyLogLike (const Double_t *x, const Double_t *p) |
static Double_t | HuberLogLike (const Double_t *x, const Double_t *p) |
static Double_t | PseudoHuberLogLike (const Double_t *x, const Double_t *p) |
static Double_t | GausKurtosisSkewness (const Double_t *x, const Double_t *p) |
Gaus convoluted with Erf to emulate kurtosis and skewness. More... | |
static Double_t | RndmGaus (Double_t mean=0, Double_t sigma=1) |
Random gaus generator as static function (to use in root TFormula, TTreeFormula) More... | |
static Double_t | RndmLandau (Double_t mean=0, Double_t sigma=1) |
Random gaus generator as static function (to use in root TFormula, TTreeFormula) More... | |
static void | SetFunctionDrawOption (TF1 *fun, Option_t *option) |
Private Member Functions | |
AliTMinuitToolkit (const AliTMinuitToolkit &) | |
AliTMinuitToolkit & | operator= (const AliTMinuitToolkit &) |
Private Attributes | |
TTreeSRedirector * | fStreamer |
Int_t | fVerbose |
TF1 * | fFormula |
TF1 * | fLogLikelihoodFunction |
TString | fFitAlgorithm |
TMatrixD * | fPoints |
TMatrixD * | fValues |
TObjArray * | fVarNames |
TObjArray * | fValueNames |
TArrayI | fPointIndex |
TMatrixD * | fInitialParam |
TVectorD * | fParam |
TVectorD * | fRMSEstimator |
TMatrixD * | fMISACEstimators |
TMatrixD * | fCovar |
Int_t | fStatus |
Double_t | fChi2 |
Int_t | fMaxCalls |
Double_t | fPrecision |
Bool_t | fIsHuberCost |
Int_t | fFCNCounter |
Static Private Attributes | |
static std::map< std::string, AliTMinuitToolkit * > | fPredefinedFitters |
The AliTMinuitToolkit serves as an easy to use interface for the TMinuit.
Even a few outliers can lead to wrong results of a least-squares fitting procedure. In this case the use of robust(resistant) methods can be helpful, but a stronger dependence on starting values or convergence to local minimum can occur.
Example usage in the $AliRoot_SRC/STAT/test/
Definition at line 104 of file AliTMinuitToolkit.h.
Enumerator | |
---|---|
kCovarFailed |
Definition at line 115 of file AliTMinuitToolkit.h.
Enumerator | |
---|---|
kPrintMinuit | |
kSysInfo | |
kPrintAll | |
kStreamFcn | |
kStreamFcnPoint | |
kFCNInfo |
Definition at line 106 of file AliTMinuitToolkit.h.
AliTMinuitToolkit::AliTMinuitToolkit | ( | const char * | streamerName = nullptr , |
const char * | mode = "recreate" |
||
) |
Definition at line 43 of file AliTMinuitToolkit.cxx.
Referenced by MISAC(), RegisterPlaneFitter(), and SetPredefinedFitter().
|
virtual |
Definition at line 78 of file AliTMinuitToolkit.cxx.
|
private |
void AliTMinuitToolkit::Bootstrap | ( | ULong_t | nIter, |
const char * | reportName, | ||
Option_t * | option = nullptr |
||
) |
Bootstrap minimization.
nIter | - number of iteration |
reportName | - report name - in case results streamed |
option | - fit option |
Definition at line 582 of file AliTMinuitToolkit.cxx.
Referenced by Fit(), GetStreamer(), and Test1D().
void AliTMinuitToolkit::ClearData | ( | ) |
Definition at line 107 of file AliTMinuitToolkit.cxx.
Referenced by FillFitter(), FitGraph(), FitHistogram(), and GetStreamer().
|
inline |
Definition at line 144 of file AliTMinuitToolkit.h.
Referenced by AliTMinuitToolkit_TestHistogram(), and Test1D().
Long64_t AliTMinuitToolkit::FillFitter | ( | TTree * | inputTree, |
TString | values, | ||
TString | variables, | ||
TString | selection, | ||
Int_t | firstEntry, | ||
Int_t | nEntries, | ||
Bool_t | doReset = kTRUE |
||
) |
Fill input matrices of the fitter.
inputTree | - pointer to input tree |
values | - values description (parameter:<error>) |
variables | - variables array (var0:<var1>) |
selection | - tree selection |
firstEntry | - first entry to query |
nEntries | - number of entries to query |
doReset | - switch reset previous/or append |
Definition at line 415 of file AliTMinuitToolkit.cxx.
Referenced by GetStreamer(), and Test1D().
void AliTMinuitToolkit::Fit | ( | Option_t * | option = "" | ) |
Internal function that calls the fitter.
option |
|
Definition at line 275 of file AliTMinuitToolkit.cxx.
Referenced by AliTMinuitToolkit_TestHistogram(), AliTMinuitToolkitTestLinear(), Bootstrap(), FitGraph(), FitHistogram(), GetStreamer(), IsVerbose(), MISAC(), Test1D(), and TwoFoldCrossValidation().
|
static |
Make histogram fit.
his | - pointer to histogram |
fitterName | - fitter name () |
fitOption | - fit Option |
hisOption | - graphic option for histogram |
funOption | - graphic option for function drawing |
xMin | - fit range min |
xMax | - fit range max |
Definition at line 943 of file AliTMinuitToolkit.cxx.
|
static |
Make graph fit.
graph | - pointer to graph |
fitterName | - fitter name () |
fitOption | - fitter setting option |
graphOption | - graphic option |
funOption | - graphic option for function |
xMin | - fit range min |
xMax | - fit range max |
Definition at line 979 of file AliTMinuitToolkit.cxx.
void AliTMinuitToolkit::FitGraph | ( | const TGraph * | gr, |
Option_t * | option = "" |
||
) |
Fit a TGraph object.
gr | - const TGraph object |
option | - fit option e.g (default, bootstrap20,misac .... )
|
Definition at line 169 of file AliTMinuitToolkit.cxx.
Referenced by Fit(), and GetStreamer().
void AliTMinuitToolkit::FitHistogram | ( | const TH1 * | his, |
Option_t * | option = "" |
||
) |
Fit a one dimensional histogram.
his | - const TH1 object |
option | - fit option e.g (default, bootstrap20,misac ... )
|
Definition at line 150 of file AliTMinuitToolkit.cxx.
Referenced by AliTMinuitToolkit_TestHistogram(), Fit(), and GetStreamer().
|
static |
Internal function used in the TMinuit minimization.
nPar | - not used in our interface |
grad | - gradient of cost function |
chi2(fVal)) | - cost function to minimize (chi2/log-likelihood/) |
gin(par) | - function parameters |
iflag(flag) | - flag used by root (1.) print/ 2.) derivative/ 3. final fit)
|
Definition at line 205 of file AliTMinuitToolkit.cxx.
Referenced by Bootstrap(), Fit(), MISAC(), SetValues(), and TwoFoldCrossValidation().
|
static |
Gaus convoluted with Erf to emulate kurtosis and skewness.
x | |
p |
|
Definition at line 562 of file AliTMinuitToolkit.cxx.
Referenced by SetPredefinedFitter().
|
static |
Log likelihood for the gaus+cauchy - used for Log likelihood minimization of distribution with tails
x | - |
p | - p[0] - gaus fraction
|
Definition at line 524 of file AliTMinuitToolkit.cxx.
Referenced by AliTMinuitToolkit_TestHistogram(), MISAC(), RegisterGaussFitters(), RegisterPlaneFitter(), and SetPredefinedFitter().
|
inline |
Definition at line 155 of file AliTMinuitToolkit.h.
|
inline |
Definition at line 152 of file AliTMinuitToolkit.h.
Referenced by AliTMinuitToolkitTestLinear().
TString AliTMinuitToolkit::GetFitFunctionAsAlias | ( | Option_t * | option = "" , |
TTree * | tree = nullptr |
||
) |
Format alias string for the for tree alias or for fit title.
option | - use latex in case of title request |
Definition at line 465 of file AliTMinuitToolkit.cxx.
Referenced by GetStreamer(), and Test1D().
|
inline |
Definition at line 148 of file AliTMinuitToolkit.h.
Referenced by AliTMinuitToolkit_TestHistogram(), Fit(), and FitterFCN().
|
inline |
Definition at line 136 of file AliTMinuitToolkit.h.
Referenced by AliTMinuitToolkitTestLinear().
|
inline |
Definition at line 158 of file AliTMinuitToolkit.h.
Referenced by GetFitFunctionAsAlias().
|
inline |
Definition at line 149 of file AliTMinuitToolkit.h.
Referenced by AliTMinuitToolkitTestLinear().
|
inline |
Definition at line 153 of file AliTMinuitToolkit.h.
Referenced by AliTMinuitToolkitTestLinear().
|
inline |
Definition at line 150 of file AliTMinuitToolkit.h.
Referenced by AliTMinuitToolkit_TestHistogram(), AliTMinuitToolkitTestLinear(), Fit(), and Test1D().
|
inline |
Definition at line 146 of file AliTMinuitToolkit.h.
Referenced by FitterFCN().
|
inlinestatic |
Definition at line 164 of file AliTMinuitToolkit.h.
Referenced by AliTMinuitToolkitTestLinear(), and Fit().
|
inline |
Definition at line 151 of file AliTMinuitToolkit.h.
Referenced by AliTMinuitToolkitTestLinear(), and GetFitFunctionAsAlias().
|
inline |
Definition at line 156 of file AliTMinuitToolkit.h.
|
inline |
Definition at line 122 of file AliTMinuitToolkit.h.
|
inline |
Definition at line 147 of file AliTMinuitToolkit.h.
Referenced by FitterFCN().
|
static |
Huber loss is a loss function used in robust regression - representation as in https://en.wikipedia.org/w/index.php?title=Huber_loss&diff=809252384&oldid=798150659
x | - |
p | - p[0] - gaus region boundary |
Definition at line 536 of file AliTMinuitToolkit.cxx.
Referenced by SetPredefinedFitter().
|
inline |
Definition at line 157 of file AliTMinuitToolkit.h.
Referenced by FitterFCN().
|
inline |
Definition at line 160 of file AliTMinuitToolkit.h.
Referenced by FitterFCN().
void AliTMinuitToolkit::MISAC | ( | Int_t | nFitPoints, |
UInt_t | nIter, | ||
const char * | reportName, | ||
Option_t * | option = nullptr |
||
) |
Fitting with outlier removal inspired by RANSAC - see https://en.wikipedia.org/w/index.php?title=Random_sample_consensus&oldid=767657742.
nFitPoints | - number of points used in one fit |
nIter | - number of iterations |
reportName | - report name (!=NULL- specified stream to save intermediate data) |
option | - fit option |
TODO- check failure of the fit
Definition at line 720 of file AliTMinuitToolkit.cxx.
Referenced by Fit(), and GetStreamer().
|
private |
Referenced by SetPredefinedFitter().
|
static |
Pseudo Huber loss is a loss function used in robust regression - representation as in https://en.wikipedia.org/w/index.php?title=Huber_loss&diff=809252384&oldid=798150659
x | - |
p | - p[0] - gaus region boundary |
Definition at line 547 of file AliTMinuitToolkit.cxx.
Referenced by AliTMinuitToolkit_TestHistogram(), MISAC(), RegisterPlaneFitter(), and SetPredefinedFitter().
|
static |
Referenced by AliTMinuitToolkitTestLinear(), MISAC(), and SetPredefinedFitter().
|
static |
Register default fitters.
nPlanes | - number of planes |
logLikeType | - type of the log likelihood
|
Definition at line 902 of file AliTMinuitToolkit.cxx.
Referenced by SetPredefinedFitter().
|
static |
Random gaus generator as static function (to use in root TFormula, TTreeFormula)
mean | - mean Gaus |
sigma | - sigma gaus |
Definition at line 505 of file AliTMinuitToolkit.cxx.
Referenced by SetPredefinedFitter().
|
static |
Random gaus generator as static function (to use in root TFormula, TTreeFormula)
mean | |
sigma |
Definition at line 513 of file AliTMinuitToolkit.cxx.
Referenced by SetPredefinedFitter().
|
inline |
Definition at line 138 of file AliTMinuitToolkit.h.
void AliTMinuitToolkit::SetFitFunction | ( | TF1 * | formula, |
Bool_t | doReset | ||
) |
Set formula of the fitted function and reset parameter. See doReset.
formula | - formula of the fitted function |
doReset | - if true reset parameters related to the fit function |
Definition at line 120 of file AliTMinuitToolkit.cxx.
Referenced by AliTMinuitToolkit_TestHistogram(), GetStreamer(), MISAC(), RegisterGaussFitters(), RegisterPlaneFitter(), and Test1D().
|
static |
Parse draw option for function
fun | - pointer to function to change |
option | - draw option for function
|
Definition at line 1014 of file AliTMinuitToolkit.cxx.
Referenced by Fit(), and SetPredefinedFitter().
void AliTMinuitToolkit::SetInitialParam | ( | const TMatrixD * | paramLimits | ) |
Set initial parameters matrix.
paramLimits | - Initialization matrix
|
Definition at line 140 of file AliTMinuitToolkit.cxx.
Referenced by AliTMinuitToolkit_TestHistogram(), GetStreamer(), MISAC(), RegisterGaussFitters(), RegisterPlaneFitter(), and Test1D().
|
inline |
Definition at line 137 of file AliTMinuitToolkit.h.
Referenced by AliTMinuitToolkit_TestHistogram(), MISAC(), RegisterGaussFitters(), RegisterPlaneFitter(), and Test1D().
|
inline |
Definition at line 139 of file AliTMinuitToolkit.h.
|
inline |
Definition at line 141 of file AliTMinuitToolkit.h.
|
inlinestatic |
Definition at line 165 of file AliTMinuitToolkit.h.
Referenced by AliTMinuitToolkit_TestHistogram(), MISAC(), RegisterGaussFitters(), and RegisterPlaneFitter().
void AliTMinuitToolkit::SetStreamer | ( | const char * | streamerName, |
const char * | mode = "recreate" |
||
) |
Set stream for writing info into TTree. see examples(STEER/STEERBase/TTreeStream.cxx)
streamerName | - Name of .root file |
mode | - option for .root file open
|
Definition at line 99 of file AliTMinuitToolkit.cxx.
Referenced by AliTMinuitToolkitTestLinear().
|
inline |
Definition at line 140 of file AliTMinuitToolkit.h.
|
inline |
Definition at line 142 of file AliTMinuitToolkit.h.
|
inline |
Definition at line 159 of file AliTMinuitToolkit.h.
Referenced by MISAC(), RegisterGaussFitters(), and Test1D().
void AliTMinuitToolkit::TwoFoldCrossValidation | ( | UInt_t | nIter, |
const char * | reportName, | ||
Option_t * | option = nullptr |
||
) |
Two fold cross validation - https://en.wikipedia.org/wiki/Cross-validation_(statistics)
nIter | - number of iterations |
reportName | - report name (identifier in the streamer) |
option | - fit option |
Definition at line 637 of file AliTMinuitToolkit.cxx.
Referenced by Fit(), GetStreamer(), and Test1D().
|
private |
Definition at line 197 of file AliTMinuitToolkit.h.
Referenced by Bootstrap(), Fit(), GetChisquare(), and MISAC().
|
private |
Definition at line 195 of file AliTMinuitToolkit.h.
Referenced by Bootstrap(), Fit(), GetCovarianceMatrix(), MISAC(), SetFitFunction(), TwoFoldCrossValidation(), and ~AliTMinuitToolkit().
|
private |
Definition at line 201 of file AliTMinuitToolkit.h.
Referenced by Bootstrap(), FitterFCN(), and TwoFoldCrossValidation().
|
private |
Definition at line 185 of file AliTMinuitToolkit.h.
Referenced by Fit(), and SetFitAlgorithm().
|
private |
Definition at line 183 of file AliTMinuitToolkit.h.
Referenced by Bootstrap(), Fit(), FitterFCN(), GetFitFunctionAsAlias(), GetFormula(), MISAC(), SetFitFunction(), and TwoFoldCrossValidation().
|
private |
Definition at line 191 of file AliTMinuitToolkit.h.
Referenced by Fit(), GetInitialParam(), SetInitialParam(), and ~AliTMinuitToolkit().
|
private |
Definition at line 200 of file AliTMinuitToolkit.h.
Referenced by EnableRobust(), and IsHuberCost().
|
private |
Definition at line 184 of file AliTMinuitToolkit.h.
Referenced by FitterFCN(), GetLogLikelihoodFunction(), and SetLogLikelihoodFunction().
|
private |
Definition at line 198 of file AliTMinuitToolkit.h.
Referenced by Fit(), and SetMaxCalls().
|
private |
Definition at line 194 of file AliTMinuitToolkit.h.
Referenced by GetMISACEstimators(), MISAC(), and ~AliTMinuitToolkit().
|
private |
Definition at line 192 of file AliTMinuitToolkit.h.
Referenced by Bootstrap(), Fit(), GetFitFunctionAsAlias(), GetParameters(), MISAC(), SetFitFunction(), TwoFoldCrossValidation(), and ~AliTMinuitToolkit().
|
private |
Definition at line 190 of file AliTMinuitToolkit.h.
Referenced by Bootstrap(), FitterFCN(), MISAC(), and TwoFoldCrossValidation().
|
private |
Definition at line 186 of file AliTMinuitToolkit.h.
Referenced by Bootstrap(), ClearData(), FillFitter(), FitGraph(), FitHistogram(), GetPoints(), MISAC(), SetPoints(), TwoFoldCrossValidation(), and ~AliTMinuitToolkit().
|
private |
Definition at line 199 of file AliTMinuitToolkit.h.
Referenced by Fit(), and SetTolerance().
|
staticprivate |
Definition at line 202 of file AliTMinuitToolkit.h.
Referenced by GetPredefinedFitter(), and SetPredefinedFitter().
|
private |
Definition at line 193 of file AliTMinuitToolkit.h.
Referenced by Bootstrap(), GetRMSEstimator(), MISAC(), SetFitFunction(), and ~AliTMinuitToolkit().
|
private |
Definition at line 196 of file AliTMinuitToolkit.h.
Referenced by Bootstrap(), Fit(), GetStatus(), MISAC(), and TwoFoldCrossValidation().
|
private |
Definition at line 181 of file AliTMinuitToolkit.h.
Referenced by Bootstrap(), FitterFCN(), GetStreamer(), MISAC(), SetStreamer(), TwoFoldCrossValidation(), and ~AliTMinuitToolkit().
|
private |
Definition at line 189 of file AliTMinuitToolkit.h.
Referenced by FillFitter().
|
private |
Definition at line 187 of file AliTMinuitToolkit.h.
Referenced by ClearData(), FillFitter(), Fit(), FitGraph(), FitHistogram(), GetValues(), SetValues(), and ~AliTMinuitToolkit().
|
private |
Definition at line 188 of file AliTMinuitToolkit.h.
Referenced by FillFitter(), GetFitFunctionAsAlias(), and GetListOfVariables().
|
private |
Definition at line 182 of file AliTMinuitToolkit.h.
Referenced by Fit(), IsVerbose(), and SetVerbose().