00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019 #ifndef __CINT__
00020
00021 #include <iostream.h>
00022
00023 #include "TH1.h"
00024 #include "TH2.h"
00025 #include "TCanvas.h"
00026 #include "TApplication.h"
00027 #include "TRandom.h"
00028
00029 #include "../Go4Fit/TGo4Fitter.h"
00030 #include "../Go4Fit/TGo4FitDataHistogram.h"
00031 #include "../Go4Fit/TGo4FitDataRidge.h"
00032 #include "../Go4Fit/TGo4FitModelFormula.h"
00033
00034 void Example13();
00035
00036 int main(int argc, char **argv) {
00037
00038 TApplication theApp("Application", 0, 0);
00039
00040 Example13();
00041
00042 theApp.Run();
00043
00044 return 0;
00045 }
00046
00047 #endif
00048
00049 void AddRidge(TH2D* histo, double a, double k, double b,
00050 double xstart, int numsteps, int numpoints, double deviation) {
00051 TRandom rnd(1234);
00052
00053 double x = xstart;
00054
00055 for(int i=0;i<numsteps;i++) {
00056 double y = a*exp(k*x) + b;
00057
00058 double alf = TMath::ATan(k*a*exp(k*x)) + TMath::Pi()/2.;
00059
00060 for(int j=0;j<numpoints;j++) {
00061 double radius = rnd.Gaus(0.,deviation);
00062 histo->Fill(x + radius*cos(alf), y + radius*sin(alf));
00063 }
00064
00065 x = x + 0.001*sin(alf);
00066 }
00067 }
00068
00069 TH2D* MakeRidgeHistogram() {
00070 TH2D* histo = new TH2D("Histo","Ridge y = a*exp(k*x)+b",100,0.,10.,100,0.,10.);
00071
00072 cout << "Generating histogram " << endl;
00073
00074 AddRidge(histo, 10, -0.3, 2, 1.5, 10000, 100, 0.3);
00075
00076 AddRidge(histo, 10, -0.5, 0, 0.5, 10000, 50, 0.2);
00077
00078 AddRidge(histo, 10, -0.1, +2, 3.5, 10000, 50, 0.2);
00079
00080 cout << " Done " << endl;
00081
00082 return histo;
00083 }
00084
00085 TCutG* MakeCut() {
00086 Float_t xx[7] = { 1.0, 4.0, 9.0, 9.5, 5.0, 2.0, 1.0 };
00087 Float_t yy[7] = { 8.0, 2.5, 2.0, 4.0, 6.0, 9.5, 8.0 };
00088 return new TCutG("RangeCut",7,xx,yy);
00089 }
00090
00091 void Example13() {
00092
00093 TH2D* histo = MakeRidgeHistogram();
00094 new TCanvas("Canvas1");
00095 histo->Draw("LEGO");
00096
00097
00098 TGo4Fitter fitter("Fitter", TGo4Fitter::ff_chi_square, kFALSE);
00099
00100
00101
00102
00103 TGo4FitDataHistogram* datah = new TGo4FitDataHistogram("histogram",histo);
00104 datah->SetExcludeLessThen(3);
00105 datah->AddRangeCut(MakeCut());
00106
00107
00108 TGo4FitDataRidge* ridge = new TGo4FitDataRidge("ridge",datah,1);
00109 fitter.AddData(ridge);
00110
00111
00112
00113 TGo4FitModelFormula* model = new TGo4FitModelFormula("RidgeFunc","parA*exp(parK*x)+parB",3,kFALSE);
00114 model->SetRange(0,2,8.);
00115 model->SetParsNames("parA","parK","parB");
00116 model->SetParsValues(9.,-0.2,1.2);
00117 fitter.AddModel("ridge", model);
00118
00119
00120 fitter.AddSimpleMinuit();
00121
00122
00123 fitter.SetMemoryUsage(1);
00124
00125
00126 fitter.DoActions();
00127
00128
00129 fitter.Print("Pars");
00130
00131
00132 fitter.Draw("#ridge");
00133 }
00134
00135