GSI Object Oriented Online Offline (Go4)  GO4-5.3.2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
Example13.cxx
Go to the documentation of this file.
1 // $Id: Example13.cxx 934 2013-01-29 15:59:24Z linev $
2 //-----------------------------------------------------------------------
3 // The GSI Online Offline Object Oriented (Go4) Project
4 // Experiment Data Processing at EE department, GSI
5 //-----------------------------------------------------------------------
6 // Copyright (C) 2000- GSI Helmholtzzentrum für Schwerionenforschung GmbH
7 // Planckstr. 1, 64291 Darmstadt, Germany
8 // Contact: http://go4.gsi.de
9 //-----------------------------------------------------------------------
10 // This software can be used under the license agreements as stated
11 // in Go4License.txt file which is part of the distribution.
12 //-----------------------------------------------------------------------
13 
14 // Go4Fit Example N13
15 
16 #ifndef __CINT__
17 
18 #include "Riostream.h"
19 
20 #include "TMath.h"
21 #include "TH1.h"
22 #include "TH2.h"
23 #include "TCanvas.h"
24 #include "TCutG.h"
25 #include "TApplication.h"
26 #include "TRandom.h"
27 
28 #include "TGo4Fitter.h"
29 #include "TGo4FitDataHistogram.h"
30 #include "TGo4FitDataRidge.h"
31 #include "TGo4FitModelFormula.h"
32 
33 void Example13();
34 
35 int main(int argc, char **argv)
36 {
37  TApplication theApp("Application", 0, 0);
38 
39  Example13();
40 
41  theApp.Run();
42 
43  return 0;
44 }
45 
46 #endif
47 
48 void AddRidge(TH2D* histo, double a, double k, double b,
49  double xstart, int numsteps, int numpoints, double deviation)
50 {
51  TRandom rnd(1234);
52 
53  double x = xstart;
54 
55  for(int i=0;i<numsteps;i++) {
56  double y = a*exp(k*x) + b;
57 
58  double alf = TMath::ATan(k*a*exp(k*x)) + TMath::Pi()/2.;
59 
60  for(int j=0;j<numpoints;j++) {
61  double radius = rnd.Gaus(0.,deviation);
62  histo->Fill(x + radius*cos(alf), y + radius*sin(alf));
63  }
64 
65  x = x + 0.001*sin(alf);
66  }
67 }
68 
70 {
71  TH2D* histo = new TH2D("Histo","Ridge y = a*exp(k*x)+b",100,0.,10.,100,0.,10.);
72 
73  std::cout << "Generating histogram " << std::endl;
74 
75  AddRidge(histo, 10, -0.3, 2, 1.5, 10000, 100, 0.3);
76 
77  AddRidge(histo, 10, -0.5, 0, 0.5, 10000, 50, 0.2);
78 
79  AddRidge(histo, 10, -0.1, +2, 3.5, 10000, 50, 0.2);
80 
81  std::cout << " Done " << std::endl;
82 
83  return histo;
84 }
85 
86 TCutG* MakeCut()
87 {
88  Float_t xx[7] = { 1.0, 4.0, 9.0, 9.5, 5.0, 2.0, 1.0 };
89  Float_t yy[7] = { 8.0, 2.5, 2.0, 4.0, 6.0, 9.5, 8.0 };
90  return new TCutG("RangeCut", 7, xx, yy);
91 }
92 
93 void Example13()
94 {
95 // model histogram with ridge with exponential function
96  TH2D* histo = MakeRidgeHistogram();
97  new TCanvas("Canvas1");
98  histo->Draw("LEGO2");
99 
100 // create fitter, select fit function and not add standard actions list
101  TGo4Fitter fitter("Fitter", TGo4Fitter::ff_chi_square, kFALSE);
102 
103 // create data object to represent of two-dimensional histogram
104 // bins, which are less then 3 are excluded
105 // select middle ridge by TCutG object
106  TGo4FitDataHistogram* datah = new TGo4FitDataHistogram("histogram",histo);
107  datah->SetExcludeLessThen(3);
108  datah->AddRangeCut(MakeCut());
109 
110 // create ridge object, which select axis number 1 (y axis) as function of the rest coordinates
111  TGo4FitDataRidge* ridge = new TGo4FitDataRidge("ridge",datah,1);
112  fitter.AddData(ridge);
113 
114 // create formula function to approximate ridge function
115 // specific range and initial values for parameters need to be set
116  TGo4FitModelFormula* model = new TGo4FitModelFormula("RidgeFunc","parA*exp(parK*x)+parB",3,kFALSE);
117  model->SetRange(0,2,8.);
118  model->SetParsNames("parA","parK","parB");
119  model->SetParsValues(9.,-0.2,1.2);
120  fitter.AddModel("ridge", model);
121 
122 // add minuit action
123  fitter.AddSimpleMinuit();
124 
125 // use memory buffers for data objects
126  fitter.SetMemoryUsage(1);
127 
128 // execute actions
129  fitter.DoActions();
130 
131 // print parameters values
132  fitter.Print("Pars");
133 
134 // draw results of fit
135  fitter.Draw("#ridge");
136 }
void SetParsNames(const char *name0="Par0", const char *name1="Par1", const char *name2="Par2", const char *name3="Par3", const char *name4="Par4", const char *name5="Par5", const char *name6="Par6", const char *name7="Par7", const char *name8="Par8", const char *name9="Par9")
virtual void Print(Option_t *option) const
Definition: TGo4Fitter.cxx:707
void SetMemoryUsage(Int_t iMemoryUsage)
Definition: TGo4Fitter.cxx:72
virtual void Draw(Option_t *option)
Definition: TGo4Fitter.cxx:983
void SetExcludeLessThen(Double_t limit=0.)
Definition: TGo4FitData.h:107
int main(int argc, char **argv)
Definition: Example13.cxx:35
tuple a
Definition: go4init.py:12
TGo4FitData * AddData(TGo4FitData *d)
Definition: TGo4Fitter.cxx:118
void AddRidge(TH2D *histo, double a, double k, double b, double xstart, int numsteps, int numpoints, double deviation)
Definition: Example13.cxx:48
void SetRange(Int_t naxis, Double_t min, Double_t max)
void DoActions(Bool_t AllowFitterChange=kFALSE, TObjArray *Actions=0)
TH2D * MakeRidgeHistogram()
Definition: Example13.cxx:69
void AddRangeCut(TCutG *cut, Bool_t exclude=kFALSE)
TGo4FitModel * AddModel(TGo4FitModel *m)
Definition: TGo4Fitter.cxx:203
TCutG * MakeCut()
Definition: Example13.cxx:86
void Example13()
Definition: Example13.cxx:93
void SetParsValues(Double_t *pars)