Logo ROOT  
Reference Guide
TMultiDimFit.h
Go to the documentation of this file.
1// @(#)root/hist:$Id$
2// Author: Christian Holm Christensen 07/11/2000
3
4#ifndef ROOT_TMultiDimFit
5#define ROOT_TMultiDimFit
6
7#include "TNamed.h"
8#include "TVectorD.h"
9#include "TMatrixD.h"
10#include "TList.h"
11#include "TVirtualFitter.h"
12
13class TBrowser;
14
15class TMultiDimFit : public TNamed {
16
17public:
22 };
23
24private:
25 static TMultiDimFit* fgInstance; // Static instance
26protected:
27
28 TVectorD fQuantity; // Training sample, dependent quantity
29 TVectorD fSqError; // Training sample, error in quantity
30 Double_t fMeanQuantity; // Mean of dependent quantity
31 Double_t fMaxQuantity; // Max value of dependent quantity
32 Double_t fMinQuantity; // Min value of dependent quantity
33 Double_t fSumSqQuantity; // SumSquare of dependent quantity
34 Double_t fSumSqAvgQuantity; // Sum of squares away from mean
35
36 TVectorD fVariables; // Training sample, independent variables
37 Int_t fNVariables; // Number of independent variables
38 TVectorD fMeanVariables; // mean value of independent variables
39 TVectorD fMaxVariables; // max value of independent variables
40 TVectorD fMinVariables; // min value of independent variables
41
42 Int_t fSampleSize; // Size of training sample
43
44 TVectorD fTestQuantity; // Test sample, dependent quantity
45 TVectorD fTestSqError; // Test sample, Error in quantity
46 TVectorD fTestVariables; // Test sample, independent variables
47
48 Int_t fTestSampleSize; // Size of test sample
49
50 Double_t fMinAngle; // Min angle for acepting new function
51 Double_t fMaxAngle; // Max angle for acepting new function
52 Int_t fMaxTerms; // Max terms expected in final expr.
53 Double_t fMinRelativeError; // Min relative error accepted
54 Int_t *fMaxPowers; // [fNVariables] maximum powers
55 Double_t fPowerLimit; // Control parameter
56
57
58 TMatrixD fFunctions; // Functions evaluated over sample
59 Int_t fMaxFunctions; // max number of functions
60 Int_t *fFunctionCodes; // [fMaxFunctions] acceptance code
61 Int_t fMaxStudy; // max functions to study
62 Int_t fMaxFuncNV; // fMaxFunctions*fNVariables
63
64 TMatrixD fOrthFunctions; // As above, but orthogonalised
65 TVectorD fOrthFunctionNorms; // Norm of the evaluated functions
66
67
68 Int_t *fMaxPowersFinal; // [fNVariables] maximum powers from fit;
69 Int_t *fPowers; // [fMaxFuncNV] where fMaxFuncNV = fMaxFunctions*fNVariables
70 Int_t *fPowerIndex; // [fMaxTerms] Index of accepted powers
71
72 TVectorD fResiduals; // Vector of the final residuals
73 Double_t fMaxResidual; // Max redsidual value
74 Double_t fMinResidual; // Min redsidual value
75 Int_t fMaxResidualRow; // Row giving max residual
76 Int_t fMinResidualRow; // Row giving min residual
77 Double_t fSumSqResidual; // Sum of Square residuals
78
79 Int_t fNCoefficients; // Dimension of model coefficients
80 TVectorD fOrthCoefficients; // The model coefficients
82 TVectorD fCoefficients; // Vector of the final coefficients
83 TVectorD fCoefficientsRMS; // Vector of RMS of coefficients
84 Double_t fRMS; // Root mean square of fit
85 Double_t fChi2; // Chi square of fit
86 Int_t fParameterisationCode; // Exit code of parameterisation
87
88 Double_t fError; // Error from parameterization
89 Double_t fTestError; // Error from test
90 Double_t fPrecision; // Relative precision of param
91 Double_t fTestPrecision; // Relative precision of test
92 Double_t fCorrelationCoeff; // Multi Correlation coefficient
93 TMatrixD fCorrelationMatrix; // Correlation matrix
94 Double_t fTestCorrelationCoeff; // Multi Correlation coefficient
95
96 TList* fHistograms; // List of histograms
97 Byte_t fHistogramMask; // Bit pattern of hisograms used
98 Int_t fBinVarX; // Number of bin in independent variables
99 Int_t fBinVarY; // Number of bin in dependent variables
100
101 TVirtualFitter* fFitter; //! Fit object (MINUIT)
102
103 EMDFPolyType fPolyType; // Type of polynomials to use
104 Bool_t fShowCorrelation; // print correlation matrix
105 Bool_t fIsUserFunction; // Flag for user defined function
107
108 virtual Double_t EvalFactor(Int_t p, Double_t x) const;
109 virtual Double_t EvalControl(const Int_t *powers) const;
110 virtual void MakeCoefficientErrors();
111 virtual void MakeCorrelation();
113 virtual void MakeCoefficients();
114 virtual void MakeCandidates();
115 virtual void MakeNormalized();
116 virtual void MakeParameterization();
117 virtual void MakeRealCode(const char *filename,
118 const char *classname,
119 Option_t *option="");
120 virtual Bool_t Select(const Int_t *iv);
121 virtual Bool_t TestFunction(Double_t squareResidual,
122 Double_t dResidur);
123public:
124 TMultiDimFit();
125 TMultiDimFit(Int_t dimension,
127 Option_t *option="");
128 virtual ~TMultiDimFit();
129
130 virtual void AddRow(const Double_t *x, Double_t D, Double_t E=0);
131 virtual void AddTestRow(const Double_t *x, Double_t D, Double_t E=0);
132 virtual void Browse(TBrowser* b);
133 virtual void Clear(Option_t *option=""); // *MENU*
134 virtual void Draw(Option_t * ="d") { }
135 virtual Double_t Eval(const Double_t *x, const Double_t *coeff=0) const;
136 virtual Double_t EvalError(const Double_t *x, const Double_t *coeff=0) const;
137 virtual void FindParameterization(Option_t* option=""); // *MENU*
138 virtual void Fit(Option_t *option=""); // *MENU*
139
140 Double_t GetChi2() const { return fChi2; }
142 const TVectorD* GetCoefficients() const { return &fCoefficients; }
143 const TVectorD* GetCoefficientsRMS() const { return &fCoefficientsRMS; }
144 Double_t GetError() const { return fError; }
146 const TMatrixD* GetFunctions() const { return &fFunctions; }
147 virtual TList* GetHistograms() const { return fHistograms; }
148 Double_t GetMaxAngle() const { return fMaxAngle; }
150 Int_t* GetMaxPowers() const { return fMaxPowers; }
152 Int_t GetMaxStudy() const { return fMaxStudy; }
153 Int_t GetMaxTerms() const { return fMaxTerms; }
154 const TVectorD* GetMaxVariables() const { return &fMaxVariables; }
156 const TVectorD* GetMeanVariables() const { return &fMeanVariables; }
157 Double_t GetMinAngle() const { return fMinAngle; }
160 const TVectorD* GetMinVariables() const { return &fMinVariables; }
161 Int_t GetNVariables() const { return fNVariables; }
163 Int_t GetPolyType() const { return fPolyType; }
164 Int_t* GetPowerIndex() const { return fPowerIndex; }
166 const Int_t* GetPowers() const { return fPowers; }
167 Double_t GetPrecision() const { return fPrecision; }
168 const TVectorD* GetQuantity() const { return &fQuantity; }
174 Double_t GetRMS() const { return fRMS; }
175 Int_t GetSampleSize() const { return fSampleSize; }
176 const TVectorD* GetSqError() const { return &fSqError; }
179 Double_t GetTestError() const { return fTestError; }
181 const TVectorD* GetTestQuantity() const { return &fTestQuantity; }
183 const TVectorD* GetTestSqError() const { return &fTestSqError; }
184 const TVectorD* GetTestVariables() const { return &fTestVariables; }
185 const TVectorD* GetVariables() const { return &fVariables; }
186
187 static TMultiDimFit* Instance();
188 virtual Bool_t IsFolder() const { return kTRUE; }
189 virtual Double_t MakeChi2(const Double_t* coeff=0);
190 virtual void MakeCode(const char *functionName="MDF", Option_t *option=""); // *MENU*
191 virtual void MakeHistograms(Option_t* option="A"); // *MENU*
192 virtual void MakeMethod(const Char_t* className="MDF", Option_t* option=""); // *MENU*
193 virtual void Print(Option_t *option="ps") const; // *MENU*
194
195 void SetBinVarX(Int_t nbbinvarx) {fBinVarX = nbbinvarx;}
196 void SetBinVarY(Int_t nbbinvary) {fBinVarY = nbbinvary;}
197 void SetMaxAngle(Double_t angle=0);
199 void SetMaxPowers(const Int_t *powers);
201 void SetMaxTerms(Int_t terms) { fMaxTerms = terms; }
202 void SetMinRelativeError(Double_t error);
203 void SetMinAngle(Double_t angle=1);
204 void SetPowerLimit(Double_t limit=1e-3);
205 virtual void SetPowers(const Int_t *powers, Int_t terms);
206
207 ClassDef(TMultiDimFit,2) // Multi dimensional fit class
208}
209;
210#endif
#define b(i)
Definition: RSha256.hxx:100
#define e(i)
Definition: RSha256.hxx:103
unsigned char Byte_t
Definition: RtypesCore.h:62
int Int_t
Definition: RtypesCore.h:43
char Char_t
Definition: RtypesCore.h:31
double Double_t
Definition: RtypesCore.h:57
const Bool_t kTRUE
Definition: RtypesCore.h:89
const char Option_t
Definition: RtypesCore.h:64
#define ClassDef(name, id)
Definition: Rtypes.h:322
int type
Definition: TGX11.cxx:120
Using a TBrowser one can browse all ROOT objects.
Definition: TBrowser.h:37
A doubly linked list.
Definition: TList.h:44
Multidimensional Fits in ROOT.
Definition: TMultiDimFit.h:15
const TVectorD * GetQuantity() const
Definition: TMultiDimFit.h:168
Double_t GetSumSqAvgQuantity() const
Definition: TMultiDimFit.h:177
virtual void MakeCorrelation()
PRIVATE METHOD: Compute the correlation matrix.
Int_t GetNCoefficients() const
Definition: TMultiDimFit.h:162
Double_t GetMinQuantity() const
Definition: TMultiDimFit.h:158
void SetMaxStudy(Int_t n)
Definition: TMultiDimFit.h:200
TVectorD fCoefficients
Definition: TMultiDimFit.h:82
Double_t fPrecision
Definition: TMultiDimFit.h:90
Double_t GetMaxAngle() const
Definition: TMultiDimFit.h:148
Byte_t fHistogramMask
Definition: TMultiDimFit.h:97
void SetBinVarX(Int_t nbbinvarx)
Definition: TMultiDimFit.h:195
virtual void MakeCoefficientErrors()
PRIVATE METHOD: Compute the errors on the coefficients.
void SetMaxFunctions(Int_t n)
Definition: TMultiDimFit.h:198
Double_t GetTestPrecision() const
Definition: TMultiDimFit.h:180
TMatrixD fOrthCurvatureMatrix
Definition: TMultiDimFit.h:81
Bool_t fShowCorrelation
Definition: TMultiDimFit.h:104
Bool_t fIsUserFunction
Definition: TMultiDimFit.h:105
Double_t GetError() const
Definition: TMultiDimFit.h:144
Double_t GetRMS() const
Definition: TMultiDimFit.h:174
Int_t fNVariables
Definition: TMultiDimFit.h:37
Int_t GetSampleSize() const
Definition: TMultiDimFit.h:175
Int_t fMaxResidualRow
Definition: TMultiDimFit.h:75
Double_t fMinQuantity
Definition: TMultiDimFit.h:32
virtual Double_t Eval(const Double_t *x, const Double_t *coeff=0) const
Evaluate parameterization at point x.
Double_t fSumSqQuantity
Definition: TMultiDimFit.h:33
TVectorD fMaxVariables
Definition: TMultiDimFit.h:39
Double_t GetPrecision() const
Definition: TMultiDimFit.h:167
const TVectorD * GetVariables() const
Definition: TMultiDimFit.h:185
Int_t GetResidualMinRow() const
Definition: TMultiDimFit.h:172
const TVectorD * GetMeanVariables() const
Definition: TMultiDimFit.h:156
void SetMaxTerms(Int_t terms)
Definition: TMultiDimFit.h:201
Double_t GetSumSqQuantity() const
Definition: TMultiDimFit.h:178
Int_t fSampleSize
Definition: TMultiDimFit.h:42
Double_t GetResidualSumSq() const
Definition: TMultiDimFit.h:173
Double_t fSumSqAvgQuantity
Definition: TMultiDimFit.h:34
Int_t fMaxTerms
Definition: TMultiDimFit.h:52
Int_t * fPowers
Definition: TMultiDimFit.h:69
Double_t GetPowerLimit() const
Definition: TMultiDimFit.h:165
virtual void MakeNormalized()
PRIVATE METHOD: Normalize data to the interval [-1;1].
Int_t * fMaxPowers
Definition: TMultiDimFit.h:54
Int_t GetNVariables() const
Definition: TMultiDimFit.h:161
Int_t fMaxStudy
Definition: TMultiDimFit.h:61
TVectorD fTestQuantity
Definition: TMultiDimFit.h:44
TVectorD fTestSqError
Definition: TMultiDimFit.h:45
Int_t fMaxFunctions
Definition: TMultiDimFit.h:59
Int_t fBinVarY
Definition: TMultiDimFit.h:99
static TMultiDimFit * Instance()
Return the static instance.
virtual ~TMultiDimFit()
Destructor.
virtual void MakeCode(const char *functionName="MDF", Option_t *option="")
Generate the file <filename> with .C appended if argument doesn't end in .cxx or ....
virtual void MakeCandidates()
PRIVATE METHOD: Create list of candidate functions for the parameterisation.
const TVectorD * GetSqError() const
Definition: TMultiDimFit.h:176
const TVectorD * GetMaxVariables() const
Definition: TMultiDimFit.h:154
Double_t fTestError
Definition: TMultiDimFit.h:89
Double_t fMinResidual
Definition: TMultiDimFit.h:74
void SetPowerLimit(Double_t limit=1e-3)
Set the user parameter for the function selection.
const TMatrixD * GetFunctions() const
Definition: TMultiDimFit.h:146
virtual void FindParameterization(Option_t *option="")
Find the parameterization.
virtual Double_t MakeGramSchmidt(Int_t function)
PRIVATE METHOD: Make Gram-Schmidt orthogonalisation.
Int_t * GetMaxPowers() const
Definition: TMultiDimFit.h:150
const TVectorD * GetMinVariables() const
Definition: TMultiDimFit.h:160
TVectorD fSqError
Definition: TMultiDimFit.h:29
virtual void AddRow(const Double_t *x, Double_t D, Double_t E=0)
Add a row consisting of fNVariables independent variables, the known, dependent quantity,...
Int_t GetTestSampleSize() const
Definition: TMultiDimFit.h:182
Int_t GetResidualMaxRow() const
Definition: TMultiDimFit.h:171
TVectorD fQuantity
Definition: TMultiDimFit.h:28
Double_t GetMinAngle() const
Definition: TMultiDimFit.h:157
TVectorD fTestVariables
Definition: TMultiDimFit.h:46
TVectorD fOrthCoefficients
Definition: TMultiDimFit.h:80
Int_t fTestSampleSize
Definition: TMultiDimFit.h:48
TVectorD fResiduals
Definition: TMultiDimFit.h:72
const Int_t * GetPowers() const
Definition: TMultiDimFit.h:166
virtual Bool_t TestFunction(Double_t squareResidual, Double_t dResidur)
PRIVATE METHOD: Test whether the currently considered function contributes to the fit.
TVectorD fVariables
Definition: TMultiDimFit.h:36
Double_t GetChi2() const
Definition: TMultiDimFit.h:140
Double_t fCorrelationCoeff
Definition: TMultiDimFit.h:92
void SetMaxPowers(const Int_t *powers)
Set the maximum power to be considered in the fit for each variable.
Int_t * fMaxPowersFinal
Definition: TMultiDimFit.h:68
Int_t fMinResidualRow
Definition: TMultiDimFit.h:76
virtual Bool_t IsFolder() const
Returns kTRUE in case object contains browsable objects (like containers or lists of other objects).
Definition: TMultiDimFit.h:188
Int_t GetMaxTerms() const
Definition: TMultiDimFit.h:153
virtual void SetPowers(const Int_t *powers, Int_t terms)
Define a user function.
Double_t fMaxQuantity
Definition: TMultiDimFit.h:31
virtual void MakeRealCode(const char *filename, const char *classname, Option_t *option="")
PRIVATE METHOD: This is the method that actually generates the code for the evaluation the parameteri...
Double_t fRMS
Definition: TMultiDimFit.h:84
TVirtualFitter * fFitter
Definition: TMultiDimFit.h:101
TMatrixD fFunctions
Definition: TMultiDimFit.h:58
Double_t GetResidualMin() const
Definition: TMultiDimFit.h:170
virtual void AddTestRow(const Double_t *x, Double_t D, Double_t E=0)
Add a row consisting of fNVariables independent variables, the known, dependent quantity,...
Int_t fParameterisationCode
Definition: TMultiDimFit.h:86
const TVectorD * GetTestVariables() const
Definition: TMultiDimFit.h:184
Double_t fMaxAngle
Definition: TMultiDimFit.h:51
Int_t fNCoefficients
Definition: TMultiDimFit.h:79
virtual void Fit(Option_t *option="")
Try to fit the found parameterisation to the test sample.
TMultiDimFit()
Empty CTOR. Do not use.
virtual TList * GetHistograms() const
Definition: TMultiDimFit.h:147
Double_t fTestCorrelationCoeff
Definition: TMultiDimFit.h:94
Int_t GetMaxStudy() const
Definition: TMultiDimFit.h:152
Double_t GetMeanQuantity() const
Definition: TMultiDimFit.h:155
Double_t fMinAngle
Definition: TMultiDimFit.h:50
const TMatrixD * GetCorrelationMatrix() const
Definition: TMultiDimFit.h:141
void SetMinRelativeError(Double_t error)
Set the acceptable relative error for when sum of square residuals is considered minimized.
Double_t fSumSqResidual
Definition: TMultiDimFit.h:77
Double_t GetResidualMax() const
Definition: TMultiDimFit.h:169
TVectorD fMeanVariables
Definition: TMultiDimFit.h:38
void SetBinVarY(Int_t nbbinvary)
Definition: TMultiDimFit.h:196
virtual void MakeMethod(const Char_t *className="MDF", Option_t *option="")
Generate the file <classname>MDF.cxx which contains the implementation of the method:
Double_t GetMinRelativeError() const
Definition: TMultiDimFit.h:159
void SetMaxAngle(Double_t angle=0)
Set the max angle (in degrees) between the initial data vector to be fitted, and the new candidate fu...
Double_t fTestPrecision
Definition: TMultiDimFit.h:91
const TVectorD * GetTestSqError() const
Definition: TMultiDimFit.h:183
virtual void MakeCoefficients()
PRIVATE METHOD: Invert the model matrix B, and compute final coefficients.
TMatrixD fCorrelationMatrix
Definition: TMultiDimFit.h:93
Int_t * GetPowerIndex() const
Definition: TMultiDimFit.h:164
Int_t * fPowerIndex
Definition: TMultiDimFit.h:70
TMatrixD fOrthFunctions
Definition: TMultiDimFit.h:64
virtual void Clear(Option_t *option="")
Clear internal structures and variables.
virtual Double_t EvalControl(const Int_t *powers) const
PRIVATE METHOD: Calculate the control parameter from the passed powers.
virtual void Browse(TBrowser *b)
Browse the TMultiDimFit object in the TBrowser.
virtual void MakeParameterization()
PRIVATE METHOD: Find the parameterization over the training sample.
EMDFPolyType fPolyType
Fit object (MINUIT)
Definition: TMultiDimFit.h:103
virtual Bool_t Select(const Int_t *iv)
Selection method.
Int_t GetMaxFunctions() const
Definition: TMultiDimFit.h:149
virtual Double_t EvalError(const Double_t *x, const Double_t *coeff=0) const
Evaluate parameterization error at point x.
TList * fHistograms
Definition: TMultiDimFit.h:96
Double_t GetTestError() const
Definition: TMultiDimFit.h:179
Double_t fChi2
Definition: TMultiDimFit.h:85
Double_t fMaxResidual
Definition: TMultiDimFit.h:73
Double_t fError
Definition: TMultiDimFit.h:88
Int_t fBinVarX
Definition: TMultiDimFit.h:98
Int_t * GetFunctionCodes() const
Definition: TMultiDimFit.h:145
virtual Double_t EvalFactor(Int_t p, Double_t x) const
PRIVATE METHOD: Evaluate function with power p at variable value x.
const TVectorD * GetTestQuantity() const
Definition: TMultiDimFit.h:181
void SetMinAngle(Double_t angle=1)
Set the min angle (in degrees) between a new candidate function and the subspace spanned by the previ...
virtual void Print(Option_t *option="ps") const
Print statistics etc.
TVectorD fMinVariables
Definition: TMultiDimFit.h:40
TVectorD fOrthFunctionNorms
Definition: TMultiDimFit.h:65
const TVectorD * GetCoefficientsRMS() const
Definition: TMultiDimFit.h:143
Double_t fMeanQuantity
Definition: TMultiDimFit.h:30
const TVectorD * GetCoefficients() const
Definition: TMultiDimFit.h:142
TVectorD fCoefficientsRMS
Definition: TMultiDimFit.h:83
Double_t fPowerLimit
Definition: TMultiDimFit.h:55
static TMultiDimFit * fgInstance
Definition: TMultiDimFit.h:25
Int_t GetPolyType() const
Definition: TMultiDimFit.h:163
virtual void MakeHistograms(Option_t *option="A")
Make histograms of the result of the analysis.
Bool_t fIsVerbose
Definition: TMultiDimFit.h:106
Int_t fMaxFuncNV
Definition: TMultiDimFit.h:62
Double_t GetMaxQuantity() const
Definition: TMultiDimFit.h:151
virtual void Draw(Option_t *="d")
Default Draw method for all objects.
Definition: TMultiDimFit.h:134
Double_t fMinRelativeError
Definition: TMultiDimFit.h:53
virtual Double_t MakeChi2(const Double_t *coeff=0)
Calculate Chi square over either the test sample.
Int_t * fFunctionCodes
Definition: TMultiDimFit.h:60
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:29
Abstract Base Class for Fitting.
Double_t x[n]
Definition: legend1.C:17
const Int_t n
Definition: legend1.C:16
void function(const Char_t *name_, T fun, const Char_t *docstring=0)
Definition: RExports.h:151
constexpr Double_t E()
Base of natural log:
Definition: TMath.h:97