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.

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