00001
00002
00003
00004 #include "TH2D.h"
00005 #include "TF2.h"
00006 #include "TCanvas.h"
00007 #include "TStyle.h"
00008 #include "TRandom3.h"
00009 #include "TVirtualFitter.h"
00010 #include "TList.h"
00011
00012 #include <vector>
00013 #include <map>
00014 #include <iostream>
00015
00016 double gauss2D(double *x, double *par) {
00017 double z1 = double((x[0]-par[1])/par[2]);
00018 double z2 = double((x[1]-par[3])/par[4]);
00019 return par[0]*exp(-0.5*(z1*z1+z2*z2));
00020 }
00021 double my2Dfunc(double *x, double *par) {
00022 double *p1 = &par[0];
00023 double *p2 = &par[5];
00024 return gauss2D(x,p1) + gauss2D(x,p2);
00025 }
00026
00027
00028
00029
00030
00031 std::vector<std::pair<double, double> > coords;
00032 std::vector<double > values;
00033 std::vector<double > errors;
00034
00035 void myFcn(Int_t & , Double_t * , Double_t &fval, Double_t *p, Int_t )
00036 {
00037 int n = coords.size();
00038 double chi2 = 0;
00039 double tmp,x[2];
00040 for (int i = 0; i <n; ++i ) {
00041 x[0] = coords[i].first;
00042 x[1] = coords[i].second;
00043 tmp = ( values[i] - my2Dfunc(x,p))/errors[i];
00044 chi2 += tmp*tmp;
00045 }
00046 fval = chi2;
00047 }
00048 TRandom3 rndm;
00049 void FillHisto(TH2D * h, int n, double * p) {
00050
00051
00052 const double mx1 = p[1];
00053 const double my1 = p[3];
00054 const double sx1 = p[2];
00055 const double sy1 = p[4];
00056 const double mx2 = p[6];
00057 const double my2 = p[8];
00058 const double sx2 = p[7];
00059 const double sy2 = p[9];
00060
00061 const double w1 = 0.5;
00062
00063 double x, y;
00064 for (int i = 0; i < n; ++i) {
00065
00066 rndm.Rannor(x,y);
00067
00068 double r = rndm.Rndm(1);
00069 if (r < w1) {
00070 x = x*sx1 + mx1;
00071 y = y*sy1 + my1;
00072 }
00073 else {
00074 x = x*sx2 + mx2;
00075 y = y*sy2 + my2;
00076 }
00077 h->Fill(x,y);
00078
00079 }
00080 }
00081
00082
00083
00084
00085 int TwoHistoFit2D(bool global = true) {
00086
00087
00088
00089 int nbx1 = 50;
00090 int nby1 = 50;
00091 int nbx2 = 50;
00092 int nby2 = 50;
00093 double xlow1 = 0.;
00094 double ylow1 = 0.;
00095 double xup1 = 10.;
00096 double yup1 = 10.;
00097 double xlow2 = 5.;
00098 double ylow2 = 5.;
00099 double xup2 = 20.;
00100 double yup2 = 20.;
00101
00102 TH2D * h1 = new TH2D("h1","core",nbx1,xlow1,xup1,nby1,ylow1,yup1);
00103 TH2D * h2 = new TH2D("h2","tails",nbx2,xlow2,xup2,nby2,ylow2,yup2);
00104
00105 double iniParams[10] = { 100, 6., 2., 7., 3, 100, 12., 3., 11., 2. };
00106
00107 TF2 * func = new TF2("func",my2Dfunc,xlow2,xup2,ylow2,yup2, 10);
00108 func->SetParameters(iniParams);
00109
00110
00111 int n1 = 1000000;
00112 int n2 = 1000000;
00113
00114
00115 FillHisto(h1,n1,iniParams);
00116 FillHisto(h2,n2,iniParams);
00117
00118
00119 double dx1 = (xup1-xlow1)/double(nbx1);
00120 double dy1 = (yup1-ylow1)/double(nby1);
00121 double dx2 = (xup2-xlow2)/double(nbx2);
00122 double dy2 = (yup2-ylow2)/double(nby2);
00123
00124
00125
00126 h2->Sumw2();
00127 h2->Scale( ( double(n1) * dx1 * dy1 ) / ( double(n2) * dx2 * dy2 ) );
00128
00129
00130 if (global) {
00131
00132 std::cout << "Do global fit" << std::endl;
00133
00134
00135
00136 TAxis *xaxis1 = h1->GetXaxis();
00137 TAxis *yaxis1 = h1->GetYaxis();
00138 TAxis *xaxis2 = h2->GetXaxis();
00139 TAxis *yaxis2 = h2->GetYaxis();
00140
00141 int nbinX1 = h1->GetNbinsX();
00142 int nbinY1 = h1->GetNbinsY();
00143 int nbinX2 = h2->GetNbinsX();
00144 int nbinY2 = h2->GetNbinsY();
00145
00146
00147 coords = std::vector<std::pair<double,double> >();
00148 values = std::vector<double>();
00149 errors = std::vector<double>();
00150
00151
00152 for (int ix = 1; ix <= nbinX1; ++ix) {
00153 for (int iy = 1; iy <= nbinY1; ++iy) {
00154 if ( h1->GetBinContent(ix,iy) > 0 ) {
00155 coords.push_back( std::make_pair(xaxis1->GetBinCenter(ix), yaxis1->GetBinCenter(iy) ) );
00156 values.push_back( h1->GetBinContent(ix,iy) );
00157 errors.push_back( h1->GetBinError(ix,iy) );
00158 }
00159 }
00160 }
00161 for (int ix = 1; ix <= nbinX2; ++ix) {
00162 for (int iy = 1; iy <= nbinY2; ++iy) {
00163 if ( h2->GetBinContent(ix,iy) > 0 ) {
00164 coords.push_back( std::make_pair(xaxis2->GetBinCenter(ix), yaxis2->GetBinCenter(iy) ) );
00165 values.push_back( h2->GetBinContent(ix,iy) );
00166 errors.push_back( h2->GetBinError(ix,iy) );
00167 }
00168 }
00169 }
00170
00171 TVirtualFitter::SetDefaultFitter("Minuit");
00172 TVirtualFitter * minuit = TVirtualFitter::Fitter(0,10);
00173 for (int i = 0; i < 10; ++i) {
00174 minuit->SetParameter(i, func->GetParName(i), func->GetParameter(i), 0.01, 0,0);
00175 }
00176 minuit->SetFCN(myFcn);
00177
00178 double arglist[100];
00179 arglist[0] = 0;
00180
00181 minuit->ExecuteCommand("SET PRINT",arglist,2);
00182
00183
00184 arglist[0] = 5000;
00185 arglist[1] = 0.01;
00186 minuit->ExecuteCommand("MIGRAD",arglist,2);
00187
00188
00189 double minParams[10];
00190 double parErrors[10];
00191 for (int i = 0; i < 10; ++i) {
00192 minParams[i] = minuit->GetParameter(i);
00193 parErrors[i] = minuit->GetParError(i);
00194 }
00195 double chi2, edm, errdef;
00196 int nvpar, nparx;
00197 minuit->GetStats(chi2,edm,errdef,nvpar,nparx);
00198
00199 func->SetParameters(minParams);
00200 func->SetParErrors(parErrors);
00201 func->SetChisquare(chi2);
00202 int ndf = coords.size()-nvpar;
00203 func->SetNDF(ndf);
00204
00205 std::cout << "Chi2 Fit = " << chi2 << " ndf = " << ndf << " " << func->GetNDF() << std::endl;
00206
00207
00208 h1->GetListOfFunctions()->Add(func);
00209 h2->GetListOfFunctions()->Add(func);
00210 }
00211 else {
00212
00213 h1->Fit(func);
00214 h2->Fit(func);
00215 }
00216
00217
00218
00219
00220 TCanvas * c1 = new TCanvas("c1","Two HIstogram Fit example",100,10,900,800);
00221 c1->Divide(2,2);
00222 gStyle->SetOptFit();
00223 gStyle->SetStatY(0.6);
00224
00225 c1->cd(1);
00226 h1->Draw();
00227 func->SetRange(xlow1,ylow1,xup1,yup1);
00228 func->DrawCopy("cont1 same");
00229 c1->cd(2);
00230 h1->Draw("lego");
00231 func->DrawCopy("surf1 same");
00232 c1->cd(3);
00233 func->SetRange(xlow2,ylow2,xup2,yup2);
00234 h2->Draw();
00235 func->DrawCopy("cont1 same");
00236 c1->cd(4);
00237 h2->Draw("lego");
00238 gPad->SetLogz();
00239 func->Draw("surf1 same");
00240
00241 return 0;
00242 }