00001 #include "TH1.h"
00002 #include "TH2.h"
00003 #include "TF1.h"
00004 #include "TF2.h"
00005 #include "TGraphErrors.h"
00006 #include "TGraphAsymmErrors.h"
00007 #include "TGraph2DErrors.h"
00008 #include "TSystem.h"
00009 #include "TRandom3.h"
00010 #include "TROOT.h"
00011 #include "TVirtualFitter.h"
00012
00013 #include "Fit/BinData.h"
00014 #include "Fit/UnBinData.h"
00015 #include "HFitInterface.h"
00016 #include "Fit/Fitter.h"
00017
00018 #include "Math/WrappedMultiTF1.h"
00019 #include "Math/WrappedParamFunction.h"
00020 #include "Math/WrappedTF1.h"
00021
00022 #include "RConfigure.h"
00023
00024 #include <string>
00025 #include <iostream>
00026 #include <cmath>
00027
00028
00029 void printData(const ROOT::Fit::BinData & data) {
00030 std::cout << "Bin data, point size is " << data.PointSize() << " data dimension is " << data.NDim() << " type is " << data.GetErrorType() << std::endl;
00031 for (unsigned int i = 0; i < data.Size(); ++i) {
00032 if (data.GetErrorType() == ROOT::Fit::BinData::kNoError)
00033 std::cout << data.Coords(i)[0] << " " << data.Value(i) << std::endl;
00034 else if (data.GetErrorType() == ROOT::Fit::BinData::kValueError)
00035 std::cout << data.Coords(i)[0] << " " << data.Value(i) << " +/- " << data.Error(i) << std::endl;
00036 else if (data.GetErrorType() == ROOT::Fit::BinData::kCoordError) {
00037 double ey = 0;
00038 const double * ex = data.GetPointError(i,ey);
00039 std::cout << data.Coords(i)[0] << " +/- " << ex[0] << " " << data.Value(i) << " +/- " << ey << std::endl;
00040 }
00041 else if (data.GetErrorType() == ROOT::Fit::BinData::kAsymError) {
00042 double eyl = 0; double eyh = 0;
00043 const double * ex = data.GetPointError(i,eyl,eyh);
00044 std::cout << data.Coords(i)[0] << " +/- " << ex[0] << " " << data.Value(i) << " - " << eyl << " + " << eyh << std::endl;
00045 }
00046
00047 }
00048 std::cout << "\ndata size is " << data.Size() << std::endl;
00049 }
00050 void printData(const ROOT::Fit::UnBinData & data) {
00051 for (unsigned int i = 0; i < data.Size(); ++i) {
00052 std::cout << data.Coords(i)[0] << "\t";
00053 }
00054 std::cout << "\ndata size is " << data.Size() << std::endl;
00055 }
00056
00057 int compareResult(double v1, double v2, std::string s = "", double tol = 0.01) {
00058
00059
00060 if (std::abs(v1-v2) < tol * std::abs(v2) ) return 0;
00061 std::cerr << s << " Failed comparison of fit results \t chi2 = " << v1 << " it should be = " << v2 << std::endl;
00062 return -1;
00063 }
00064
00065 double chi2FromFit(const TF1 * func ) {
00066
00067 assert(TVirtualFitter::GetFitter() != 0 );
00068 return (TVirtualFitter::GetFitter()->Chisquare(func->GetNpar(), func->GetParameters() ) );
00069 }
00070
00071 int testHisto1DFit() {
00072
00073
00074 std::string fname("gaus");
00075 TF1 * func = (TF1*)gROOT->GetFunction(fname.c_str());
00076 func->SetParameter(0,10);
00077 func->SetParameter(1,0);
00078 func->SetParameter(2,3.0);
00079
00080 TRandom3 rndm;
00081 int iret = 0;
00082 double chi2ref = 0;
00083
00084
00085 TH1D * h1 = new TH1D("h1","h1",30,-5.,5.);
00086
00087 for (int i = 0; i <1000; ++i)
00088 h1->Fill( rndm.Gaus(0,1) );
00089
00090 h1->Print();
00091
00092
00093
00094
00095
00096
00097
00098
00099 ROOT::Fit::BinData d;
00100 ROOT::Fit::FillData(d,h1,func);
00101
00102
00103 printData(d);
00104
00105
00106 ROOT::Math::WrappedMultiTF1 wf(*func);
00107
00108 ROOT::Math::IParamMultiFunction & f = wf;
00109
00110 double p[3] = {100,0,3.};
00111 f.SetParameters(p);
00112
00113
00114
00115 ROOT::Fit::Fitter fitter;
00116
00117 bool ret = fitter.Fit(d, f);
00118 if (ret)
00119 fitter.Result().Print(std::cout);
00120 else {
00121 std::cout << "Chi2 Fit Failed " << std::endl;
00122 return -1;
00123 }
00124 chi2ref = fitter.Result().Chi2();
00125
00126
00127 TVirtualFitter::SetDefaultFitter("Minuit2");
00128 std::cout << "\n******************************\n\t TH1::Fit Result \n" << std::endl;
00129 func->SetParameters(p);
00130 h1->Fit(func);
00131
00132 iret |= compareResult( chi2FromFit(func), chi2ref,"1D histogram chi2 fit");
00133
00134
00135
00136 std::cout << "\n\nTest Binned Likelihood Fit" << std::endl;
00137
00138 ROOT::Fit::DataOptions opt;
00139 opt.fUseEmpty = true;
00140 ROOT::Fit::BinData dl(opt);
00141 ROOT::Fit::FillData(dl,h1,func);
00142
00143 ret = fitter.LikelihoodFit(dl, f);
00144 f.SetParameters(p);
00145 if (ret)
00146 fitter.Result().Print(std::cout);
00147 else {
00148 std::cout << "Binned Likelihood Fit Failed " << std::endl;
00149 return -1;
00150 }
00151 iret |= compareResult(fitter.Result().Chi2(), chi2ref,"1D histogram likelihood fit",0.3);
00152
00153
00154 std::cout << "\n******************************\n\t TH1::Fit Result \n" << std::endl;
00155 func->SetParameters(p);
00156 h1->Fit(func,"L");
00157
00158 iret |= compareResult(func->GetChisquare(),fitter.Result().Chi2(),"TH1::Fit likelihood ",0.001);
00159
00160
00161 std::cout << "\n\nTest Chi2 Fit using integral option" << std::endl;
00162
00163
00164 opt = ROOT::Fit::DataOptions();
00165 opt.fIntegral = true;
00166 ROOT::Fit::BinData d2(opt);
00167 ROOT::Fit::FillData(d2,h1,func);
00168
00169 f.SetParameters(p);
00170 ret = fitter.Fit(d2, f);
00171 if (ret)
00172 fitter.Result().Print(std::cout);
00173 else {
00174 std::cout << "Integral Chi2 Fit Failed " << std::endl;
00175 return -1;
00176 }
00177 iret |= compareResult(fitter.Result().Chi2(), chi2ref,"1D histogram integral chi2 fit",0.2);
00178
00179
00180 std::cout << "\n******************************\n\t TH1::Fit Result \n" << std::endl;
00181 func->SetParameters(p);
00182 h1->Fit(func,"I");
00183 iret |= compareResult(func->GetChisquare(),fitter.Result().Chi2(),"TH1::Fit integral ",0.001);
00184
00185
00186 opt = ROOT::Fit::DataOptions();
00187 opt.fIntegral = true;
00188 opt.fUseEmpty = true;
00189 ROOT::Fit::BinData dl2(opt);
00190 ROOT::Fit::FillData(dl2,h1,func);
00191 f.SetParameters(p);
00192 ret = fitter.LikelihoodFit(dl2, f);
00193 if (ret)
00194 fitter.Result().Print(std::cout);
00195 else {
00196 std::cout << "Integral Likelihood Fit Failed " << std::endl;
00197 return -1;
00198 }
00199 iret |= compareResult(fitter.Result().Chi2(), chi2ref,"1D histogram integral likelihood fit",0.3);
00200
00201
00202 std::cout << "\n******************************\n\t TH1::Fit Result \n" << std::endl;
00203 func->SetParameters(p);
00204 h1->Fit(func,"IL");
00205
00206 iret |= compareResult(func->GetChisquare(),fitter.Result().Chi2(),"TH1::Fit likelihood integral ",0.001);
00207
00208
00209
00210
00211 std::cout << "\n\nRedo Chi2 Hist Fit" << std::endl;
00212 f.SetParameters(p);
00213 ret = fitter.Fit(d, f);
00214 if (ret)
00215 fitter.Result().Print(std::cout);
00216 else {
00217 std::cout << "Chi2 Fit Failed " << std::endl;
00218 return -1;
00219 }
00220 iret |= compareResult(fitter.Result().Chi2(), chi2ref,"1D histogram chi2 fit (2)",0.001);
00221
00222
00223
00224
00225 std::cout << "\n\nTest Same Fit from a TGraphErrors - no coord errors" << std::endl;
00226 TGraphErrors gr(h1);
00227 ROOT::Fit::BinData dg;
00228 dg.Opt().fCoordErrors = false;
00229 ROOT::Fit::FillData(dg,&gr);
00230
00231 f.SetParameters(p);
00232 ret = fitter.Fit(dg, f);
00233 if (ret)
00234 fitter.Result().Print(std::cout);
00235 else {
00236 std::cout << "Chi2 Graph Errors Fit Failed " << std::endl;
00237 return -1;
00238 }
00239 iret |= compareResult(fitter.Result().Chi2(), chi2ref,"TGraphErrors chi2 fit",0.001);
00240
00241
00242
00243 std::cout << "\n\nTest Same Fit from a TGraphErrors - use coord errors" << std::endl;
00244 ROOT::Fit::BinData dger;
00245
00246
00247 dger.Opt().fUseEmpty = true;
00248 ROOT::Fit::FillData(dger,&gr);
00249
00250 f.SetParameters(p);
00251 ret = fitter.Fit(dger, f);
00252 if (ret)
00253 fitter.Result().Print(std::cout);
00254 else {
00255 std::cout << "Chi2 Graph Errors Fit Failed " << std::endl;
00256 return -1;
00257 }
00258 iret |= compareResult(fitter.Result().Chi2(), chi2ref,"TGraphErrors effective chi2 fit ",0.7);
00259
00260
00261 std::cout << "\n******************************\n\t TGraphErrors::Fit Result \n" << std::endl;
00262 func->SetParameters(p);
00263
00264 for (int ip = 0; ip < gr.GetN(); ++ip)
00265 if (gr.GetErrorY(ip) <= 0) gr.SetPointError(ip, gr.GetErrorX(ip), 1.);
00266
00267 gr.Fit(func);
00268 std::cout << "Ndf of TGraphErrors::Fit = " << func->GetNDF() << std::endl;
00269 iret |= compareResult(func->GetChisquare(),fitter.Result().Chi2(),"TGraphErrors::Fit ",0.001);
00270
00271
00272
00273 std::cout << "\n\nTest Same Fit from a TGraph" << std::endl;
00274 fitter.Config().SetNormErrors(true);
00275 TGraph gr2(h1);
00276 ROOT::Fit::BinData dg2;
00277 ROOT::Fit::FillData(dg2,&gr2);
00278
00279 f.SetParameters(p);
00280 ret = fitter.Fit(dg2, f);
00281 if (ret)
00282 fitter.Result().Print(std::cout);
00283 else {
00284 std::cout << "Chi2 Graph Fit Failed " << std::endl;
00285 return -1;
00286 }
00287
00288
00289
00290
00291 std::cout << "\n******************************\n\t TGraph::Fit Result \n" << std::endl;
00292 func->SetParameters(p);
00293 gr2.Fit(func);
00294 std::cout << "Ndf of TGraph::Fit = " << func->GetNDF() << std::endl;
00295
00296 iret |= compareResult(func->GetChisquare(),fitter.Result().Chi2(),"TGraph::Fit ",0.001);
00297
00298
00299
00300 std::cout << "\n\nRedo Chi2 Hist Fit using FUMILI2" << std::endl;
00301 f.SetParameters(p);
00302 fitter.Config().SetMinimizer("Minuit2","Fumili");
00303 ret = fitter.Fit(d, f);
00304 if (ret)
00305 fitter.Result().Print(std::cout);
00306 else {
00307 std::cout << "Chi2 Fit Failed " << std::endl;
00308 return -1;
00309 }
00310 iret |= compareResult(fitter.Result().Chi2(), chi2ref,"1D Histo Fumili2 fit ");
00311
00312
00313 std::cout << "\n\nRedo Chi2 Hist Fit using FUMILI" << std::endl;
00314 f.SetParameters(p);
00315 fitter.Config().SetMinimizer("Fumili");
00316 ret = fitter.Fit(d, f);
00317 if (ret)
00318 fitter.Result().Print(std::cout);
00319 else {
00320 std::cout << "Chi2 Fit Failed " << std::endl;
00321 return -1;
00322 }
00323 iret |= compareResult(fitter.Result().Chi2(), chi2ref,"1D Histo Fumili fit ");
00324
00325
00326 std::cout << "\n\nRedo Chi2 Hist Fit using GSLMultiFit" << std::endl;
00327 f.SetParameters(p);
00328 fitter.Config().SetMinimizer("GSLMultiFit");
00329 ret = fitter.Fit(d, f);
00330 if (ret)
00331 fitter.Result().Print(std::cout);
00332 else {
00333 std::cout << "Chi2 Fit Failed " << std::endl;
00334 return -1;
00335 }
00336 iret |= compareResult(fitter.Result().Chi2(), chi2ref,"1D Histo GSL NLS fit ");
00337
00338
00339 std::cout << "\n\nRedo Chi2 Hist Fit using GSLMultiMin" << std::endl;
00340 f.SetParameters(p);
00341 fitter.Config().SetMinimizer("GSLMultiMin","BFGS2");
00342 ret = fitter.Fit(d, f);
00343 if (ret)
00344 fitter.Result().Print(std::cout);
00345 else {
00346 std::cout << "Chi2 Fit Failed " << std::endl;
00347 return -1;
00348 }
00349 iret |= compareResult(fitter.Result().Chi2(), chi2ref,"1D Histo GSL Minimizer fit ");
00350
00351 return iret;
00352 }
00353
00354
00355 class Func1D : public ROOT::Math::IParamFunction {
00356 public:
00357 void SetParameters(const double *p) { std::copy(p,p+NPar(),fp);}
00358 const double * Parameters() const { return fp; }
00359 ROOT::Math::IGenFunction * Clone() const {
00360 Func1D * f = new Func1D();
00361 f->SetParameters(fp);
00362 return f;
00363 };
00364 unsigned int NPar() const { return 3; }
00365 private:
00366 double DoEvalPar( double x, const double *p) const {
00367 return p[0]*x*x + p[1]*x + p[2];
00368 }
00369 double fp[3];
00370
00371 };
00372
00373
00374 class GradFunc2D : public ROOT::Math::IParamMultiGradFunction {
00375 public:
00376 void SetParameters(const double *p) { std::copy(p,p+NPar(),fp);}
00377 const double * Parameters() const { return fp; }
00378 ROOT::Math::IMultiGenFunction * Clone() const {
00379 GradFunc2D * f = new GradFunc2D();
00380 f->SetParameters(fp);
00381 return f;
00382 };
00383 unsigned int NDim() const { return 2; }
00384 unsigned int NPar() const { return 5; }
00385
00386 void ParameterGradient( const double * x, const double * , double * grad) const {
00387 grad[0] = x[0]*x[0];
00388 grad[1] = x[0];
00389 grad[2] = x[1]*x[1];
00390 grad[3] = x[1];
00391 grad[4] = 1;
00392 }
00393
00394 private:
00395
00396 double DoEvalPar( const double *x, const double * p) const {
00397 return p[0]*x[0]*x[0] + p[1]*x[0] + p[2]*x[1]*x[1] + p[3]*x[1] + p[4];
00398 }
00399
00400
00401
00402
00403
00404
00405
00406
00407 double DoParameterDerivative(const double * x, const double * p, unsigned int ipar) const {
00408 std::vector<double> grad(NPar());
00409 ParameterGradient(x, p, &grad[0] );
00410 return grad[ipar];
00411 }
00412
00413 double fp[5];
00414
00415 };
00416
00417 int testHisto1DPolFit() {
00418
00419
00420
00421 std::string fname("pol2");
00422 TF1 * func = (TF1*)gROOT->GetFunction(fname.c_str());
00423 func->SetParameter(0,1.);
00424 func->SetParameter(1,2.);
00425 func->SetParameter(2,3.0);
00426
00427 TRandom3 rndm;
00428 int iret = 0;
00429 double chi2ref = 0;
00430
00431
00432 TH1D * h2 = new TH1D("h2","h2",30,-5.,5.);
00433
00434 for (int i = 0; i <1000; ++i)
00435 h2->Fill( func->GetRandom() );
00436
00437
00438 ROOT::Fit::BinData d;
00439 ROOT::Fit::FillData(d,h2,func);
00440
00441
00442 printData(d);
00443
00444
00445 Func1D f;
00446
00447 double p[3] = {100,0,3.};
00448 f.SetParameters(p);
00449
00450
00451
00452
00453 std::cout << "\n\nTest histo polynomial fit using Fitter" << std::endl;
00454
00455 ROOT::Fit::Fitter fitter;
00456 bool ret = fitter.Fit(d, f);
00457 if (ret)
00458 fitter.Result().Print(std::cout,true);
00459 else {
00460 std::cout << " Fit Failed " << std::endl;
00461 return -1;
00462 }
00463 chi2ref = fitter.Result().Chi2();
00464
00465
00466 std::cout << "\n******************************\n\t TH1::Fit(pol2) Result \n" << std::endl;
00467 func->SetParameters(p);
00468 h2->Fit(func,"F");
00469 iret |= compareResult(func->GetChisquare(),chi2ref,"TH1::Fit ",0.001);
00470
00471
00472 std::cout << "\n\nTest histo polynomial linear fit " << std::endl;
00473
00474 ROOT::Math::WrappedTF1 pf(*func);
00475
00476 pf.SetParameters(p);
00477
00478 fitter.Config().SetMinimizer("Linear");
00479 ret = fitter.Fit(d, pf);
00480 if (ret) {
00481 fitter.Result().Print(std::cout);
00482 fitter.Result().PrintCovMatrix(std::cout);
00483 }
00484 else {
00485 std::cout << " Fit Failed " << std::endl;
00486 return -1;
00487 }
00488 iret |= compareResult(fitter.Result().Chi2(),chi2ref,"1D histo linear Fit ");
00489
00490
00491 std::cout << "\n******************************\n\t TH1::Fit(pol2) Result with TLinearFitter \n" << std::endl;
00492 func->SetParameters(p);
00493 h2->Fit(func);
00494 iret |= compareResult(func->GetChisquare(),fitter.Result().Chi2(),"TH1::Fit linear",0.001);
00495
00496
00497
00498 return iret;
00499 }
00500
00501 int testHisto2DFit() {
00502
00503
00504
00505
00506 std::string fname("pol2");
00507 TF2 * func = new TF2("f2d",ROOT::Math::ParamFunctor(GradFunc2D() ), 0.,10.,0,10,5);
00508 double p0[5] = { 1.,2.,0.5,1.,3. };
00509 func->SetParameters(p0);
00510 assert(func->GetNpar() == 5);
00511
00512 TRandom3 rndm;
00513 double chi2ref = 0;
00514 int iret = 0;
00515
00516
00517 TH2D * h2 = new TH2D("h2d","h2d",30,0,10.,30,0.,10.);
00518
00519 for (int i = 0; i <10000; ++i) {
00520 double x,y = 0;
00521 func->GetRandom2(x,y);
00522 h2->Fill(x,y);
00523 }
00524
00525 ROOT::Fit::BinData d;
00526 ROOT::Fit::FillData(d,h2,func);
00527
00528
00529
00530
00531
00532 GradFunc2D f;
00533
00534 double p[5] = { 2.,1.,1,2.,100. };
00535 f.SetParameters(p);
00536
00537
00538
00539 ROOT::Fit::Fitter fitter;
00540
00541 fitter.Config().SetMinimizer("Minuit2");
00542
00543 std::cout <<"\ntest 2D histo fit using gradient" << std::endl;
00544 bool ret = fitter.Fit(d, f);
00545 if (ret)
00546 fitter.Result().Print(std::cout);
00547 else {
00548 std::cout << "Gradient Fit Failed " << std::endl;
00549 return -1;
00550 }
00551 chi2ref = fitter.Result().Chi2();
00552
00553
00554 std::cout <<"\ntest result without using gradient" << std::endl;
00555 ROOT::Math::WrappedParamFunction<GradFunc2D *> f2(&f,2,5,p);
00556 ret = fitter.Fit(d, f2);
00557 if (ret)
00558 fitter.Result().Print(std::cout);
00559 else {
00560 std::cout << " Chi2 Fit Failed " << std::endl;
00561 return -1;
00562 }
00563 iret |= compareResult(fitter.Result().Chi2(), chi2ref,"2D histogram chi2 fit");
00564
00565
00566 std::cout <<"\ntest result without gradient and binned likelihood" << std::endl;
00567 f.SetParameters(p);
00568 fitter.SetFunction(static_cast<const ROOT::Math::IParamMultiFunction &>(f) );
00569 fitter.Config().ParSettings(0).SetLimits(0,100);
00570 fitter.Config().ParSettings(1).SetLimits(0,100);
00571 fitter.Config().ParSettings(2).SetLimits(0,100);
00572 fitter.Config().ParSettings(3).SetLimits(0,100);
00573 fitter.Config().ParSettings(4).SetLowerLimit(0);
00574
00575 ret = fitter.LikelihoodFit(d);
00576 if (ret)
00577 fitter.Result().Print(std::cout);
00578 else {
00579 std::cout << "Poisson 2D Bin Likelihood Fit Failed " << std::endl;
00580 return -1;
00581 }
00582
00583
00584 std::cout <<"\ntest result using gradient and binned likelihood" << std::endl;
00585 f.SetParameters(p);
00586 fitter.SetFunction(f);
00587
00588 fitter.Config().ParSettings(0).SetLimits(0,100);
00589 fitter.Config().ParSettings(1).SetLimits(0,100);
00590 fitter.Config().ParSettings(2).SetLimits(0,100);
00591 fitter.Config().ParSettings(3).SetLimits(0,100);
00592 fitter.Config().ParSettings(4).SetLowerLimit(0);
00593 ret = fitter.LikelihoodFit(d);
00594 if (ret) {
00595
00596 f.SetParameters(&(fitter.Result().Parameters().front()) );
00597 ret = fitter.LikelihoodFit(d,f);
00598 }
00599 if (ret)
00600 fitter.Result().Print(std::cout);
00601 else {
00602 std::cout << "Gradient Bin Likelihood Fit Failed " << std::endl;
00603 return -1;
00604 }
00605
00606
00607 std::cout <<"\ntest result using linear fitter" << std::endl;
00608 fitter.Config().SetMinimizer("Linear");
00609 f.SetParameters(p);
00610 ret = fitter.Fit(d, f);
00611 if (ret)
00612 fitter.Result().Print(std::cout);
00613 else {
00614 std::cout << "Linear 2D Fit Failed " << std::endl;
00615 return -1;
00616 }
00617 iret |= compareResult(fitter.Result().Chi2(), chi2ref,"2D histogram linear fit");
00618
00619
00620
00621 TGraph2D g2(h2);
00622
00623 std::cout <<"\ntest using TGraph2D" << std::endl;
00624 ROOT::Fit::BinData d2;
00625 ROOT::Fit::FillData(d2,&g2,func);
00626
00627 std::cout << "data size from graph " << d2.Size() << std::endl;
00628
00629 f2.SetParameters(p);
00630 fitter.Config().SetMinimizer("Minuit2");
00631 ret = fitter.Fit(d2, f2);
00632 if (ret)
00633 fitter.Result().Print(std::cout);
00634 else {
00635 std::cout << " TGraph2D Fit Failed " << std::endl;
00636 return -1;
00637 }
00638 double chi2ref2 = fitter.Result().Chi2();
00639
00640
00641 std::cout << "\n******************************\n\t TGraph::Fit Result \n" << std::endl;
00642 func->SetParameters(p);
00643 g2.Fit(func);
00644
00645 std::cout <<"\ntest using TGraph2D and gradient function" << std::endl;
00646 f.SetParameters(p);
00647 ret = fitter.Fit(d2, f);
00648 if (ret)
00649 fitter.Result().Print(std::cout);
00650 else {
00651 std::cout << " TGraph2D Grad Fit Failed " << std::endl;
00652 return -1;
00653 }
00654 iret |= compareResult(fitter.Result().Chi2(), chi2ref2,"TGraph2D chi2 fit");
00655
00656 std::cout <<"\ntest using TGraph2DErrors - error only in Z" << std::endl;
00657 TGraph2DErrors g2err(g2.GetN() );
00658
00659 for (int i = 0; i < g2.GetN(); ++i) {
00660 double x = g2.GetX()[i];
00661 double y = g2.GetY()[i];
00662 g2err.SetPoint(i,x,y,g2.GetZ()[i]);
00663 g2err.SetPointError(i,0,0,h2->GetBinError(h2->FindBin(x,y) ) );
00664 }
00665 func->SetParameters(p);
00666
00667 f.SetParameters(p);
00668 ROOT::Fit::BinData d3;
00669 ROOT::Fit::FillData(d3,&g2err,func);
00670 ret = fitter.Fit(d3, f);
00671 if (ret)
00672 fitter.Result().Print(std::cout);
00673 else {
00674 std::cout << " TGraph2DErrors Fit Failed " << std::endl;
00675 return -1;
00676 }
00677
00678 iret |= compareResult(fitter.Result().Chi2(), chi2ref,"TGraph2DErrors chi2 fit");
00679
00680
00681
00682 std::cout <<"\ntest using TGraph2DErrors - with error in X,Y,Z" << std::endl;
00683 for (int i = 0; i < g2err.GetN(); ++i) {
00684 double x = g2.GetX()[i];
00685 double y = g2.GetY()[i];
00686 g2err.SetPointError(i,0.5* h2->GetXaxis()->GetBinWidth(1),0.5*h2->GetXaxis()->GetBinWidth(1),h2->GetBinError(h2->FindBin(x,y) ) );
00687 }
00688 std::cout << "\n******************************\n\t TGraph2DErrors::Fit Result \n" << std::endl;
00689 func->SetParameters(p);
00690 g2err.Fit(func);
00691
00692
00693 return iret;
00694 }
00695
00696
00697 int testUnBin1DFit() {
00698
00699 int iret = 0;
00700
00701 std::string fname("gausn");
00702 TF1 * func = (TF1*)gROOT->GetFunction(fname.c_str());
00703
00704 TRandom3 rndm;
00705
00706 int n = 100;
00707 ROOT::Fit::UnBinData d(n);
00708
00709 for (int i = 0; i <n; ++i)
00710 d.Add( rndm.Gaus(0,1) );
00711
00712
00713
00714
00715
00716 ROOT::Math::WrappedMultiTF1 wf(*func);
00717
00718 ROOT::Math::IParamMultiFunction & f = wf;
00719 double p[3] = {1,2,10.};
00720 f.SetParameters(p);
00721
00722
00723
00724
00725 ROOT::Fit::Fitter fitter;
00726 fitter.SetFunction(f);
00727 std::cout << "fix parameter 0 " << " to value " << f.Parameters()[0] << std::endl;
00728 fitter.Config().ParSettings(0).Fix();
00729
00730 std::cout << "set lower range to 0 for sigma " << std::endl;
00731 fitter.Config().ParSettings(2).SetLowerLimit(0);
00732
00733 #ifdef DEBUG
00734 fitter.Config().MinimizerOptions().SetPrintLevel(3);
00735 #endif
00736
00737
00738
00739
00740
00741
00742
00743 fitter.Config().SetMinimizer("Minuit2");
00744
00745
00746 bool ret = fitter.Fit(d);
00747 if (ret)
00748 fitter.Result().Print(std::cout);
00749 else {
00750 std::cout << "Unbinned Likelihood Fit Failed " << std::endl;
00751 iret |= 1;
00752 }
00753 double lref = fitter.Result().MinFcnValue();
00754
00755 std::cout << "\n\nRedo Fit using FUMILI2" << std::endl;
00756 f.SetParameters(p);
00757 fitter.Config().SetMinimizer("Fumili2");
00758
00759 fitter.SetFunction(f);
00760 fitter.Config().ParSettings(0).Fix();
00761
00762 fitter.Config().ParSettings(2).SetLowerLimit(0);
00763
00764 ret = fitter.Fit(d);
00765 if (ret)
00766 fitter.Result().Print(std::cout);
00767 else {
00768 std::cout << "Unbinned Likelihood Fit using FUMILI2 Failed " << std::endl;
00769 iret |= 1;
00770 }
00771
00772 iret |= compareResult(fitter.Result().MinFcnValue(), lref,"1D unbin FUMILI2 fit");
00773
00774 std::cout << "\n\nRedo Fit using FUMILI" << std::endl;
00775 f.SetParameters(p);
00776 fitter.Config().SetMinimizer("Fumili");
00777
00778
00779 fitter.SetFunction(f);
00780 fitter.Config().ParSettings(0).Fix();
00781
00782 fitter.Config().ParSettings(2).SetLowerLimit(0);
00783
00784 ret = fitter.Fit(d);
00785 if (ret)
00786 fitter.Result().Print(std::cout);
00787 else {
00788 std::cout << "Unbinned Likelihood Fit using FUMILI Failed " << std::endl;
00789 iret |= 1;
00790 }
00791
00792 iret |= compareResult(fitter.Result().MinFcnValue(), lref,"1D unbin FUMILI fit");
00793
00794
00795 return iret;
00796 }
00797
00798 int testGraphFit() {
00799
00800 int iret = 0;
00801
00802
00803
00804 double x[5] = {1,2,3,4,5};
00805 double y[5] = {2.1, 3.5, 6.5, 8.8, 9.5};
00806 double ex[5] = {.3,.3,.3,.3,.3};
00807 double ey[5] = {.5,.5,.5,.5,.5};
00808 double eyl[5] = {.2,.2,.2,.2,.2};
00809 double eyh[5] = {.8,.8,.8,.8,.8};
00810
00811 std::cout << "\n********************************************************\n";
00812 std::cout << "Test simple fit of Tgraph of 5 points" << std::endl;
00813 std::cout << "\n********************************************************\n";
00814
00815
00816 double p[2] = {1,1};
00817 TF1 * func = new TF1("f","pol1",0,10);
00818 func->SetParameters(p);
00819
00820 ROOT::Math::WrappedMultiTF1 wf(*func);
00821
00822 ROOT::Math::IParamMultiFunction & f = wf;
00823 f.SetParameters(p);
00824
00825 ROOT::Fit::Fitter fitter;
00826 fitter.SetFunction(f);
00827
00828
00829 std::cout <<"\ntest TGraph (no errors) " << std::endl;
00830 TGraph gr(5, x,y);
00831
00832 ROOT::Fit::BinData dgr;
00833 ROOT::Fit::FillData(dgr,&gr);
00834
00835
00836
00837 f.SetParameters(p);
00838 bool ret = fitter.Fit(dgr, f);
00839 if (ret)
00840 fitter.Result().Print(std::cout);
00841 else {
00842 std::cout << "Chi2 Graph Fit Failed " << std::endl;
00843 return -1;
00844 }
00845 double chi2ref = fitter.Result().Chi2();
00846
00847
00848
00849 std::cout << "\n******************************\n\t TGraph::Fit Result \n" << std::endl;
00850 func->SetParameters(p);
00851 gr.Fit(func,"F");
00852
00853 iret |= compareResult(func->GetChisquare(),fitter.Result().Chi2(),"TGraph::Fit ",0.001);
00854
00855
00856 std::cout <<"\ntest TGraphErrors " << std::endl;
00857 TGraphErrors grer(5, x,y,ex,ey);
00858
00859 ROOT::Fit::BinData dgrer;
00860 dgrer.Opt().fCoordErrors = true;
00861 ROOT::Fit::FillData(dgrer,&grer);
00862
00863
00864
00865 f.SetParameters(p);
00866 ret = fitter.Fit(dgrer, f);
00867 if (ret)
00868 fitter.Result().Print(std::cout);
00869 else {
00870 std::cout << "Chi2 Graph Fit Failed " << std::endl;
00871 return -1;
00872 }
00873
00874 iret |= compareResult(fitter.Result().Chi2(),chi2ref,"TGraphErrors fit with coord errors",0.8);
00875
00876
00877
00878 std::cout << "\n******************************\n\t TGraphErrors::Fit Result \n" << std::endl;
00879 func->SetParameters(p);
00880 grer.Fit(func,"F");
00881 iret |= compareResult(func->GetChisquare(),fitter.Result().Chi2(),"TGraphErrors::Fit ",0.001);
00882
00883 chi2ref = fitter.Result().Chi2();
00884
00885 std::cout <<"\ntest TGraphAsymmErrors " << std::endl;
00886 TGraphAsymmErrors graer(5, x,y,ex,ex,eyl, eyh);
00887
00888 ROOT::Fit::BinData dgraer;
00889
00890 dgraer.Opt().fCoordErrors = true;
00891 dgraer.Opt().fAsymErrors = true;
00892 ROOT::Fit::FillData(dgraer,&graer);
00893
00894
00895 f.SetParameters(p);
00896 ret = fitter.Fit(dgraer, f);
00897 if (ret)
00898 fitter.Result().Print(std::cout);
00899 else {
00900 std::cout << "Chi2 Graph Fit Failed " << std::endl;
00901 return -1;
00902 }
00903
00904
00905
00906 std::cout << "\n******************************\n\t TGraphAsymmErrors::Fit Result \n" << std::endl;
00907 func->SetParameters(p);
00908 graer.Fit(func,"F");
00909 iret |= compareResult(func->GetChisquare(),fitter.Result().Chi2(),"TGraphAsymmErrors::Fit ",0.001);
00910
00911
00912
00913 return iret;
00914 }
00915
00916
00917 template<typename Test>
00918 int testFit(Test t, std::string name) {
00919 std::cout << name << "\n\t\t";
00920 int iret = t();
00921 std::cout << "\n" << name << ":\t\t";
00922 if (iret == 0)
00923 std::cout << "OK" << std::endl;
00924 else
00925 std::cout << "Failed" << std::endl;
00926 return iret;
00927 }
00928
00929 int main() {
00930
00931 int iret = 0;
00932 iret |= testFit( testHisto1DFit, "Histogram1D Fit");
00933 iret |= testFit( testHisto1DPolFit, "Histogram1D Polynomial Fit");
00934 iret |= testFit( testHisto2DFit, "Histogram2D Gradient Fit");
00935 iret |= testFit( testUnBin1DFit, "Unbin 1D Fit");
00936 iret |= testFit( testGraphFit, "Graph 1D Fit");
00937
00938 std::cout << "\n******************************\n";
00939 if (iret) std::cerr << "\n\t testFit FAILED !!!!!!!!!!!!!!!! \n";
00940 else std::cout << "\n\t testFit all OK \n";
00941 return iret;
00942 }
00943