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