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