AliRoot Core  3dc7879 (3dc7879)
AliTMinuitToolkit.h
Go to the documentation of this file.
1 #ifndef ALITMINUITTOOLKIT_H
2 #define ALITMINUITTOOLKIT_H
3 
85 
86 //--------------------------------------------------------------------------------------
87 
88 
89 #include <TNamed.h>
90 class TH1;
91 class TH1F;
92 class TH2F;
93 class TFormula;
94 class TF1;
95 class TTree;
96 class TGraph;
97 class TTreeSRedirector;
98 #include <TVectorDfwd.h>
99 #include <TMatrixDfwd.h>
100 #include <TArrayI.h>
101 #include <map>
102 
103 
104 class AliTMinuitToolkit: public TNamed{
105 public:
106  enum EVerboseFlags{ // flags for verbosity
107  kPrintMinuit=0x00001, // flag:print Minuit messages
108  kSysInfo=0x00002, // flag:print system Information for profiling (Memory, CPU...)
109  kPrintAll=0x00004, // flag:print detailed information in minuit
110  //
111  kStreamFcn=0x00010, // flag: stream fcn values
112  kStreamFcnPoint=0x00020, // flag: stream fcn points
113  kFCNInfo=0x000400 // flag print FCN information
114  };
115  enum EFitFlags{
116  kCovarFailed=0x0001 // flag - problem to extract covariance matrix
117  };
118  AliTMinuitToolkit(const char *streamerName=nullptr, const char *mode="recreate");
119  virtual ~AliTMinuitToolkit();
120  // streaming intermediate results
121  void SetStreamer(const char *streamerName, const char *mode="recreate");
123  //
124  void ClearData();
125  void Fit(Option_t* option = "");
126  void FitHistogram(const TH1 *his, Option_t* option = "");
127  void FitGraph(const TGraph * gr, Option_t* option = "");
128  Long64_t FillFitter(TTree * inputTree, TString values, TString variables, TString selection, Int_t firstEntry, Int_t nEntries, Bool_t doReset=kTRUE);
129  TString GetFitFunctionAsAlias(Option_t *option="", TTree * tree= nullptr);
130  void Bootstrap(ULong_t nIter, const char* reportName, Option_t *option= nullptr);
131  void TwoFoldCrossValidation(UInt_t nIter, const char*reportName, Option_t *option= nullptr);
132  void MISAC(Int_t nFitPoints, UInt_t nIter, const char*reportName, Option_t *option= nullptr);
133  //
134  void SetFitFunction(TF1 * formula, Bool_t doReset);
135  void SetInitialParam(const TMatrixD * initialParam);
136  TMatrixD * GetInitialParam() {return fInitialParam;}
137  void SetLogLikelihoodFunction(TF1 *logLikelihood){fLogLikelihoodFunction=logLikelihood;}
138  void SetFitAlgorithm(Char_t *const name) {fFitAlgorithm=name;};
139  void SetMaxCalls(Int_t calls) {fMaxCalls = calls;};
140  void SetTolerance(Double_t tol) {fPrecision = tol;};
141  void SetPoints(TMatrixD * points) {fPoints = points;};
142  void SetValues(TMatrixD * values) {fValues = values;};
143  static void FitterFCN(int &nParam, double *info, double &chi2, double *gin, int iflag);
144  void EnableRobust(Bool_t b) {fIsHuberCost = b;};
145  //
146  const TMatrixD * GetPoints() const {return fPoints;};
147  const TMatrixD * GetValues() const {return fValues;};
148  TF1 * GetFormula() const {return fFormula;};
149  const TF1 *GetLogLikelihoodFunction() const { return fLogLikelihoodFunction;}
150  const TVectorD * GetParameters() const {return fParam;};
151  const TVectorD * GetRMSEstimator() const {return fRMSEstimator;};
152  const TMatrixD * GetCovarianceMatrix() const {return fCovar;};
153  const TMatrixD * GetMISACEstimators() const {return fMISACEstimators;};
154 
155  Double_t GetChisquare(){return fChi2;}
156  Int_t GetStatus(){return fStatus;}
157  Bool_t IsHuberCost() const {return fIsHuberCost;};
159  void SetVerbose(Int_t verbose){fVerbose=verbose;}
160  Bool_t IsVerbose(Int_t flag){ return (fVerbose&flag)>0;}
161  // static fitting
162  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);
163  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);
164  static AliTMinuitToolkit * GetPredefinedFitter(const std::string & fitterName){return fPredefinedFitters[fitterName];}
165  static void SetPredefinedFitter(const std::string & fitterName, AliTMinuitToolkit * fitter){fPredefinedFitters[fitterName]=fitter;}
166  static void RegisterDefaultFitters();
167  static AliTMinuitToolkit * RegisterPlaneFitter(Int_t nPlanes, Int_t logLikeType);
168  // helper functions for TFormula and TTreeFormula drawing
169  static Double_t GaussCachyLogLike(const Double_t *x, const Double_t *p);
170  static Double_t HuberLogLike(const Double_t *x, const Double_t *p);
171  static Double_t PseudoHuberLogLike(const Double_t *x, const Double_t *p);
172  static Double_t GausKurtosisSkewness(const Double_t *x, const Double_t *p);
173  static Double_t RndmGaus(Double_t mean=0, Double_t sigma=1);
174  static Double_t RndmLandau(Double_t mean=0, Double_t sigma=1);
175  static void SetFunctionDrawOption(TF1 *fun, Option_t *option); // parse draw option and set function attributes
176 private:
177  //
178  AliTMinuitToolkit(const AliTMinuitToolkit&); // private copy constructor -- suppress warnings
179  AliTMinuitToolkit& operator=(const AliTMinuitToolkit&); // private assignment operator
180  //
181  TTreeSRedirector *fStreamer; // !pointer to the streamer
182  Int_t fVerbose; // verbosity flag
183  TF1 * fFormula; // formula of the fitted function
184  TF1 * fLogLikelihoodFunction; // Log likelihood (const) function
185  TString fFitAlgorithm; // fit algorithm for TMinuit: migrad, simplex, ...
186  TMatrixD * fPoints; // !points - IsOwner
187  TMatrixD * fValues; // !matrix of values and their errors - IsOwner
188  TObjArray * fVarNames; // !variable names
189  TObjArray * fValueNames; // value names
190  TArrayI fPointIndex; // point index - to specify points for likelihood calculation (used for Bootstrap and CrossValidation)
191  TMatrixD * fInitialParam; // limits of the parameters (up, down)
192  TVectorD * fParam; // parameter values
193  TVectorD * fRMSEstimator; // parameter spread as estimated in Bootstrap resp. TwoFold cross validation
194  TMatrixD * fMISACEstimators; // MISAC estimators
195  TMatrixD * fCovar; // covariance matrix
196  Int_t fStatus; // fStatus
197  Double_t fChi2; // chi-square
198  Int_t fMaxCalls; // maximum number of calls, default value depends on fit algorithm
199  Double_t fPrecision; // normalized tolerance, default value depends on fit algorithm
200  Bool_t fIsHuberCost; // switch on/off robust option, default: kFalse
201  Int_t fFCNCounter; // function call counter
202  static std::map<std::string, AliTMinuitToolkit*> fPredefinedFitters;
203 ClassDef(AliTMinuitToolkit,2);
204 };
205 
206 
207 
208 #endif
AliTMinuitToolkit(const char *streamerName=nullptr, const char *mode="recreate")
TBrowser b
Definition: RunAnaESD.C:12
void FitHistogram(const TH1 *his, Option_t *option="")
Fit a one dimensional histogram.
void SetTolerance(Double_t tol)
void EnableRobust(Bool_t b)
void SetInitialParam(const TMatrixD *initialParam)
Set initial parameters matrix.
static Double_t RndmGaus(Double_t mean=0, Double_t sigma=1)
Random gaus generator as static function (to use in root TFormula, TTreeFormula)
const TVectorD * GetRMSEstimator() const
const TMatrixD * GetPoints() const
void Fit(Option_t *option="")
Internal function that calls the fitter.
TTreeSRedirector * fStreamer
#define TObjArray
TMatrixD * GetInitialParam()
TObjArray * GetListOfVariables()
const TMatrixD * GetMISACEstimators() const
const TMatrixD * GetCovarianceMatrix() const
static AliTMinuitToolkit * GetPredefinedFitter(const std::string &fitterName)
Float_t p[]
Definition: kNNTest.C:133
static Double_t RndmLandau(Double_t mean=0, Double_t sigma=1)
Random gaus generator as static function (to use in root TFormula, TTreeFormula)
static Double_t PseudoHuberLogLike(const Double_t *x, const Double_t *p)
Bool_t IsVerbose(Int_t flag)
Bool_t IsHuberCost() const
const TVectorD * GetParameters() const
AliTMinuitToolkit & operator=(const AliTMinuitToolkit &)
AliExternalInfo info
void SetValues(TMatrixD *values)
static Double_t HuberLogLike(const Double_t *x, const Double_t *p)
TTree * tree
Double_t chi2
Definition: AnalyzeLaser.C:7
TString GetFitFunctionAsAlias(Option_t *option="", TTree *tree=nullptr)
Format alias string for the for tree alias or for fit title.
void SetFitAlgorithm(Char_t *const name)
static void SetFunctionDrawOption(TF1 *fun, Option_t *option)
TGraph * gr
Definition: CalibTime.C:25
static void SetPredefinedFitter(const std::string &fitterName, AliTMinuitToolkit *fitter)
void SetVerbose(Int_t verbose)
void Bootstrap(ULong_t nIter, const char *reportName, Option_t *option=nullptr)
Bootstrap minimization.
void FitGraph(const TGraph *gr, Option_t *option="")
Fit a TGraph object.
void SetMaxCalls(Int_t calls)
TF1 * GetFormula() const
static Double_t GaussCachyLogLike(const Double_t *x, const Double_t *p)
static std::map< std::string, AliTMinuitToolkit * > fPredefinedFitters
const TMatrixD * GetValues() const
void SetStreamer(const char *streamerName, const char *mode="recreate")
Set stream for writing info into TTree. see examples(STEER/STEERBase/TTreeStream.cxx) ...
const TF1 * GetLogLikelihoodFunction() const
The AliTMinuitToolkit serves as an easy to use interface for the TMinuit.
void SetLogLikelihoodFunction(TF1 *logLikelihood)
static void RegisterDefaultFitters()
TTree * inputTree
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.
void TwoFoldCrossValidation(UInt_t nIter, const char *reportName, Option_t *option=nullptr)
void SetFitFunction(TF1 *formula, Bool_t doReset)
Set formula of the fitted function and reset parameter. See doReset.
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.
TMatrixD * fMISACEstimators
TTreeSRedirector * GetStreamer()
static AliTMinuitToolkit * RegisterPlaneFitter(Int_t nPlanes, Int_t logLikeType)
Register default fitters.
static Double_t GausKurtosisSkewness(const Double_t *x, const Double_t *p)
Gaus convoluted with Erf to emulate kurtosis and skewness.
void SetPoints(TMatrixD *points)
static void FitterFCN(int &nParam, double *info, double &chi2, double *gin, int iflag)
Internal function used in the TMinuit minimization.