GSI Object Oriented Online Offline (Go4)  GO4-6.3.0
Example13.cxx
Go to the documentation of this file.
1 // $Id$
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 fuer 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 <iostream>
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", nullptr, nullptr);
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, double xstart, int numsteps, int numpoints, double deviation)
49 {
50  TRandom rnd(1234);
51 
52  double x = xstart;
53 
54  for (int i = 0; i < numsteps; i++) {
55  double y = a * exp(k * x) + b;
56 
57  double alf = TMath::ATan(k * a * exp(k * x)) + TMath::Pi() / 2.;
58 
59  for (int j = 0; j < numpoints; j++) {
60  double radius = rnd.Gaus(0., deviation);
61  histo->Fill(x + radius * cos(alf), y + radius * sin(alf));
62  }
63 
64  x = x + 0.001 * sin(alf);
65  }
66 }
67 
69 {
70  TH2D *histo = new TH2D("Histo", "Ridge y = a*exp(k*x)+b", 100, 0., 10., 100, 0., 10.);
71 
72  std::cout << "Generating histogram " << std::endl;
73 
74  AddRidge(histo, 10, -0.3, 2, 1.5, 10000, 100, 0.3);
75 
76  AddRidge(histo, 10, -0.5, 0, 0.5, 10000, 50, 0.2);
77 
78  AddRidge(histo, 10, -0.1, +2, 3.5, 10000, 50, 0.2);
79 
80  std::cout << " Done " << std::endl;
81 
82  return histo;
83 }
84 
85 TCutG *MakeCut()
86 {
87  Float_t xx[7] = {1.0, 4.0, 9.0, 9.5, 5.0, 2.0, 1.0};
88  Float_t yy[7] = {8.0, 2.5, 2.0, 4.0, 6.0, 9.5, 8.0};
89  return new TCutG("RangeCut", 7, xx, yy);
90 }
91 
92 void Example13()
93 {
94  // model histogram with ridge with exponential function
95  TH2D *histo = MakeRidgeHistogram();
96  new TCanvas("Canvas1");
97  histo->Draw("LEGO2");
98 
99  // create fitter, select fit function and not add standard actions list
100  TGo4Fitter fitter("Fitter", TGo4Fitter::ff_chi_square, kFALSE);
101 
102  // create data object to represent of two-dimensional histogram
103  // bins, which are less then 3 are excluded
104  // select middle ridge by TCutG object
105  TGo4FitDataHistogram *datah = new TGo4FitDataHistogram("histogram", histo);
106  datah->SetExcludeLessThen(3);
107  datah->AddRangeCut(MakeCut());
108 
109  // create ridge object, which select axis number 1 (y axis) as function of the rest coordinates
110  TGo4FitDataRidge *ridge = new TGo4FitDataRidge("ridge", datah, 1);
111  fitter.AddData(ridge);
112 
113  // create formula function to approximate ridge function
114  // specific range and initial values for parameters need to be set
115  TGo4FitModelFormula *model = new TGo4FitModelFormula("RidgeFunc", "parA*exp(parK*x)+parB", 3, kFALSE);
116  model->SetRange(0, 2, 8.);
117  model->SetParsNames("parA", "parK", "parB");
118  model->SetParsValues(9., -0.2, 1.2);
119  fitter.AddModel("ridge", model);
120 
121  // add minuit action
122  fitter.AddSimpleMinuit();
123 
124  // use memory buffers for data objects
125  fitter.SetMemoryUsage(1);
126 
127  // execute actions
128  fitter.DoActions();
129 
130  // print parameters values
131  fitter.Print("Pars");
132 
133  // draw results of fit
134  fitter.Draw("#ridge");
135 }
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")
void SetMemoryUsage(Int_t iMemoryUsage)
Definition: TGo4Fitter.cxx:70
void SetExcludeLessThen(Double_t limit=0.)
Definition: TGo4FitData.h:117
int main(int argc, char **argv)
Definition: Example13.cxx:35
TGo4FitData * AddData(TGo4FitData *d)
Definition: TGo4Fitter.cxx:119
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)
TH2D * MakeRidgeHistogram()
Definition: Example13.cxx:68
void AddRangeCut(TCutG *cut, Bool_t exclude=kFALSE)
TGo4FitModel * AddModel(TGo4FitModel *m)
Definition: TGo4Fitter.cxx:210
TCutG * MakeCut()
Definition: Example13.cxx:85
void Example13()
Definition: Example13.cxx:92
void DoActions(Bool_t AllowFitterChange=kFALSE, TObjArray *Actions=nullptr)
void Draw(Option_t *option) override
void SetParsValues(Double_t *pars)
void Print(Option_t *option="") const override
Definition: TGo4Fitter.cxx:779