00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014 #ifndef __CINT__
00015
00016 #include "Riostream.h"
00017 #include "TH1.h"
00018 #include "TH2.h"
00019 #include "TFile.h"
00020 #include "TCanvas.h"
00021 #include "TApplication.h"
00022
00023 #include "TGo4FitMinuit.h"
00024 #include "TGo4Fitter.h"
00025 #include "TGo4FitDataHistogram.h"
00026 #include "TGo4FitModelPolynom.h"
00027 #include "TGo4FitModelGauss2.h"
00028 #include "TGo4FitModelGaussN.h"
00029
00030
00031 void Example9();
00032
00033 int main(int argc, char **argv)
00034 {
00035 TApplication theApp("Application", 0, 0);
00036
00037 Example9();
00038
00039 theApp.Run();
00040
00041 return 0;
00042 }
00043
00044 #endif
00045
00046 void DrawHistogram(TH1* histo, const char* CanvasName, const char* DrawOption)
00047 {
00048 TCanvas *fCanvas = new TCanvas(CanvasName,"Draw of histogram",3);
00049 fCanvas->cd();
00050 histo->Draw(DrawOption);
00051 fCanvas->Update();
00052 }
00053
00054 void Example9()
00055 {
00056
00057 TGo4Fitter *fitter = new TGo4Fitter("Fitter","Example fitter object");
00058
00059
00060 TH2D *histo = new TH2D("histo","dummy histogram",1000,0.,10.,1000,0.,10.);
00061 TGo4FitData* data = fitter->AddData(new TGo4FitDataHistogram("data",histo,kTRUE));
00062 data->SetRange(0,1.,3.); data->SetRange(0,7.,9.);
00063 data->SetRange(1,1.,3.); data->SetRange(1,7.,9.);
00064
00065
00066 TGo4FitModel* model1 = fitter->AddModel( "data", new TGo4FitModelGauss2("Gauss1",2.,8.,.5,.5,0.) );
00067 model1->SetAmplValue(1000.);
00068
00069 TGo4FitModel* model2 = fitter->AddModel( "data", new TGo4FitModelGauss2("Gauss2",8.,8.,.5,.5,0.) );
00070 model2->SetAmplValue(1000.);
00071
00072 TGo4FitModel* model3 = fitter->AddModel( "data", new TGo4FitModelGauss2("Gauss3",8.,2.,.5,.5,0.) );
00073 model3->SetAmplValue(1000.);
00074
00075 TGo4FitModel* model4 = fitter->AddModel( "data", new TGo4FitModelGauss2("Gauss4",2.,2.,.5,.5,0.) );
00076 model4->SetAmplValue(1000.);
00077
00078
00079
00080 TH1* res1 = (TH1*) fitter->CreateDrawObject("Large", "data", kTRUE);
00081
00082 histo = new TH2D("histo2","dummy histogram",10,0.,10.,10,0.,10.);
00083 fitter->SetObject("data", histo, kTRUE);
00084
00085 TH1* res2 = (TH1*) fitter->CreateDrawObject("Small", "data", kTRUE);
00086
00087 model1->SetIntegrationsProperty(5);
00088 model2->SetIntegrationsProperty(5);
00089 model3->SetIntegrationsProperty(5);
00090 model4->SetIntegrationsProperty(5);
00091
00092 TH1* res3 = (TH1*) fitter->CreateDrawObject("SmallI", "data", kTRUE);
00093
00094 delete fitter;
00095
00096
00097
00098
00099 DrawHistogram(res2,"Canvas2","SURF1");
00100 DrawHistogram(res3,"Canvas3","SURF1");
00101
00102 Double_t i1 = res1->Integral()/1000000.;
00103 Double_t i2 = res2->Integral()/100.;
00104 Double_t i3 = res3->Integral()/100.;
00105
00106 std::cout << "Integral over 1000x1000 points = " << i1 << std::endl;
00107 std::cout << "Integral over 10x10 points = " << i2 << std::endl;
00108 std::cout << "Integral over 10x10 with model integr. = " << i3 << std::endl << std::endl;
00109
00110 std::cout << "Integral2/Integral1 = " << i2/i1 << std::endl;
00111 std::cout << "Integral3/Integral1 = " << i3/i1 << std::endl;
00112 }