ROOT Analysis for the Inclusive Detachted Dilepton Trigger Lines
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.

213 lines
8.6 KiB

  1. #include "TH1D.h"
  2. #include "TH2D.h"
  3. #include "THStack.h"
  4. #include "TGraph.h"
  5. #include "TTree.h"
  6. #include "TChain.h"
  7. #include "TFile.h"
  8. #include "TCanvas.h"
  9. #include "TROOT.h"
  10. #include "TStyle.h"
  11. #include "TColor.h"
  12. #include "TLorentzVector.h"
  13. #include "TRandom3.h"
  14. #include "TLorentzVector.h"
  15. #include "RooDataHist.h"
  16. #include "RooRealVar.h"
  17. #include "RooPlot.h"
  18. #include "RooGaussian.h"
  19. #include "RooExponential.h"
  20. #include "RooRealConstant.h"
  21. #include "RooAddPdf.h"
  22. #include "RooFitResult.h"
  23. #include "RooProduct.h"
  24. #include "RooCrystalBall.h"
  25. #include "RooBreitWigner.h"
  26. #include "RooArgSet.h"
  27. #include "RooFFTConvPdf.h"
  28. #include "RooNovosibirsk.h"
  29. #include <string>
  30. #include <iostream>
  31. #include <cmath>
  32. const int nBins = 70;
  33. const Double_t MASS_HIST_MIN = 5150.;
  34. const Double_t MASS_HIST_MAX = 6000.;
  35. const char *B0_DECAY = "Hlt2RD_B0ToKpPimMuMu_woHtl1cut";
  36. const char *Bu_DECAY = "Hlt2RD_BuToKpMuMu_woHtl1cut";
  37. const char *TITLE = Bu_DECAY;
  38. const char *HIST_TITLE = "Hlt2RD_BuToKpMuMu w/o Hlt1 Decision Cut";
  39. const char *X_AXIS = "m(K^{+}#mu^{+}#mu^{-})";
  40. void CreateRooFitAndDraw(TH1D *hist, int fitting_entries);
  41. int analysis_turbo()
  42. {
  43. std::string tree = "";
  44. if (TITLE == B0_DECAY)
  45. {
  46. tree = "Hlt2RD_B0ToKpPimMuMu";
  47. }
  48. else
  49. {
  50. tree = "Hlt2RD_BuToKpMuMu";
  51. }
  52. TChain *data_chain = new TChain(TString::Format("%s/DecayTree", tree.c_str()).Data());
  53. data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/Collision23_Beam6800GeV-VeloClosed-MagDown-Excl-UT_RealData_SprucingPass23r1_94000000_RD.root");
  54. Double_t B_M, Jpsi_M, B0_M, Kst_M;
  55. Bool_t Hlt1TrackMVADecision, Hlt1TwoTrackMVADecision;
  56. if (TITLE == B0_DECAY)
  57. {
  58. data_chain->SetBranchAddress("B0_M", &B0_M);
  59. data_chain->SetBranchAddress("Jpsi_M", &Jpsi_M);
  60. data_chain->SetBranchAddress("Kst0_M", &Kst_M);
  61. data_chain->SetBranchAddress("Hlt1TrackMVADecision", &Hlt1TrackMVADecision);
  62. data_chain->SetBranchAddress("Hlt1TwoTrackMVADecision", &Hlt1TwoTrackMVADecision);
  63. }
  64. else
  65. {
  66. data_chain->SetBranchAddress("B_M", &B_M);
  67. data_chain->SetBranchAddress("Jpsi_M", &Jpsi_M);
  68. data_chain->SetBranchAddress("Hlt1TrackMVADecision", &Hlt1TrackMVADecision);
  69. data_chain->SetBranchAddress("Hlt1TwoTrackMVADecision", &Hlt1TwoTrackMVADecision);
  70. }
  71. TH1D *h1_B_M = new TH1D("h1_B_M", HIST_TITLE, nBins, MASS_HIST_MIN, MASS_HIST_MAX);
  72. h1_B_M->GetXaxis()->SetTitle(X_AXIS);
  73. int fitting_entries = 0;
  74. unsigned int entries = data_chain->GetEntries();
  75. for (unsigned int i = 0; i < entries; i++)
  76. {
  77. data_chain->GetEntry(i);
  78. if (TITLE == B0_DECAY)
  79. {
  80. if ((TMath::Abs(Jpsi_M - 3096.9) < 100) && (B0_M > 4500) && (B0_M < 6000) && (TMath::Abs(Kst_M - 895.55) < 50) /*&& ((Hlt1TrackMVADecision) || (Hlt1TwoTrackMVADecision))*/)
  81. {
  82. h1_B_M->Fill(B0_M);
  83. if (MASS_HIST_MIN <= B0_M && B0_M <= MASS_HIST_MAX)
  84. {
  85. fitting_entries++;
  86. }
  87. }
  88. }
  89. else
  90. {
  91. if ((TMath::Abs(Jpsi_M - 3096.9) < 100.) && (B_M > 4500.) && (B_M < 6000.) /*&& ((Hlt1TrackMVADecision) || (Hlt1TwoTrackMVADecision))*/)
  92. {
  93. h1_B_M->Fill(B_M);
  94. if (MASS_HIST_MIN <= B_M && B_M <= MASS_HIST_MAX)
  95. {
  96. fitting_entries++;
  97. }
  98. }
  99. }
  100. if (i % 10000 == 0)
  101. {
  102. std::cout << "["
  103. << TITLE
  104. << "] Processed event: " << i << " (" << TString::Format("%.2f", ((double)i / (double)entries) * 100.) << "%)" << std::endl;
  105. }
  106. }
  107. std::cout << "["
  108. << TITLE
  109. << "] Processed event: " << entries << " (" << TString::Format("%.2f", 100.) << "%)" << std::endl;
  110. TCanvas *c1 = new TCanvas("c1", "c1", 0, 0, 800, 600);
  111. h1_B_M->Draw();
  112. c1->Draw();
  113. c1->SaveAs(TString::Format("images/root_hist_%s_bmass.pdf", TITLE).Data());
  114. CreateRooFitAndDraw(h1_B_M, fitting_entries);
  115. return 0;
  116. }
  117. void CreateRooFitAndDraw(TH1D *hist, int fitting_entries)
  118. {
  119. RooRealVar roo_var_mass("var_mass", "B Mass Variable", MASS_HIST_MIN, MASS_HIST_MAX);
  120. roo_var_mass.setRange("fitting_range", MASS_HIST_MIN, MASS_HIST_MAX);
  121. RooDataHist roohist_B_M("roohist_B_M", "B Mass Histogram", roo_var_mass, RooFit::Import(*hist));
  122. RooPlot *roo_frame_mass = roo_var_mass.frame(RooFit::Title(HIST_TITLE));
  123. roohist_B_M.plotOn(roo_frame_mass, RooFit::Binning(nBins), RooFit::Name("B Mass Distribution"));
  124. roo_frame_mass->GetXaxis()->SetTitle(X_AXIS);
  125. // Crystal Ball for Signal
  126. RooRealVar var_mass_x0("var_mass_x0", "#mu", 5278., 5170., 5500.);
  127. RooRealVar var_mass_sigmaLR("var_mass_sigmaLR", "#sigma_{LR}", 16., 5., 40.);
  128. // Same Variables for Left and Right Tail
  129. // RooRealVar var_mass_alphaL("var_mass_alphaL", "#alpha_{L}", 2., 0., 4.);
  130. // RooRealVar var_mass_nL("var_mass_nL", "n_{L}", 5., 0., 15.);
  131. // RooRealVar var_mass_alphaR("var_mass_alphaR", "#alpha_{R}", 2., 0., 4.);
  132. // RooRealVar var_mass_nR("var_mass_nR", "n_{R}", 5., 0., 15.);
  133. // B+/-
  134. auto var_mass_alphaL = RooRealConstant::value(1.912);
  135. auto var_mass_nL = RooRealConstant::value(1.125);
  136. auto var_mass_alphaR = RooRealConstant::value(2.447);
  137. auto var_mass_nR = RooRealConstant::value(1.491);
  138. // B0
  139. // auto var_mass_alphaL = RooRealConstant::value(1.948);
  140. // auto var_mass_nL = RooRealConstant::value(0.618);
  141. // auto var_mass_alphaR = RooRealConstant::value(2.320);
  142. // auto var_mass_nR = RooRealConstant::value(0.473);
  143. RooCrystalBall sig_cb("sig_cb", "Signal Crystal Ball", roo_var_mass, var_mass_x0, var_mass_sigmaLR, var_mass_alphaL, var_mass_nL, var_mass_alphaR, var_mass_nR);
  144. // Exponential for Background
  145. RooRealVar var_mass_bkg_c("var_mass_bkg_c", "#lambda", -0.0014, -0.004, -0.000);
  146. RooExponential bkg_exp("bkg_exp", "Exp Background", roo_var_mass, var_mass_bkg_c);
  147. RooRealVar var_mass_nsig("nsig", "Mass N Signal", 0., hist->GetEntries());
  148. RooRealVar var_mass_nbkg("nbkg", "Mass N Background", 0., hist->GetEntries());
  149. // TString pdf_name = TString::Format("%s_sigplusbkg", var_id.c_str());
  150. RooAddPdf sigplusbkg("sigplusbkg", "Sig and Bkg PDF", RooArgList(sig_cb, bkg_exp), RooArgList(var_mass_nsig, var_mass_nbkg));
  151. RooFitResult *fitres = sigplusbkg.fitTo(roohist_B_M, RooFit::Save(), RooFit::PrintLevel(1), RooFit::Range("fitting_range"));
  152. sigplusbkg.plotOn(roo_frame_mass, RooFit::VisualizeError(*fitres, 1), RooFit::FillColor(kYellow -7), RooFit::FillStyle(3144));
  153. sigplusbkg.plotOn(roo_frame_mass, RooFit::LineColor(kRed), RooFit::LineStyle(kSolid), RooFit::Range("fitting_range"), RooFit::Name("sigplusbkg"));
  154. sigplusbkg.plotOn(roo_frame_mass, RooFit::Components(RooArgSet(bkg_exp)), RooFit::LineColor(kBlue - 7), RooFit::LineStyle(kDashed), RooFit::Range("fitting_range"), RooFit::Name("bkg_exp"));
  155. sigplusbkg.plotOn(roo_frame_mass, RooFit::Components(RooArgSet(sig_cb)), RooFit::LineColor(kRed - 7), RooFit::LineStyle(kDashed), RooFit::Range("fitting_range"), RooFit::Name("sig_cb"));
  156. TCanvas *c = new TCanvas("roofit_c", "roofit_c", 0, 0, 800, 600);
  157. roo_frame_mass->Draw();
  158. TLegend *leg1 = new TLegend(0.50, 0.58, 0.87, 0.89);
  159. leg1->SetFillColor(kWhite);
  160. leg1->SetLineColor(kBlack);
  161. leg1->AddEntry(roo_frame_mass->findObject("sigplusbkg"), "Signal + Background", "LP");
  162. leg1->AddEntry(roo_frame_mass->findObject("sig_cb"), "Signal", "LP");
  163. leg1->AddEntry(roo_frame_mass->findObject("bkg_exp"), "Background", "LP");
  164. leg1->AddEntry((TObject *)0, "", "");
  165. leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", var_mass_x0.getTitle().Data(), var_mass_x0.getVal(), var_mass_x0.getError()).Data(), "");
  166. leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", var_mass_sigmaLR.getTitle().Data(), var_mass_sigmaLR.getVal(), var_mass_sigmaLR.getError()).Data(), "");
  167. leg1->AddEntry((TObject *)0, TString::Format("%s = %.6f #pm %.6f", var_mass_bkg_c.getTitle().Data(), var_mass_bkg_c.getVal(), var_mass_bkg_c.getError()).Data(), "");
  168. leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", "S", var_mass_nsig.getVal(), var_mass_nsig.getError()).Data(), "");
  169. leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", "B", var_mass_nbkg.getVal(), var_mass_nbkg.getError()).Data(), "");
  170. leg1->AddEntry((TObject *)0, TString::Format("Entries: %d", fitting_entries), "");
  171. leg1->Draw();
  172. c->Draw();
  173. c->SaveAs(TString::Format("images/data/roofit_hist_%s_fitted_bmass.pdf", TITLE).Data());
  174. }