ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
multidimfit.C
Go to the documentation of this file.
1 //Multi-Dimensional Parametrisation and Fitting
2 //Authors: Rene Brun, Christian Holm Christensen
3 
4 #include "Riostream.h"
5 #include "TROOT.h"
6 #include "TApplication.h"
7 #include "TCanvas.h"
8 #include "TH1.h"
9 #include "TSystem.h"
10 #include "TBrowser.h"
11 #include "TFile.h"
12 #include "TRandom.h"
13 #include "TMultiDimFit.h"
14 #include "TVectorD.h"
15 #include "TMath.h"
16 
17 
18 //____________________________________________________________________
20 {
21  // Make data points
22  Double_t upp[5] = { 10, 10, 10, 10, 1 };
23  Double_t low[5] = { 0, 0, 0, 0, .1 };
24  for (int i = 0; i < 4; i++)
25  x[i] = (upp[i] - low[i]) * gRandom->Rndm() + low[i];
26 
27  d = x[0] * TMath::Sqrt(x[1] * x[1] + x[2] * x[2] + x[3] * x[3]);
28 
29  e = gRandom->Gaus(upp[4],low[4]);
30 }
31 
32 //____________________________________________________________________
34 {
35  //Compare results with reference run
36 
37 
38  // the right coefficients (before fit)
39  double GoodCoeffsNoFit[] = {
40  -4.37056,
41  43.1468,
42  13.432,
43  13.4632,
44  13.3964,
45  13.328,
46  13.3016,
47  13.3519,
48  4.49724,
49  4.63876,
50  4.89036,
51  -3.69982,
52  -3.98618,
53  -3.86195,
54  4.36054,
55  -4.02597,
56  4.57037,
57  4.69845,
58  2.83819,
59  -3.48855,
60  -3.97612
61 };
62 
63  // the right coefficients (after fit)
64  double GoodCoeffs[] = {
65  -4.399,
66  43.15,
67  13.41,
68  13.49,
69  13.4,
70  13.23,
71  13.34,
72  13.29,
73  4.523,
74  4.659,
75  4.948,
76  -4.026,
77  -4.045,
78  -3.939,
79  4.421,
80  -4.006,
81  4.626,
82  4.378,
83  3.516,
84  -4.111,
85  -3.823,
86 };
87 
88 // Good Powers
89  int GoodPower[] = {
90  1, 1, 1, 1,
91  2, 1, 1, 1,
92  1, 1, 1, 2,
93  1, 1, 2, 1,
94  1, 2, 1, 1,
95  2, 2, 1, 1,
96  2, 1, 1, 2,
97  2, 1, 2, 1,
98  1, 1, 1, 3,
99  1, 3, 1, 1,
100  1, 1, 5, 1,
101  1, 1, 2, 2,
102  1, 2, 1, 2,
103  1, 2, 2, 1,
104  2, 1, 1, 3,
105  2, 2, 1, 2,
106  2, 1, 3, 1,
107  2, 3, 1, 1,
108  1, 2, 2, 2,
109  2, 1, 2, 2,
110  2, 2, 2, 1
111 };
112 
113  Int_t nc = fit->GetNCoefficients();
114  Int_t nv = fit->GetNVariables();
115  const Int_t *powers = fit->GetPowers();
116  const Int_t *pindex = fit->GetPowerIndex();
117  if (nc != 21) return 1;
118  const TVectorD *coeffs = fit->GetCoefficients();
119  int k = 0;
120  for (Int_t i=0;i<nc;i++) {
121  if (doFit) {
122  if (!TMath::AreEqualRel((*coeffs)[i],GoodCoeffs[i],1e-3)) return 2;
123  }
124  else {
125  if (TMath::Abs((*coeffs)[i] - GoodCoeffsNoFit[i]) > 5e-5) return 2;
126  }
127  for (Int_t j=0;j<nv;j++) {
128  if (powers[pindex[i]*nv+j] != GoodPower[k]) return 3;
129  k++;
130  }
131  }
132 
133  // now test the result of the generated function
134  gROOT->ProcessLine(".L MDF.C");
135 
136  Double_t refMDF = (doFit) ? 43.95 : 43.98;
137  // this does not work in CLing since the function is not defined
138  //Double_t x[] = {5,5,5,5};
139  //Double_t rMDF = MDF(x);
140  //LM: need to return the address of the result since it is casted to a long (this should not be in a tutorial !)
141  Long_t iret = gROOT->ProcessLine(" Double_t x[] = {5,5,5,5}; double result=MDF(x); &result;");
142  Double_t rMDF = * ( (Double_t*)iret);
143  //printf("%f\n",rMDF);
144  if (TMath::Abs(rMDF -refMDF) > 1e-2) return 4;
145  return 0;
146 }
147 
148 //____________________________________________________________________
149 Int_t multidimfit(bool doFit = true)
150 {
151 
152  cout << "*************************************************" << endl;
153  cout << "* Multidimensional Fit *" << endl;
154  cout << "* *" << endl;
155  cout << "* By Christian Holm <cholm@nbi.dk> 14/10/00 *" << endl;
156  cout << "*************************************************" << endl;
157  cout << endl;
158 
159  // Initialize global TRannom object.
160  gRandom = new TRandom();
161 
162  // Open output file
163  TFile* output = new TFile("mdf.root", "RECREATE");
164 
165  // Global data parameters
166  Int_t nVars = 4;
167  Int_t nData = 500;
168  Double_t x[4];
169 
170  // make fit object and set parameters on it.
171  TMultiDimFit* fit = new TMultiDimFit(nVars, TMultiDimFit::kMonomials,"v");
172 
173  Int_t mPowers[] = { 6 , 6, 6, 6 };
174  fit->SetMaxPowers(mPowers);
175  fit->SetMaxFunctions(1000);
176  fit->SetMaxStudy(1000);
177  fit->SetMaxTerms(30);
178  fit->SetPowerLimit(1);
179  fit->SetMinAngle(10);
180  fit->SetMaxAngle(10);
181  fit->SetMinRelativeError(.01);
182 
183  // variables to hold the temporary input data
184  Double_t d;
185  Double_t e;
186 
187  // Print out the start parameters
188  fit->Print("p");
189 
190  printf("======================================\n");
191 
192  // Create training sample
193  Int_t i;
194  for (i = 0; i < nData ; i++) {
195 
196  // Make some data
197  makeData(x,d,e);
198 
199  // Add the row to the fit object
200  fit->AddRow(x,d,e);
201  }
202 
203  // Print out the statistics
204  fit->Print("s");
205 
206  // Book histograms
207  fit->MakeHistograms();
208 
209  // Find the parameterization
210  fit->FindParameterization();
211 
212  // Print coefficents
213  fit->Print("rc");
214 
215  // Get the min and max of variables from the training sample, used
216  // for cuts in test sample.
217  Double_t *xMax = new Double_t[nVars];
218  Double_t *xMin = new Double_t[nVars];
219  for (i = 0; i < nVars; i++) {
220  xMax[i] = (*fit->GetMaxVariables())(i);
221  xMin[i] = (*fit->GetMinVariables())(i);
222  }
223 
224  nData = fit->GetNCoefficients() * 100;
225  Int_t j;
226 
227  // Create test sample
228  for (i = 0; i < nData ; i++) {
229  // Make some data
230  makeData(x,d,e);
231 
232  for (j = 0; j < nVars; j++)
233  if (x[j] < xMin[j] || x[j] > xMax[j])
234  break;
235 
236  // If we get through the loop above, all variables are in range
237  if (j == nVars)
238  // Add the row to the fit object
239  fit->AddTestRow(x,d,e);
240  else
241  i--;
242  }
243  //delete gRandom;
244 
245  // Test the parameterizatio and coefficents using the test sample.
246  if (doFit)
247  fit->Fit("M");
248 
249  // Print result
250  fit->Print("fc v");
251 
252  // Write code to file
253  fit->MakeCode();
254 
255  // Write histograms to disk, and close file
256  output->Write();
257  output->Close();
258  delete output;
259 
260  // Compare results with reference run
262  if (!compare) {
263  printf("\nmultidimfit .............................................. OK\n");
264  } else {
265  printf("\nmultidimfit .............................................. fails case %d\n",compare);
266  }
267 
268  // We're done
269  delete fit;
270  return compare;
271 }
const TVectorD * GetCoefficients() const
Definition: TMultiDimFit.h:152
Int_t multidimfit(bool doFit=true)
Definition: multidimfit.C:149
void doFit(int n, const char *fitter)
const TVectorD * GetMaxVariables() const
Definition: TMultiDimFit.h:164
virtual Double_t Rndm(Int_t i=0)
Machine independent random number generator.
Definition: TRandom.cxx:512
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
Samples a random number from the standard Normal (Gaussian) Distribution with the given mean and sigm...
Definition: TRandom.cxx:235
virtual void Fit(Option_t *option="")
Try to fit the found parameterisation to the test sample.
void makeData(Double_t *x, Double_t &d, Double_t &e)
Definition: multidimfit.C:19
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format...
Definition: TFile.h:45
virtual void Print(Option_t *option="ps") const
Print statistics etc.
void SetMaxPowers(const Int_t *powers)
Set the maximum power to be considered in the fit for each variable.
int compare(double v1, double v2, const std::string &name="", double scale=1.0)
#define gROOT
Definition: TROOT.h:344
int Int_t
Definition: RtypesCore.h:41
Int_t GetNVariables() const
Definition: TMultiDimFit.h:171
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
Int_t GetNCoefficients() const
Definition: TMultiDimFit.h:172
Double_t x[n]
Definition: legend1.C:17
void SetPowerLimit(Double_t limit=1e-3)
Set the user parameter for the function selection.
This is the base class for the ROOT Random number generators.
Definition: TRandom.h:29
int d
Definition: tornado.py:11
void SetMinAngle(Double_t angle=1)
Set the min angle (in degrees) between a new candidate function and the subspace spanned by the previ...
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...
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, and optionally, the square error in the dependent quantity, to the training sample to be used for the parameterization.
void SetMaxTerms(Int_t terms)
Definition: TMultiDimFit.h:211
const TVectorD * GetMinVariables() const
Definition: TMultiDimFit.h:170
Bool_t AreEqualRel(Double_t af, Double_t bf, Double_t relPrec)
Definition: TMath.h:196
virtual Int_t Write(const char *name=0, Int_t opt=0, Int_t bufsiz=0)
Write memory objects to this file.
Definition: TFile.cxx:2248
Multidimensional Fits in ROOT.
Definition: TMultiDimFit.h:25
virtual void FindParameterization(Option_t *option="")
Find the parameterization.
R__EXTERN TRandom * gRandom
Definition: TRandom.h:62
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, and optionally, the square error in the dependent quantity, to the test sample to be used for the test of the parameterization.
virtual void MakeHistograms(Option_t *option="A")
Make histograms of the result of the analysis.
long Long_t
Definition: RtypesCore.h:50
const Int_t * GetPowers() const
Definition: TMultiDimFit.h:176
double Double_t
Definition: RtypesCore.h:55
void SetMinRelativeError(Double_t error)
Set the acceptable relative error for when sum of square residuals is considered minimized.
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
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 .C.
int CompareResults(TMultiDimFit *fit, bool doFit)
Definition: multidimfit.C:33
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
static void output(int code)
Definition: gifencode.c:226
void SetMaxStudy(Int_t n)
Definition: TMultiDimFit.h:210
void SetMaxFunctions(Int_t n)
Definition: TMultiDimFit.h:208
Int_t * GetPowerIndex() const
Definition: TMultiDimFit.h:174
virtual void Close(Option_t *option="")
Close a file.
Definition: TFile.cxx:898