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.

332 lines
14 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 "TLegend.h"
  15. #include "TLatex.h"
  16. #include "RooDataHist.h"
  17. #include "RooRealVar.h"
  18. #include "RooPlot.h"
  19. #include "RooGaussian.h"
  20. #include "RooExponential.h"
  21. #include "RooRealConstant.h"
  22. #include "RooAddPdf.h"
  23. #include "RooFitResult.h"
  24. #include "RooProduct.h"
  25. #include "RooCrystalBall.h"
  26. #include "RooBreitWigner.h"
  27. #include "RooArgSet.h"
  28. #include "RooFFTConvPdf.h"
  29. #include "RooNovosibirsk.h"
  30. #include <string>
  31. #include <iostream>
  32. #include <cmath>
  33. const int nBins = 70;
  34. const Double_t MASS_HIST_MIN = 5100.;
  35. const Double_t MASS_HIST_MAX = 5700.;
  36. const Double_t MASS_HIST_FIT_MIN = 5100.;
  37. const Double_t MASS_HIST_FIT_MAX = 5700.;
  38. const char *TITLE = "SpruceRD_B0ToHpHmMuMu (#pi^{-} #rightarrow K^{-})";
  39. const char *FILE_NAME = "SpruceRD_B0ToHpHmMuMu_Pip2Kp";
  40. const char *MASS_LITERAL = "m(#pi^{+}#pi^{-}_{(#rightarrow K^{-})}#mu^{+}#mu^{-})";
  41. const Double_t K_MASS = 493.677;
  42. const Double_t PSI2S_MASS = 3686.093;
  43. const Double_t JPSI_MASS = 3096.9;
  44. const Double_t KSTAR_MASS = 891.76;
  45. void CreateRooFitAndDraw(TH1D *hist);
  46. int analysis_fullstream_b02hphmmumu()
  47. {
  48. TChain *data_chain = new TChain("B0ToHpHmMuMu/DecayTree");
  49. data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/spruce_magdown_2023_v0r1_tuple_90000000_2023_v0r0p6288631.root");
  50. data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/spruce_magdown_2023_v0_tuple_90000000_v0r0p6288631.root");
  51. Float_t B0_BPVFDCHI2,
  52. B0_BPVIPCHI2,
  53. L1_BPVIPCHI2,
  54. L2_BPVIPCHI2,
  55. L1_PT,
  56. L2_PT,
  57. Jpsi_BPVFDCHI2,
  58. Hp_PT,
  59. Hp_BPVIPCHI2,
  60. Hp_P,
  61. Hm_PT,
  62. Hm_BPVIPCHI2,
  63. Hm_P,
  64. Res_PT,
  65. Res_BPVIPCHI2,
  66. Res_P;
  67. Double_t L1_PID_MU,
  68. L2_PID_MU,
  69. B0_CHI2VXNDOF,
  70. Jpsi_MAXDOCACHI2,
  71. Jpsi_CHI2DOF,
  72. Res_CHI2DOF,
  73. Hp_PID_K,
  74. Hm_PID_K,
  75. Jpsi_M,
  76. B0_M,
  77. Res_MAXDOCACHI2,
  78. Res_M;
  79. Bool_t L1_ISMUON,
  80. L2_ISMUON,
  81. Hlt2RD_B0ToKpPimMuMuDecision,
  82. Hlt2_InclDetDiMuon_4BodyDecision,
  83. Hlt2_InclDetDiMuon_3BodyDecision,
  84. Hlt2_InclDetDiMuonDecision,
  85. Hlt1TrackMVADecision,
  86. Hlt1TwoTrackMVADecision;
  87. data_chain->SetBranchAddress("B0_M", &B0_M);
  88. data_chain->SetBranchAddress("B0_BPVFDCHI2", &B0_BPVFDCHI2);
  89. data_chain->SetBranchAddress("B0_BPVIPCHI2", &B0_BPVIPCHI2);
  90. data_chain->SetBranchAddress("L1_BPVIPCHI2", &L1_BPVIPCHI2);
  91. data_chain->SetBranchAddress("L2_BPVIPCHI2", &L2_BPVIPCHI2);
  92. data_chain->SetBranchAddress("L1_PID_MU", &L1_PID_MU);
  93. data_chain->SetBranchAddress("L2_PID_MU", &L2_PID_MU);
  94. data_chain->SetBranchAddress("L1_ISMUON", &L1_ISMUON);
  95. data_chain->SetBranchAddress("L2_ISMUON", &L2_ISMUON);
  96. data_chain->SetBranchAddress("L1_PT", &L1_PT);
  97. data_chain->SetBranchAddress("L2_PT", &L2_PT);
  98. data_chain->SetBranchAddress("B0_CHI2VXNDOF", &B0_CHI2VXNDOF);
  99. data_chain->SetBranchAddress("Jpsi_MAXDOCACHI2", &Jpsi_MAXDOCACHI2);
  100. data_chain->SetBranchAddress("Jpsi_CHI2DOF", &Jpsi_CHI2DOF);
  101. data_chain->SetBranchAddress("Res_CHI2DOF", &Res_CHI2DOF);
  102. data_chain->SetBranchAddress("Hp_PT", &Hp_PT);
  103. data_chain->SetBranchAddress("Hp_BPVIPCHI2", &Hp_BPVIPCHI2);
  104. data_chain->SetBranchAddress("Hp_P", &Hp_P);
  105. data_chain->SetBranchAddress("Hp_PID_K", &Hp_PID_K);
  106. data_chain->SetBranchAddress("Hm_PT", &Hm_PT);
  107. data_chain->SetBranchAddress("Hm_BPVIPCHI2", &Hm_BPVIPCHI2);
  108. data_chain->SetBranchAddress("Hm_P", &Hm_P);
  109. data_chain->SetBranchAddress("Hm_PID_K", &Hm_PID_K);
  110. data_chain->SetBranchAddress("Res_PT", &Res_PT);
  111. data_chain->SetBranchAddress("Res_BPVIPCHI2", &Res_BPVIPCHI2);
  112. data_chain->SetBranchAddress("Res_P", &Res_P);
  113. data_chain->SetBranchAddress("Res_MAXDOCACHI2", &Res_MAXDOCACHI2);
  114. data_chain->SetBranchAddress("Res_M", &Res_M);
  115. data_chain->SetBranchAddress("Jpsi_BPVFDCHI2", &Jpsi_BPVFDCHI2);
  116. data_chain->SetBranchAddress("Jpsi_M", &Jpsi_M);
  117. data_chain->SetBranchAddress("Hlt2RD_B0ToKpPimMuMuDecision", &Hlt2RD_B0ToKpPimMuMuDecision);
  118. data_chain->SetBranchAddress("Hlt2_InclDetDiMuon_4BodyDecision", &Hlt2_InclDetDiMuon_4BodyDecision);
  119. data_chain->SetBranchAddress("Hlt2_InclDetDiMuon_3BodyDecision", &Hlt2_InclDetDiMuon_3BodyDecision);
  120. data_chain->SetBranchAddress("Hlt2_InclDetDiMuonDecision", &Hlt2_InclDetDiMuonDecision);
  121. data_chain->SetBranchAddress("Hlt1TrackMVADecision", &Hlt1TrackMVADecision);
  122. data_chain->SetBranchAddress("Hlt1TwoTrackMVADecision", &Hlt1TwoTrackMVADecision);
  123. Float_t
  124. L1_PX,
  125. L1_PY, L1_PZ, L1_ENERGY,
  126. L2_PX, L2_PY, L2_PZ, L2_ENERGY,
  127. Hp_PX, Hp_PY, Hp_PZ, Hp_ENERGY,
  128. Hm_PX, Hm_PY, Hm_PZ, Hm_ENERGY;
  129. data_chain->SetBranchAddress("L1_PX", &L1_PX);
  130. data_chain->SetBranchAddress("L1_PY", &L1_PY);
  131. data_chain->SetBranchAddress("L1_PZ", &L1_PZ);
  132. data_chain->SetBranchAddress("L1_ENERGY", &L1_ENERGY);
  133. data_chain->SetBranchAddress("L2_PX", &L2_PX);
  134. data_chain->SetBranchAddress("L2_PY", &L2_PY);
  135. data_chain->SetBranchAddress("L2_PZ", &L2_PZ);
  136. data_chain->SetBranchAddress("L2_ENERGY", &L2_ENERGY);
  137. data_chain->SetBranchAddress("Hp_PX", &Hp_PX);
  138. data_chain->SetBranchAddress("Hp_PY", &Hp_PY);
  139. data_chain->SetBranchAddress("Hp_PZ", &Hp_PZ);
  140. data_chain->SetBranchAddress("Hp_ENERGY", &Hp_ENERGY);
  141. data_chain->SetBranchAddress("Hm_PX", &Hm_PX);
  142. data_chain->SetBranchAddress("Hm_PY", &Hm_PY);
  143. data_chain->SetBranchAddress("Hm_PZ", &Hm_PZ);
  144. data_chain->SetBranchAddress("Hm_ENERGY", &Hm_ENERGY);
  145. TH1D *h1_B0_M = new TH1D("h1_B0_M", TITLE, nBins, MASS_HIST_MIN, MASS_HIST_MAX);
  146. TH1D *h1_Res_M = new TH1D("h1_Res_M", "K^{*0} Mass", nBins, 0., 3500.);
  147. TH1D *h1_Hm_PID_K = new TH1D("h1_Hm_PID_K", "H^{-} PID K", nBins, -10., 10.);
  148. TH1D *h1_Hp_PID_K = new TH1D("h1_Hp_PID_K", "H^{+} PID K", nBins, -10., 10.);
  149. unsigned int entries = data_chain->GetEntries();
  150. unsigned int fitting_entries = 0;
  151. for (unsigned int i = 0; i < entries; i++)
  152. {
  153. data_chain->GetEntry(i);
  154. TLorentzVector l1_4v(L1_PX, L1_PY, L1_PZ, L1_ENERGY);
  155. TLorentzVector l2_4v(L2_PX, L2_PY, L2_PZ, L2_ENERGY);
  156. // Pi+ -> K+
  157. TVector3 Kp_momentum(Hp_PX, Hp_PY, Hp_PZ);
  158. double Kp_energy = TMath::Sqrt(TMath::Sq(K_MASS) + Kp_momentum.Mag2());
  159. TLorentzVector Pim_4v(Hm_PX, Hm_PY, Hm_PZ, Hm_ENERGY);
  160. TLorentzVector Kp_4v(Kp_momentum, Kp_energy);
  161. Double_t reconstructed_Kp_Pim_Mass = (Kp_4v + Pim_4v).M();
  162. Double_t reconstructed_B0_Kp_Pim_Mass = (Kp_4v + Pim_4v + l1_4v + l2_4v).M();
  163. // Pi- -> K-
  164. TVector3 Km_momentum(Hm_PX, Hm_PY, Hm_PZ);
  165. double Km_energy = TMath::Sqrt(TMath::Sq(K_MASS) + Km_momentum.Mag2());
  166. TLorentzVector Pip_4v(Hp_PX, Hp_PY, Hp_PZ, Hp_ENERGY);
  167. TLorentzVector Km_4v(Km_momentum, Km_energy);
  168. Double_t reconstructed_Km_Pip_Mass = (Km_4v + Pip_4v).M();
  169. Double_t reconstructed_B0_Km_Pip_Mass = (Km_4v + Pip_4v + l1_4v + l2_4v).M();
  170. h1_Res_M->Fill(reconstructed_Kp_Pim_Mass);
  171. h1_Hm_PID_K->Fill(Hm_PID_K);
  172. h1_Hp_PID_K->Fill(Hp_PID_K);
  173. if ((
  174. (B0_BPVFDCHI2 > 64) & (B0_CHI2VXNDOF < 9) & (B0_BPVIPCHI2 < 25) &
  175. //
  176. (Jpsi_MAXDOCACHI2 < 36) &
  177. //
  178. (L1_BPVIPCHI2 > 9) & (L2_BPVIPCHI2 > 9) & (L1_PID_MU > -3) & (L2_PID_MU > -3) & (L1_ISMUON) & (L2_ISMUON) & (L1_PT > 350) & (L2_PT > 350) &
  179. //
  180. (Jpsi_CHI2DOF < 9) & (Jpsi_BPVFDCHI2 > 16) & (Jpsi_M < 5500) &
  181. //
  182. (Res_PT > 400) & (Res_MAXDOCACHI2 < 36) & (reconstructed_Kp_Pim_Mass < 2600) & (reconstructed_Km_Pip_Mass < 2600) & (Res_CHI2DOF < 25) & (Res_BPVIPCHI2 > 4) &
  183. //
  184. (Hm_BPVIPCHI2 > 6) & (Hm_PT > 250) & (Hm_P > 2000) & (Hm_PID_K > 0) &
  185. //
  186. (Hp_BPVIPCHI2 > 6) & (Hp_PT > 250) & (Hp_P > 2000) & (Hp_PID_K < 0)
  187. ) &
  188. //
  189. (Hlt2RD_B0ToKpPimMuMuDecision) & (((Hlt2_InclDetDiMuon_4BodyDecision) | (Hlt2_InclDetDiMuon_3BodyDecision) | (Hlt2_InclDetDiMuonDecision))))
  190. {
  191. if ((TMath::Abs(Jpsi_M - JPSI_MASS) < 100) & ((Hlt1TrackMVADecision) | (Hlt1TwoTrackMVADecision)))
  192. {
  193. h1_B0_M->Fill(reconstructed_B0_Km_Pip_Mass);
  194. }
  195. }
  196. if ((i + 1) % 10000 == 0 || i + 1 == entries)
  197. {
  198. std::cout << "["
  199. << FILE_NAME
  200. << "] Processed event: " << i + 1 << " (" << TString::Format("%.2f", ((double)(i + 1) / (double)entries) * 100.) << "%)" << std::endl;
  201. }
  202. }
  203. h1_B0_M->GetXaxis()->SetTitle(MASS_LITERAL);
  204. TCanvas *c1 = new TCanvas("c1", "c1", 0, 0, 800, 600);
  205. h1_B0_M->Draw();
  206. c1->Draw();
  207. c1->SaveAs(TString::Format("images/root_hist_%s_bmass.pdf", FILE_NAME).Data());
  208. TCanvas *c2 = new TCanvas("c2", "c2", 0, 0, 800, 600);
  209. h1_Res_M->Draw();
  210. c2->Draw();
  211. TCanvas *c3 = new TCanvas("c3", "c3", 0, 0, 800, 600);
  212. h1_Hm_PID_K->SetStats(0);
  213. h1_Hm_PID_K->SetLineColor(kRed);
  214. h1_Hp_PID_K->SetStats(0);
  215. h1_Hm_PID_K->Draw();
  216. h1_Hp_PID_K->Draw("SAME");
  217. c3->BuildLegend();
  218. c3->Draw();
  219. CreateRooFitAndDraw(h1_B0_M);
  220. return 0;
  221. }
  222. void CreateRooFitAndDraw(TH1D *hist)
  223. {
  224. RooRealVar roo_var_mass("var_mass", "B0 Mass Variable", MASS_HIST_FIT_MIN, MASS_HIST_FIT_MAX);
  225. roo_var_mass.setRange("fitting_range", MASS_HIST_FIT_MIN, MASS_HIST_FIT_MAX);
  226. RooDataHist roohist_B0_M("roohist_B0_M", "B0 Mass Histogram", roo_var_mass, RooFit::Import(*hist));
  227. RooPlot *roo_frame_mass = roo_var_mass.frame(RooFit::Title(TITLE));
  228. roohist_B0_M.plotOn(roo_frame_mass, RooFit::Binning(nBins), RooFit::Name("B Mass Distribution"));
  229. roo_frame_mass->GetXaxis()->SetTitle(MASS_LITERAL);
  230. // Crystal Ball for Signal
  231. RooRealVar var_mass_x0("var_mass_x0", "#mu", 5278., 5170., 5500.);
  232. RooRealVar var_mass_sigmaLR("var_mass_sigmaLR", "#sigma_{LR}", 16., 5., 40.);
  233. // Same Variables for Left and Right Tail
  234. // RooRealVar var_mass_alphaL("var_mass_alphaL", "#alpha_{L}", 2., 0., 4.);
  235. // RooRealVar var_mass_nL("var_mass_nL", "n_{L}", 5., 0., 15.);
  236. // RooRealVar var_mass_alphaR("var_mass_alphaR", "#alpha_{R}", 2., 0., 4.);
  237. // RooRealVar var_mass_nR("var_mass_nR", "n_{R}", 5., 0., 15.);
  238. auto var_mass_alphaL = RooRealConstant::value(1.915);
  239. auto var_mass_nL = RooRealConstant::value(0.573);
  240. auto var_mass_alphaR = RooRealConstant::value(2.257);
  241. auto var_mass_nR = RooRealConstant::value(0.418);
  242. 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);
  243. // Exponential for Background
  244. RooRealVar var_mass_bkg_c("var_mass_bkg_c", "#lambda", -0.0014, -0.004, -0.000);
  245. RooExponential bkg_exp("bkg_exp", "Exp Background", roo_var_mass, var_mass_bkg_c);
  246. RooRealVar var_mass_nsig("nsig", "Mass N Signal", 0., hist->GetEntries());
  247. RooRealVar var_mass_nbkg("nbkg", "Mass N Background", 0., hist->GetEntries());
  248. // TString pdf_name = TString::Format("%s_sigplusbkg", var_id.c_str());
  249. RooAddPdf sigplusbkg("sigplusbkg", "Sig and Bkg PDF", RooArgList(sig_cb, bkg_exp), RooArgList(var_mass_nsig, var_mass_nbkg));
  250. RooFitResult *fitres = sigplusbkg.fitTo(roohist_B0_M, RooFit::Save(), RooFit::PrintLevel(1), RooFit::Range("fitting_range"));
  251. sigplusbkg.plotOn(roo_frame_mass, RooFit::VisualizeError(*fitres, 1), RooFit::FillColor(kOrange + 1), RooFit::FillStyle(3144));
  252. sigplusbkg.plotOn(roo_frame_mass, RooFit::LineColor(kRed), RooFit::LineStyle(kSolid), RooFit::Range("fitting_range"), RooFit::Name("sigplusbkg"));
  253. 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"));
  254. 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"));
  255. TCanvas *c = new TCanvas("roofit_c", "roofit_c", 0, 0, 800, 600);
  256. roo_frame_mass->Draw();
  257. TLegend *leg1 = new TLegend(0.50, 0.58, 0.87, 0.89);
  258. leg1->SetFillColor(kWhite);
  259. leg1->SetLineColor(kBlack);
  260. leg1->AddEntry(roo_frame_mass->findObject("sigplusbkg"), "Signal + Background", "LP");
  261. leg1->AddEntry(roo_frame_mass->findObject("sig_cb"), "Signal", "LP");
  262. leg1->AddEntry(roo_frame_mass->findObject("bkg_exp"), "Background", "LP");
  263. leg1->AddEntry((TObject *)0, "", "");
  264. leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", var_mass_x0.getTitle().Data(), var_mass_x0.getVal(), var_mass_x0.getError()).Data(), "");
  265. leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", var_mass_sigmaLR.getTitle().Data(), var_mass_sigmaLR.getVal(), var_mass_sigmaLR.getError()).Data(), "");
  266. 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(), "");
  267. leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", "S", var_mass_nsig.getVal(), var_mass_nsig.getError()).Data(), "");
  268. leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", "B", var_mass_nbkg.getVal(), var_mass_nbkg.getError()).Data(), "");
  269. leg1->AddEntry((TObject *)0, TString::Format("Entries: %.1f", var_mass_nsig.getVal() + var_mass_nbkg.getVal()), "");
  270. leg1->Draw();
  271. c->Draw();
  272. c->SaveAs(TString::Format("images/data/roofit_hist_%s_fitted_bmass.pdf", FILE_NAME).Data());
  273. }