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.

236 lines
10 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. #include <filesystem>
  33. const int nBins = 70;
  34. const Double_t MASS_HIST_MIN = 5100.;
  35. const Double_t MASS_HIST_MAX = 6000.;
  36. const Double_t K_MASS = 493.677;
  37. const int PID_KAON = 321;
  38. const int PID_PION = 211;
  39. const int PID_MUON = 13;
  40. const char *NAME = "BuToKpMuMu";
  41. const char *X_AXIS = "m(K^{+}#mu^{+}#mu^{-})";//"m(K^{+}#pi^{-}#mu^{+}#mu^{-})";
  42. const char *B_NAME = "B";
  43. const bool HAS_PI = false;
  44. TString quickformat_name(const char *format)
  45. {
  46. return TString::Format(format, NAME);
  47. };
  48. void CreateRooFitAndDraw(TH1D *hist, double selected_entries, double total_entries);
  49. int simulation()
  50. {
  51. TChain *data_chain = new TChain(TString::Format("%s_noPID/DecayTree", NAME));
  52. // data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/MC/rd_btoxll_simulation_fullstream_v0r0p6671378_BuToKpMuMu_12143001_magdown.root");
  53. data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/MC/rd_btoxll_simulation_turbo_v0r0p6657752_BuToKpMuMu_12143001_magdown.root");
  54. // data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/MC/rd_btoxll_simulation_turbo_v0r0p6657752_B0ToKpPimMuMu_11144002_magdown.root");
  55. /*
  56. rd_btoxll_simulation_fullstream_v0r0p6671378_B0ToKpPimMuMu_11144002_magdown.root
  57. rd_btoxll_simulation_fullstream_v0r0p6671378_BuToKpMuMu_12143001_magdown.root
  58. rd_btoxll_simulation_turbo_v0r0p6657752_B0ToKpPimMuMu_11144002_magdown.root
  59. rd_btoxll_simulation_turbo_v0r0p6657752_BuToKpMuMu_12143001_magdown.root
  60. */
  61. Int_t B_BKGCAT;
  62. data_chain->SetBranchAddress(TString::Format("%s_BKGCAT", B_NAME).Data(), &B_BKGCAT);
  63. Double_t B_M;
  64. Int_t L1_TRUEID, L2_TRUEID, K_TRUEID, Pi_TRUEID;
  65. data_chain->SetBranchAddress("L1_TRUEID", &L1_TRUEID);
  66. data_chain->SetBranchAddress("L2_TRUEID", &L2_TRUEID);
  67. data_chain->SetBranchAddress(TString::Format("%s_M", B_NAME).Data(), &B_M);
  68. data_chain->SetBranchAddress("K_TRUEID", &K_TRUEID);
  69. if (HAS_PI)
  70. {
  71. data_chain->SetBranchAddress("Pi_TRUEID", &Pi_TRUEID);
  72. }
  73. TH1D *h1_B_M = new TH1D("h1_B_M", "B Mass", nBins, MASS_HIST_MIN, MASS_HIST_MAX);
  74. TH1D *h1_B_M_trueid = new TH1D("h1_B_M_trueid", "B Mass, TrueID matched", nBins, MASS_HIST_MIN, MASS_HIST_MAX);
  75. TH1I *h1_B_BKGCAT = new TH1I("h1_B_BKGCAT", "B Background Category", nBins, 0, 120);
  76. TH2D *h2_B_M_vs_B_BKGCAT = new TH2D("h2_B_M_vs_B_BKGCAT", "B Mass vs Background Category", nBins, MASS_HIST_MIN, MASS_HIST_MAX, nBins, 0, 120);
  77. h1_B_M->GetXaxis()->SetTitle(X_AXIS);
  78. h1_B_M_trueid->GetXaxis()->SetTitle(X_AXIS);
  79. unsigned int entries = data_chain->GetEntries();
  80. unsigned int selected_entries = 0;
  81. for (size_t i = 0; i < entries; i++)
  82. {
  83. data_chain->GetEntry(i);
  84. if (B_BKGCAT == 0 && TMath::Abs(L1_TRUEID) == 13 && L2_TRUEID == -L1_TRUEID)
  85. {
  86. if (((HAS_PI) && (TMath::Abs(K_TRUEID) == PID_KAON) && (TMath::Abs(Pi_TRUEID) == PID_PION)) || ((!HAS_PI) && (TMath::Abs(K_TRUEID) == PID_KAON)))
  87. {
  88. h1_B_M_trueid->Fill(B_M);
  89. selected_entries++;
  90. }
  91. }
  92. h1_B_M->Fill(B_M);
  93. h1_B_BKGCAT->Fill(B_BKGCAT);
  94. h2_B_M_vs_B_BKGCAT->Fill(B_M, B_BKGCAT);
  95. if (i % 10000 == 0)
  96. {
  97. std::cout << "["
  98. << NAME
  99. << "] Processed event: " << i << " (" << TString::Format("%.2f", ((double)i / (double)entries) * 100.) << "%), "
  100. << "Current Efficiency: " << TString::Format("%.2f%%", ((double)selected_entries / (double)entries) * 100.)
  101. << std::endl;
  102. }
  103. }
  104. std::cout << "["
  105. << NAME
  106. << "] Processed event: " << entries << " (" << TString::Format("%.2f", 100.) << "%)" << std::endl;
  107. h1_B_M->SetLineColor(kBlue);
  108. h1_B_M_trueid->SetLineColor(kMagenta + 2);
  109. h1_B_M_trueid->SetFillColor(kMagenta + 2);
  110. h1_B_M_trueid->SetFillStyle(3244);
  111. h1_B_BKGCAT->SetLineColor(kBlue);
  112. h1_B_BKGCAT->SetFillColor(kBlue);
  113. h1_B_BKGCAT->SetFillStyle(3351);
  114. std::filesystem::create_directory(TString::Format("images/sim/%s_Turbo", NAME).Data());
  115. TCanvas *c1 = new TCanvas("c1", "c1", 0, 0, 800, 600);
  116. h1_B_M->SetStats(0);
  117. h1_B_M_trueid->SetStats(0);
  118. h1_B_M->Draw();
  119. h1_B_M_trueid->Draw("SAME");
  120. c1->BuildLegend(0.58, 0.65, 0.85, 0.87);
  121. c1->Draw();
  122. c1->SaveAs(quickformat_name("images/sim/%s/h1_B_M.pdf"));
  123. TCanvas *c2 = new TCanvas("c2", "c2", 0, 0, 800, 600);
  124. h1_B_BKGCAT->SetStats(0);
  125. h1_B_BKGCAT->Draw();
  126. c2->BuildLegend(0.58, 0.77, 0.85, 0.87);
  127. c2->Draw();
  128. c2->SaveAs(quickformat_name("images/sim/%s/h1_B_BKGCAT.pdf"));
  129. TCanvas *c3 = new TCanvas("c3", "c3", 0, 0, 800, 600);
  130. h2_B_M_vs_B_BKGCAT->SetStats(0);
  131. h2_B_M_vs_B_BKGCAT->Draw("COLZ");
  132. c3->Draw();
  133. c3->SaveAs(quickformat_name("images/sim/%s/h2_B_M_vs_B_BKGCAT.pdf"));
  134. CreateRooFitAndDraw(h1_B_M_trueid, selected_entries, entries);
  135. return 0;
  136. }
  137. void CreateRooFitAndDraw(TH1D *hist, double selected_entries, double total_entries)
  138. {
  139. if (hist->GetEntries() == 0) {
  140. std::cout << "### Histogram for RooFit has not entries." << std::endl;
  141. return;
  142. }
  143. // std::cout << "### Entries: " << hist->GetEntries() << std::endl;
  144. RooRealVar roo_var_mass("var_mass", "B Mass Variable", MASS_HIST_MIN, MASS_HIST_MAX);
  145. roo_var_mass.setRange("fitting_range", MASS_HIST_MIN, MASS_HIST_MAX);
  146. RooDataHist roohist_B_M("roohist_B_M", "B Mass Histogram", roo_var_mass, RooFit::Import(*hist));
  147. RooPlot *roo_frame_mass = roo_var_mass.frame(RooFit::Title(quickformat_name("Simulation of %s")));
  148. roohist_B_M.plotOn(roo_frame_mass, RooFit::Binning(nBins), RooFit::Name("B Mass Distribution"));
  149. roo_frame_mass->GetXaxis()->SetTitle(X_AXIS);
  150. // Crystal Ball for Signal
  151. RooRealVar var_mass_x0("var_mass_x0", "#mu", 5278., 5170., 5500.);
  152. RooRealVar var_mass_sigmaLR("var_mass_sigmaLR", "#sigma_{LR}", 16., 5., 25.);
  153. // Same Variables for Left and Right Tail
  154. RooRealVar var_mass_alphaL("var_mass_alphaL", "#alpha_{L}", 2., 0., 4.);
  155. RooRealVar var_mass_nL("var_mass_nL", "n_{L}", 5., 0., 15.);
  156. RooRealVar var_mass_alphaR("var_mass_alphaR", "#alpha_{R}", 2., 0., 4.);
  157. RooRealVar var_mass_nR("var_mass_nR", "n_{R}", 5., 0., 15.);
  158. 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);
  159. // Exponential for Background
  160. // RooRealVar var_mass_bkg_c("var_mass_bkg_c", "Background C", -0.00231049, -0.003, -0.002);
  161. // RooExponential bkg_exp("bkg_exp", "Exp Background", var_mass, var_mass_bkg_c);
  162. // RooRealVar var_mass_nsig(TString::Format("var_%s_nsig", var_id.c_str()), "Mass N Signal", 0., hist->GetEntries());
  163. // RooRealVar var_mass_nbkg(TString::Format("var_%s_nbkg", var_id.c_str()), "Mass N Background", 0., hist->GetEntries());
  164. // TString pdf_name = TString::Format("%s_sigplusbkg", var_id.c_str());
  165. // RooAddPdf sigplusbkg(pdf_name, "Sig and Bkg PDF", RooArgList(sig_cb, bkg_exp), RooArgList(var_mass_nsig, var_mass_nbkg));
  166. RooFitResult *fitres = sig_cb.fitTo(roohist_B_M, RooFit::Save(), RooFit::PrintLevel(1), RooFit::Range("fitting_range"));
  167. sig_cb.plotOn(roo_frame_mass, RooFit::VisualizeError(*fitres, 1), RooFit::FillColor(kOrange + 1), RooFit::FillStyle(3144));
  168. sig_cb.plotOn(roo_frame_mass, RooFit::LineColor(kRed), RooFit::LineStyle(kSolid), RooFit::Range("fitting_range"), RooFit::Name("sig_cb"));
  169. // sigplusbkg.plotOn(roo_frame_mass, RooFit::Components(RooArgSet(bkg_exp)), RooFit::LineColor(kBlue), RooFit::LineStyle(kDashed), RooFit::Range(5170., 5500.));
  170. TCanvas *c = new TCanvas("roofit_c", "roofit_c", 0, 0, 800, 600);
  171. roo_frame_mass->Draw();
  172. TLegend *leg1 = new TLegend(0.58, 0.50, 0.87, 0.87);
  173. // leg1->SetFillColor(kWhite);
  174. leg1->SetLineColor(kWhite);
  175. // leg1->AddEntry(roo_frame_mass->findObject(pdf_name.c_str()), "Signal + Background", "LP");
  176. leg1->AddEntry(roo_frame_mass->findObject("sig_cb"), "Signal", "LP");
  177. // leg1->AddEntry(roo_frame_mass->findObject(name_fit_func_bkg.c_str()), "Background", "LP");
  178. leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", var_mass_x0.getTitle().Data(), var_mass_x0.getVal(), var_mass_x0.getError()).Data(), "");
  179. leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", var_mass_sigmaLR.getTitle().Data(), var_mass_sigmaLR.getVal(), var_mass_sigmaLR.getError()).Data(), "");
  180. // leg1->AddEntry((TObject *)0, TString::Format("%s = %.8f #pm %.8f", roo_sig_tail.getTitle().Data(), roo_sig_tail.getVal(), roo_sig_tail.getError()).Data(), "");
  181. leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", var_mass_alphaL.getTitle().Data(), var_mass_alphaL.getVal(), var_mass_alphaL.getError()).Data(), "");
  182. leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", var_mass_nL.getTitle().Data(), var_mass_nL.getVal(), var_mass_nL.getError()).Data(), "");
  183. leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", var_mass_alphaR.getTitle().Data(), var_mass_alphaR.getVal(), var_mass_alphaR.getError()).Data(), "");
  184. leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", var_mass_nR.getTitle().Data(), var_mass_nR.getVal(), var_mass_nR.getError()).Data(), "");
  185. leg1->AddEntry((TObject *)0, TString::Format("#epsilon = %.2f%%", (selected_entries / total_entries) * 100).Data(), "");
  186. leg1->Draw();
  187. c->Draw();
  188. c->SaveAs(quickformat_name("images/sim/%s/roof_B_M.pdf"));
  189. }