00001 #include "TH2.h"
00002 #include "TH3.h"
00003 #include "TF2.h"
00004 #include "TF3.h"
00005 #include "TCanvas.h"
00006 #include "TApplication.h"
00007 #include "TMath.h"
00008
00009 #include "Fit/SparseData.h"
00010 #include "HFitInterface.h"
00011 #include "Fit/Fitter.h"
00012 #include "Math/WrappedMultiTF1.h"
00013
00014 #include "TRandom.h"
00015
00016 #include <iostream>
00017 #include <iterator>
00018 #include <algorithm>
00019
00020 #include <list>
00021 #include <vector>
00022
00023 #include <cmath>
00024 #include <cassert>
00025
00026 using namespace std;
00027
00028 const bool __DRAW__ = 1;
00029
00030 double gaus2D(double *x, double *p)
00031 {
00032 return p[0]*TMath::Gaus(x[0],p[1],p[2]) * TMath::Gaus(x[1],p[3],p[4]);
00033 }
00034
00035 double gaus3D(double *x, double *p)
00036 {
00037 return p[0] * TMath::Gaus(x[0],p[1],p[2])
00038 * TMath::Gaus(x[1],p[3],p[4])
00039 * TMath::Gaus(x[2],p[5],p[6]);
00040 }
00041
00042 double pol2D(double *x, double *p)
00043 {
00044 return p[0]*x[0]+ p[1]*x[0]*x[0] + p[2]*x[1] + p[3]*x[1]*x[1] + p[4];
00045 }
00046
00047 ostream& operator << (ostream& out, ROOT::Fit::BinData& bd)
00048 {
00049 const unsigned int ndim( bd.NDim() );
00050 const unsigned int npoints( bd.NPoints() );
00051 for ( unsigned int i = 0; i < npoints; ++i )
00052 {
00053 double value, error;
00054 const double *x = bd.GetPoint(i, value, error);
00055 for ( unsigned int j = 0; j < ndim; ++j )
00056 {
00057 out << " x[" << j << "]: " << x[j];
00058 }
00059 out << " value: " << value;
00060 out << " error: " << error;
00061 out << endl;
00062 }
00063 return out;
00064 }
00065
00066 int findBin(ROOT::Fit::BinData& bd, const double *x)
00067 {
00068 const unsigned int ndim = bd.NDim();
00069 const unsigned int npoints = bd.NPoints();
00070
00071 for ( unsigned int i = 0; i < npoints; ++i )
00072 {
00073 double value1, error1;
00074 const double *x1 = bd.GetPoint(i, value1, error1);
00075 bool thisIsIt = true;
00076 for ( unsigned int j = 0; j < ndim; ++j )
00077 {
00078 thisIsIt &= fabs(x1[j] - x[j]) < 1E-15;
00079 }
00080 if ( thisIsIt ) return i;
00081 }
00082
00083 cout << "NO ENCONTRADO!";
00084 copy(x, x+ndim, ostream_iterator<double>(cout, " " ));
00085 cout << endl;
00086
00087 return -1;
00088 }
00089
00090 bool operator ==(ROOT::Fit::BinData& bd1, ROOT::Fit::BinData& bd2)
00091 {
00092 const unsigned int ndim = bd1.NDim();
00093 const unsigned int npoints = bd1.NPoints();
00094
00095 bool equals = true;
00096
00097 cout << "Equals" << endl;
00098
00099 for ( unsigned int i = 0; i < npoints && equals; ++i )
00100 {
00101 double value1, error1;
00102 const double *x1 = bd1.GetPoint(i, value1, error1);
00103
00104 int bin = findBin(bd2, x1);
00105
00106 double value2 = 0, error2;
00107 const double *x2 = bd2.GetPoint(bin, value2, error2);
00108
00109 equals &= ( value1 == value2 );
00110 cout << " v: " << equals;
00111 equals &= ( error1 == error2 );
00112 cout << " e: " << equals;
00113 for ( unsigned int j = 0; j < ndim; ++j )
00114 {
00115 equals &= fabs(x1[j] - x2[j]) < 1E-15;
00116 cout << " x[" << j << "]: " << equals;
00117 }
00118
00119 cout << " bd1: ";
00120 std::copy(x1, &x1[ndim], ostream_iterator<double>(cout, " "));
00121 cout << " value:" << value1;
00122 cout << " error:" << error1;
00123
00124 cout << " bd2: ";
00125 std::copy(x2, &x2[ndim], ostream_iterator<double>(cout, " "));
00126 cout << " value:" << value2;
00127 cout << " error:" << error2;
00128
00129 cout << " equals: " << equals;
00130
00131 cout << endl;
00132 }
00133
00134 return equals;
00135 }
00136
00137 void fit3DHist()
00138 {
00139 vector<double> min(3); min[0] = 0.; min[1] = 0.; min[2] = 0.;
00140 vector<double> max(3); max[0] = 10.; max[1] = 10.; max[2] = 10.;
00141 vector<int> nbins(3); nbins[0] = 10; nbins[1] = 10; nbins[2] = 10;
00142
00143 TH3D* h1 = new TH3D("3D Original Hist Fit", "h1-title",
00144 nbins[0], min[0], max[0],
00145 nbins[1], min[1], max[1],
00146 nbins[2], min[2], max[2]);
00147 TH3D* h2 = new TH3D("3D Blanked Hist Fit", "h1-title",
00148 nbins[0], min[0], max[0],
00149 nbins[1], min[1], max[1],
00150 nbins[2], min[2], max[2]);
00151
00152 TF3* f1 = new TF3("func3D", gaus3D,
00153 min[0],max[0],
00154 min[1],max[1],
00155 min[2],max[2],
00156 7);
00157 double initialPars[] = {20,5,2,5,1,5,2};
00158
00159
00160
00161 f1->SetParameters(initialPars);
00162
00163
00164
00165
00166 for (int ix=1; ix <= h1->GetNbinsX(); ++ix) {
00167 for (int iy=1; iy <= h1->GetNbinsY(); ++iy) {
00168 for (int iz=1; iz <= h1->GetNbinsZ(); ++iz) {
00169 double x = h1->GetXaxis()->GetBinCenter(ix);
00170 double y = h1->GetYaxis()->GetBinCenter(iy);
00171 double z = h1->GetZaxis()->GetBinCenter(iz);
00172
00173 h1->SetBinContent( ix, iy, iz, gRandom->Poisson( f1->Eval(x,y,z) ) );
00174 }
00175 }
00176 }
00177
00178
00179 cout << "Retrieving the Sparse Data Structure" << endl;
00180 ROOT::Fit::SparseData d(min,max);
00181 ROOT::Fit::FillData(d, h1, 0);
00182 ROOT::Fit::BinData bd;
00183 d.GetBinData(bd);
00184
00185
00186 cout << "Filling second histogram" << endl;
00187 for ( unsigned int i = 0; i < bd.NPoints(); ++i )
00188 {
00189 const double* x;
00190 double value, error;
00191 x = bd.GetPoint(i, value, error);
00192 value = (value)?value:-10;
00193 h2->Fill(x[0], x[1], x[2], value);
00194 }
00195
00196
00197
00198
00199 bool ret;
00200 ROOT::Fit::Fitter fitter;
00201 ROOT::Math::WrappedMultiTF1 wf1(*f1);
00202 ROOT::Math::IParametricFunctionMultiDim & if1 = wf1;
00203 fitter.Config().SetMinimizer("Minuit2");
00204
00205
00206
00207 cout << "\n ******* Chi2Fit with Original BinData *******" << endl;
00208 ROOT::Fit::BinData bdOriginal;
00209 ROOT::Fit::FillData(bdOriginal, h1, 0);
00210 ret = fitter.LikelihoodFit(bdOriginal, if1);
00211 fitter.Result().Print(std::cout);
00212 if (!ret)
00213 std::cout << "Fit Failed " << std::endl;
00214
00215
00216 cout << "\n ******* Chi2Fit with Original BinData with Ceros*******" << endl;
00217 ROOT::Fit::DataOptions opt;
00218 opt.fUseEmpty = true;
00219 opt.fIntegral = true;
00220 ROOT::Fit::BinData bdOriginalWithCeros(opt);
00221 ROOT::Fit::FillData(bdOriginalWithCeros, h1, 0);
00222 ret = fitter.LikelihoodFit(bdOriginalWithCeros, if1);
00223 fitter.Result().Print(std::cout);
00224 if (!ret)
00225 std::cout << "Fit Failed " << std::endl;
00226
00227
00228 cout << "\n ******* Chi2Fit with BinData and NoCeros *******" << endl;
00229 ROOT::Fit::BinData bdNoCeros;
00230 d.GetBinDataNoZeros(bdNoCeros);
00231 ret = fitter.LikelihoodFit(bdNoCeros, if1);
00232 fitter.Result().Print(std::cout);
00233 if (!ret)
00234 std::cout << "Fit Failed " << std::endl;
00235
00236
00237
00238
00239
00240 cout << "\n ******* Chi2Fit with BinData with Ceros *******" << endl;
00241 ROOT::Fit::BinData bdWithCeros(opt);
00242 d.GetBinDataIntegral(bdWithCeros);
00243 ret = fitter.LikelihoodFit(bdWithCeros, if1);
00244 fitter.Result().Print(std::cout);
00245 if (!ret)
00246 std::cout << "Fit Failed " << std::endl;
00247
00248
00249
00250
00251
00252
00253
00254
00255 TCanvas* c = new TCanvas("Histogram 3D");
00256 c->Divide(1,2);
00257 c->cd(1);
00258 h1->Draw("lego");
00259 c->cd(2);
00260 h2->Draw("lego");
00261 }
00262
00263 void fit2DHist()
00264 {
00265 vector<double> min(2); min[0] = 0.; min[1] = 0.;
00266 vector<double> max(2); max[0] = 10.; max[1] = 10.;
00267 vector<int> nbins(2); nbins[0] = 10; nbins[1] = 10;
00268
00269 TH2D* h1 = new TH2D("2D Original Hist Fit", "h1-title", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
00270 TH2D* h2 = new TH2D("2D Blanked Hist Fit", "h1-title", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
00271
00272 TF2* f2 = new TF2("func2D", gaus2D,
00273 min[0],max[0], min[1], max[1], 5);
00274 double initialPars[] = {20,5,2,5,1};
00275
00276
00277
00278 f2->SetParameters(initialPars);
00279
00280
00281
00282
00283 for (int ix=1; ix <= h1->GetNbinsX(); ++ix) {
00284 for (int iy=1; iy <= h1->GetNbinsY(); ++iy) {
00285 double x= h1->GetXaxis()->GetBinCenter(ix);
00286 double y= h1->GetYaxis()->GetBinCenter(iy);
00287
00288 h1->SetBinContent( ix, iy, gRandom->Poisson( f2->Eval(x,y) ) );
00289 }
00290 }
00291
00292
00293 cout << "Retrieving the Sparse Data Structure" << endl;
00294 ROOT::Fit::SparseData d(min,max);
00295 ROOT::Fit::FillData(d, h1, 0);
00296 ROOT::Fit::BinData bd;
00297 d.GetBinData(bd);
00298
00299
00300 cout << "Filling second histogram" << endl;
00301 for ( unsigned int i = 0; i < bd.NPoints(); ++i )
00302 {
00303 const double* x;
00304 double value, error;
00305 x = bd.GetPoint(i, value, error);
00306 value = (value)?value:-10;
00307 h2->Fill(x[0], x[1], value);
00308 }
00309
00310
00311
00312 bool ret;
00313 ROOT::Fit::Fitter fitter;
00314 ROOT::Math::WrappedMultiTF1 wf2(*f2);
00315 ROOT::Math::IParametricFunctionMultiDim & if2 = wf2;
00316 fitter.Config().SetMinimizer("Minuit2");
00317
00318
00319
00320 cout << "\n ******* Chi2Fit with Original BinData *******" << endl;
00321 ROOT::Fit::BinData bdOriginal;
00322 ROOT::Fit::FillData(bdOriginal, h1, 0);
00323 ret = fitter.LikelihoodFit(bdOriginal, if2);
00324 fitter.Result().Print(std::cout);
00325 if (!ret)
00326 std::cout << "Fit Failed " << std::endl;
00327
00328
00329 cout << "\n ******* Chi2Fit with Original BinData with Ceros*******" << endl;
00330 ROOT::Fit::DataOptions opt;
00331 opt.fUseEmpty = true;
00332 opt.fIntegral = true;
00333 ROOT::Fit::BinData bdOriginalWithCeros(opt);
00334 ROOT::Fit::FillData(bdOriginalWithCeros, h1, 0);
00335 ret = fitter.LikelihoodFit(bdOriginalWithCeros, if2);
00336 fitter.Result().Print(std::cout);
00337 if (!ret)
00338 std::cout << "Fit Failed " << std::endl;
00339
00340
00341 cout << "\n ******* Chi2Fit with BinData and NoCeros *******" << endl;
00342 ROOT::Fit::BinData bdNoCeros;
00343 d.GetBinDataNoZeros(bdNoCeros);
00344 ret = fitter.LikelihoodFit(bdNoCeros, if2);
00345 fitter.Result().Print(std::cout);
00346 if (!ret)
00347 std::cout << "Fit Failed " << std::endl;
00348
00349
00350
00351
00352
00353 cout << "\n ******* Chi2Fit with BinData with Ceros *******" << endl;
00354 ROOT::Fit::BinData bdWithCeros(opt);
00355 d.GetBinDataIntegral(bdWithCeros);
00356 ret = fitter.LikelihoodFit(bdWithCeros, if2);
00357 fitter.Result().Print(std::cout);
00358 if (!ret)
00359 std::cout << "Fit Failed " << std::endl;
00360
00361
00362
00363
00364
00365
00366
00367
00368 TCanvas* c = new TCanvas("Histogram 2D");
00369 c->Divide(2,2);
00370 c->cd(1);
00371 h1->Draw("colz");
00372 c->cd(2);
00373 h1->Draw("text");
00374 c->cd(3);
00375 h2->Draw("colz");
00376 c->cd(4);
00377 h2->Draw("text");
00378 }
00379
00380
00381
00382
00383
00384
00385
00386
00387
00388
00389
00390
00391
00392
00393
00394
00395
00396
00397
00398
00399
00400
00401
00402
00403
00404
00405
00406
00407
00408
00409
00410
00411
00412
00413
00414
00415
00416
00417
00418
00419
00420
00421
00422
00423
00424
00425
00426
00427
00428
00429
00430
00431
00432
00433
00434
00435
00436
00437
00438
00439
00440
00441
00442
00443
00444
00445
00446
00447
00448
00449
00450
00451
00452 int main(int argc, char** argv)
00453 {
00454 TApplication* theApp = 0;
00455
00456 if ( __DRAW__ )
00457 theApp = new TApplication("App",&argc,argv);
00458
00459 fit3DHist();
00460
00461
00462
00463 if ( __DRAW__ ) {
00464 theApp->Run();
00465 delete theApp;
00466 theApp = 0;
00467 }
00468
00469 return 0;
00470 }
00471