data analysis scripts
You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

145 lines
4.1 KiB

  1. //+ Combined (simultaneous) fit of two histogram with separate functions
  2. // and some common parameters
  3. //
  4. // See http://root.cern.ch/phpBB3//viewtopic.php?f=3&t=11740#p50908
  5. // for a modified version working with Fumili or GSLMultiFit
  6. //
  7. // N.B. this macro must be compiled with ACliC
  8. //
  9. //Author: L. Moneta - Dec 2010
  10. #include "Fit/Fitter.h"
  11. #include "Fit/BinData.h"
  12. #include "Fit/Chi2FCN.h"
  13. #include "TH1.h"
  14. #include "TList.h"
  15. #include "Math/WrappedMultiTF1.h"
  16. #include "HFitInterface.h"
  17. #include "TCanvas.h"
  18. #include "TStyle.h"
  19. // definition of shared parameter
  20. // background function
  21. int iparB[2] = { 0, // exp amplitude in B histo
  22. 2 // exp common parameter
  23. };
  24. // signal + background function
  25. int iparSB[5] = { 1, // exp amplitude in S+B histo
  26. 2, // exp common parameter
  27. 3, // gaussian amplitude
  28. 4, // gaussian mean
  29. 5 // gaussian sigma
  30. };
  31. struct GlobalChi2 {
  32. GlobalChi2( ROOT::Math::IMultiGenFunction & f1,
  33. ROOT::Math::IMultiGenFunction & f2) :
  34. fChi2_1(&f1), fChi2_2(&f2) {}
  35. // parameter vector is first background (in common 1 and 2)
  36. // and then is signal (only in 2)
  37. double operator() (const double *par) const {
  38. double p1[2];
  39. for (int i = 0; i < 2; ++i) p1[i] = par[iparB[i] ];
  40. double p2[5];
  41. for (int i = 0; i < 5; ++i) p2[i] = par[iparSB[i] ];
  42. return (*fChi2_1)(p1) + (*fChi2_2)(p2);
  43. }
  44. const ROOT::Math::IMultiGenFunction * fChi2_1;
  45. const ROOT::Math::IMultiGenFunction * fChi2_2;
  46. };
  47. void combinedFit() {
  48. #if defined(__CINT__) && !defined(__MAKECINT__)
  49. cout << "ERROR: This tutorial can run only using ACliC, you must run it by doing: " << endl;
  50. cout << "\t .x $ROOTSYS/tutorials/fit/combinedFit.C+" << endl;
  51. return;
  52. #endif
  53. TH1D * hB = new TH1D("hB","histo B",100,0,100);
  54. TH1D * hSB = new TH1D("hSB","histo S+B",100, 0,100);
  55. TF1 * fB = new TF1("fB","expo",0,100);
  56. fB->SetParameters(1,-0.05);
  57. hB->FillRandom("fB");
  58. TF1 * fS = new TF1("fS","gaus",0,100);
  59. fS->SetParameters(1,30,5);
  60. hSB->FillRandom("fB",2000);
  61. hSB->FillRandom("fS",1000);
  62. // perform now global fit
  63. TF1 * fSB = new TF1("fSB","expo + gaus(2)",0,100);
  64. ROOT::Math::WrappedMultiTF1 wfB(*fB,1);
  65. ROOT::Math::WrappedMultiTF1 wfSB(*fSB,1);
  66. ROOT::Fit::DataOptions opt;
  67. ROOT::Fit::DataRange rangeB;
  68. // set the data range
  69. rangeB.SetRange(10,90);
  70. ROOT::Fit::BinData dataB(opt,rangeB);
  71. ROOT::Fit::FillData(dataB, hB);
  72. ROOT::Fit::DataRange rangeSB;
  73. rangeSB.SetRange(10,50);
  74. ROOT::Fit::BinData dataSB(opt,rangeSB);
  75. ROOT::Fit::FillData(dataSB, hSB);
  76. ROOT::Fit::Chi2Function chi2_B(dataB, wfB);
  77. ROOT::Fit::Chi2Function chi2_SB(dataSB, wfSB);
  78. GlobalChi2 globalChi2(chi2_B, chi2_SB);
  79. ROOT::Fit::Fitter fitter;
  80. const int Npar = 6;
  81. double par0[Npar] = { 5,5,-0.1,100, 30,10};
  82. // create before the parameter settings in order to fix or set range on them
  83. fitter.Config().SetParamsSettings(6,par0);
  84. // fix 5-th parameter
  85. fitter.Config().ParSettings(4).Fix();
  86. // set limits on the third and 4-th parameter
  87. fitter.Config().ParSettings(2).SetLimits(-10,-1.E-4);
  88. fitter.Config().ParSettings(3).SetLimits(0,10000);
  89. fitter.Config().ParSettings(3).SetStepSize(5);
  90. fitter.Config().MinimizerOptions().SetPrintLevel(0);
  91. fitter.Config().SetMinimizer("Minuit2","Migrad");
  92. // fit FCN function directly
  93. // (specify optionally data size and flag to indicate that is a chi2 fit)
  94. fitter.FitFCN(6,globalChi2,0,dataB.Size()+dataSB.Size(),true);
  95. ROOT::Fit::FitResult result = fitter.Result();
  96. result.Print(std::cout);
  97. TCanvas * c1 = new TCanvas("Simfit","Simultaneous fit of two histograms",
  98. 10,10,700,700);
  99. c1->Divide(1,2);
  100. c1->cd(1);
  101. gStyle->SetOptFit(1111);
  102. fB->SetFitResult( result, iparB);
  103. fB->SetRange(rangeB().first, rangeB().second);
  104. fB->SetLineColor(kBlue);
  105. hB->GetListOfFunctions()->Add(fB);
  106. hB->Draw();
  107. c1->cd(2);
  108. fSB->SetFitResult( result, iparSB);
  109. fSB->SetRange(rangeSB().first, rangeSB().second);
  110. fSB->SetLineColor(kRed);
  111. hSB->GetListOfFunctions()->Add(fSB);
  112. hSB->Draw();
  113. }