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.

385 lines
15 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 "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 = 5100.;
  34. const Double_t MASS_HIST_MAX = 6000.;
  35. const Double_t MASS_HIST_FIT_MIN = 5100.;
  36. const Double_t MASS_HIST_FIT_MAX = 6000.;
  37. const Double_t K_MASS = 493.677;
  38. const Double_t PSI2S_MASS = 3686.093;
  39. const Double_t JPSI_MASS = 3096.900;
  40. const Double_t KSTAR_MASS = 891.760;
  41. void CreateRooFitAndDraw(TH1D *hist, const char *name);
  42. const char *TITLE = "SpruceRD_BuToHpMuMu (#pi^{+} #rightarrow K^{+}, w/o excl Hlt2 decision cut)";
  43. const char *FILE_NAME = "SpruceRD_BuToHpMuMu_Pip2Kp_wohlt2cut";
  44. const char *MASS_LITERAL = "m(#pi^{+}_{(#rightarrow K^{+})}#mu^{+}#mu^{-})";
  45. struct Hlt1Decision
  46. {
  47. std::string name;
  48. int index;
  49. int flag;
  50. Bool_t value;
  51. std::string GetName() const
  52. {
  53. return TString::Format("Hlt1%sDecision", name.c_str()).Data();
  54. }
  55. Bool_t *GetValuePointer()
  56. {
  57. return &value;
  58. }
  59. };
  60. const int FLAG_TrackMVA = 0b0000100000000000;
  61. const int FLAG_TwoTrackMVA = 0b0001000000000000;
  62. std::vector<Hlt1Decision> Hlt1Decisions{
  63. Hlt1Decision{"DiMuonHighMass", 1, 0b0000000000000001}, // 0001
  64. Hlt1Decision{"DiMuonLowMass", 2, 0b0000000000000010}, // 0002
  65. Hlt1Decision{"DiMuonNoIP", 3, 0b0000000000000100}, // 0004
  66. Hlt1Decision{"DiMuonSoft", 4, 0b0000000000001000}, // 0008
  67. Hlt1Decision{"DisplacedLeptons", 5, 0b0000000000010000}, // 0016
  68. Hlt1Decision{"LowPtDiMuon", 6, 0b0000000000100000}, // 0032
  69. Hlt1Decision{"LowPtMuon", 7, 0b0000000001000000}, // 0064
  70. Hlt1Decision{"OneMuonTrackLine", 8, 0b0000000010000000}, // 0128
  71. Hlt1Decision{"SingleHighEt", 9, 0b0000000100000000}, // 0256
  72. Hlt1Decision{"SingleHighPtMuon", 10, 0b0000001000000000}, // 0512
  73. Hlt1Decision{"TrackMVA", 11, FLAG_TrackMVA}, // 1024
  74. Hlt1Decision{"TrackMuonMVA", 12, 0b0000100000000000}, // 2048
  75. Hlt1Decision{"TwoTrackMVA", 13, FLAG_TwoTrackMVA}, // 4096
  76. };
  77. // 0b0001000000000000
  78. // 0b0001000000110111
  79. void ConnectHlt1Decisions(TChain *chain)
  80. {
  81. for (auto &var : Hlt1Decisions)
  82. {
  83. chain->SetBranchAddress(var.GetName().c_str(), var.GetValuePointer());
  84. }
  85. }
  86. int CombineHlt1DecisionFlags()
  87. {
  88. int result = 0;
  89. for (const auto &var : Hlt1Decisions)
  90. {
  91. if (var.value)
  92. {
  93. result = result | var.flag;
  94. }
  95. }
  96. return result;
  97. }
  98. int analysis_fullstream_bu2hpmumu()
  99. {
  100. TChain *data_chain = new TChain("BuToHpMuMu/DecayTree");
  101. data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/spruce_magdown_2023_v0r1_tuple_90000000_2023_v0r0p6288631.root");
  102. data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/spruce_magdown_2023_v0_tuple_90000000_v0r0p6288631.root");
  103. Float_t B_BPVFDCHI2,
  104. B_BPVIPCHI2,
  105. L1_BPVIPCHI2,
  106. L2_BPVIPCHI2,
  107. L1_PT,
  108. L2_PT,
  109. Jpsi_BPVFDCHI2,
  110. Hp_PT,
  111. Hp_BPVIPCHI2,
  112. Hp_P;
  113. Double_t L1_PID_MU,
  114. L2_PID_MU,
  115. B_CHI2VXNDOF,
  116. Jpsi_MAXDOCACHI2,
  117. Jpsi_CHI2DOF,
  118. Hp_PID_K,
  119. Jpsi_M,
  120. B_M;
  121. Bool_t L1_ISMUON,
  122. L2_ISMUON,
  123. Hlt2RD_BuToKpMuMuDecision,
  124. Hlt2_InclDetDiMuon_4BodyDecision,
  125. Hlt2_InclDetDiMuon_3BodyDecision,
  126. Hlt2_InclDetDiMuonDecision,
  127. Hlt1TrackMVADecision,
  128. Hlt1TwoTrackMVADecision;
  129. data_chain->SetBranchAddress("B_M", &B_M);
  130. data_chain->SetBranchAddress("B_BPVFDCHI2", &B_BPVFDCHI2);
  131. data_chain->SetBranchAddress("B_BPVIPCHI2", &B_BPVIPCHI2);
  132. data_chain->SetBranchAddress("L1_BPVIPCHI2", &L1_BPVIPCHI2);
  133. data_chain->SetBranchAddress("L2_BPVIPCHI2", &L2_BPVIPCHI2);
  134. data_chain->SetBranchAddress("L1_PID_MU", &L1_PID_MU);
  135. data_chain->SetBranchAddress("L2_PID_MU", &L2_PID_MU);
  136. data_chain->SetBranchAddress("L1_ISMUON", &L1_ISMUON);
  137. data_chain->SetBranchAddress("L2_ISMUON", &L2_ISMUON);
  138. data_chain->SetBranchAddress("L1_PT", &L1_PT);
  139. data_chain->SetBranchAddress("L2_PT", &L2_PT);
  140. data_chain->SetBranchAddress("B_CHI2VXNDOF", &B_CHI2VXNDOF);
  141. data_chain->SetBranchAddress("Jpsi_MAXDOCACHI2", &Jpsi_MAXDOCACHI2);
  142. data_chain->SetBranchAddress("Jpsi_CHI2DOF", &Jpsi_CHI2DOF);
  143. data_chain->SetBranchAddress("Hp_PT", &Hp_PT);
  144. data_chain->SetBranchAddress("Hp_BPVIPCHI2", &Hp_BPVIPCHI2);
  145. data_chain->SetBranchAddress("Hp_P", &Hp_P);
  146. data_chain->SetBranchAddress("Hp_PID_K", &Hp_PID_K);
  147. data_chain->SetBranchAddress("Jpsi_BPVFDCHI2", &Jpsi_BPVFDCHI2);
  148. data_chain->SetBranchAddress("Jpsi_M", &Jpsi_M);
  149. data_chain->SetBranchAddress("Hlt2RD_BuToKpMuMuDecision", &Hlt2RD_BuToKpMuMuDecision);
  150. data_chain->SetBranchAddress("Hlt2_InclDetDiMuon_4BodyDecision", &Hlt2_InclDetDiMuon_4BodyDecision);
  151. data_chain->SetBranchAddress("Hlt2_InclDetDiMuon_3BodyDecision", &Hlt2_InclDetDiMuon_3BodyDecision);
  152. data_chain->SetBranchAddress("Hlt2_InclDetDiMuonDecision", &Hlt2_InclDetDiMuonDecision);
  153. // data_chain->SetBranchAddress("Hlt1TrackMVADecision", &Hlt1TrackMVADecision);
  154. // data_chain->SetBranchAddress("Hlt1TwoTrackMVADecision", &Hlt1TwoTrackMVADecision);
  155. Float_t L1_PX, L1_PY, L1_PZ, L1_ENERGY, L2_PX, L2_PY, L2_PZ, L2_ENERGY, Hp_PX, Hp_PY, Hp_PZ, Hp_ENERGY;
  156. data_chain->SetBranchAddress("L1_PX", &L1_PX);
  157. data_chain->SetBranchAddress("L1_PY", &L1_PY);
  158. data_chain->SetBranchAddress("L1_PZ", &L1_PZ);
  159. data_chain->SetBranchAddress("L1_ENERGY", &L1_ENERGY);
  160. data_chain->SetBranchAddress("L2_PX", &L2_PX);
  161. data_chain->SetBranchAddress("L2_PY", &L2_PY);
  162. data_chain->SetBranchAddress("L2_PZ", &L2_PZ);
  163. data_chain->SetBranchAddress("L2_ENERGY", &L2_ENERGY);
  164. data_chain->SetBranchAddress("Hp_PX", &Hp_PX);
  165. data_chain->SetBranchAddress("Hp_PY", &Hp_PY);
  166. data_chain->SetBranchAddress("Hp_PZ", &Hp_PZ);
  167. data_chain->SetBranchAddress("Hp_ENERGY", &Hp_ENERGY);
  168. ConnectHlt1Decisions(data_chain);
  169. TH1D *h1_B_M_jpsi = new TH1D("h1_B_M_jpsi", "B Mass, J/#psi Mode", nBins, MASS_HIST_MIN, MASS_HIST_MAX);
  170. TH1D *h1_B_M_psi2s = new TH1D("h1_B_M_psi2s", "B Mass, #psi(2S) Mode", nBins, MASS_HIST_MIN, MASS_HIST_MAX);
  171. TH2D *h2_Hlt1_flags_B_Mass = new TH2D("h2_Hlt1_flags_B_Mass", "Hlt1 Decision vs B Mass", 50, 5100, 5400, 13, 1., 14.);
  172. h2_Hlt1_flags_B_Mass->SetCanExtend(TH1::kYaxis);
  173. for (const auto &var : Hlt1Decisions)
  174. {
  175. h2_Hlt1_flags_B_Mass->Fill(0., var.name.c_str(), 0.);
  176. }
  177. unsigned int entries = data_chain->GetEntries();
  178. unsigned int fitting_entries = 0;
  179. for (unsigned int i = 0; i < entries; i++)
  180. {
  181. data_chain->GetEntry(i);
  182. TVector3 K_momentum(Hp_PX, Hp_PY, Hp_PZ);
  183. double K_energy = TMath::Sqrt(TMath::Sq(K_MASS) + K_momentum.Mag2());
  184. TLorentzVector K_4v(K_momentum, K_energy);
  185. TLorentzVector l1_4v(L1_PX, L1_PY, L1_PZ, L1_ENERGY);
  186. TLorentzVector l2_4v(L2_PX, L2_PY, L2_PZ, L2_ENERGY);
  187. Double_t reconstructed_B_Mass = (K_4v + l1_4v + l2_4v).M();
  188. int combinedFlags = CombineHlt1DecisionFlags();
  189. bool trackMVADec = (((combinedFlags & FLAG_TrackMVA) == FLAG_TrackMVA) | ((combinedFlags & FLAG_TwoTrackMVA) == FLAG_TwoTrackMVA));
  190. if (TMath::Abs(Jpsi_M - JPSI_MASS) < 100 && trackMVADec)
  191. {
  192. for (const auto &var : Hlt1Decisions)
  193. {
  194. if (var.value)
  195. {
  196. h2_Hlt1_flags_B_Mass->Fill(reconstructed_B_Mass, var.name.c_str(), 1);
  197. }
  198. }
  199. }
  200. // std::cout << "HLT1 Decision Flags: " << combinedFlags << std::endl;
  201. if (trackMVADec)
  202. {
  203. if (TMath::Abs(Jpsi_M - JPSI_MASS) < 100)
  204. {
  205. h1_B_M_jpsi->Fill(reconstructed_B_Mass);
  206. }
  207. else if (TMath::Abs(Jpsi_M - PSI2S_MASS) < 100)
  208. {
  209. h1_B_M_psi2s->Fill(reconstructed_B_Mass);
  210. }
  211. }
  212. // if ((((B_BPVFDCHI2 > 36) & (B_CHI2VXNDOF < 16) & (B_BPVIPCHI2 < 25) & (Jpsi_MAXDOCACHI2 < 36) & (L1_BPVIPCHI2 > 9) & (L2_BPVIPCHI2 > 9) & (L1_PID_MU > -3) & (L2_PID_MU > -3) & (L1_ISMUON == 1) & (L2_ISMUON == 1) & (L1_PT > 350) & (L2_PT > 350) & (Jpsi_CHI2DOF < 9) & (Jpsi_BPVFDCHI2 > 16) & (Jpsi_M < 5500) & (Hp_PT > 400) & (Hp_BPVIPCHI2 > 6) & (Hp_PT > 400) & (Hp_BPVIPCHI2 > 6) & (Hp_P > 2000) & (Hp_PID_K > -4)) /*& (Hlt2RD_BuToKpMuMuDecision == 1)*/ & (((Hlt2_InclDetDiMuon_4BodyDecision == 1) | (Hlt2_InclDetDiMuon_3BodyDecision == 1) | (Hlt2_InclDetDiMuonDecision == 1)))))
  213. // {
  214. // if ((TMath::Abs(Jpsi_M - 3096.9) < 100) & ((Hlt1TrackMVADecision == 1) | (Hlt1TwoTrackMVADecision == 1)))
  215. // {
  216. // h1_B_M->Fill(reconstructed_B_Mass);
  217. // if (MASS_HIST_FIT_MIN <= reconstructed_B_Mass && reconstructed_B_Mass <= MASS_HIST_FIT_MAX)
  218. // {
  219. // fitting_entries++;
  220. // }
  221. // }
  222. // }
  223. if ((i + 1) % 10000 == 0 || i + 1 == entries)
  224. {
  225. std::cout << "["
  226. << FILE_NAME
  227. << "] Processed event: " << i + 1 << " (" << TString::Format("%.2f", ((double)(i + 1) / (double)entries) * 100.) << "%)" << std::endl;
  228. }
  229. }
  230. TCanvas *c1 = new TCanvas("c1", "c1", 0, 0, 800, 600);
  231. h1_B_M_jpsi->Draw();
  232. c1->Draw();
  233. c1->SaveAs(TString::Format("images/root_hist_%s_%s_bmass.pdf", FILE_NAME, "jpsi").Data());
  234. TCanvas *c2 = new TCanvas("c2", "c2", 0, 0, 800, 600);
  235. h1_B_M_psi2s->Draw();
  236. c2->Draw();
  237. c2->SaveAs(TString::Format("images/root_hist_%s_%s_bmass.pdf", FILE_NAME, "psi2s").Data());
  238. TCanvas *c3 = new TCanvas("c3", "c3", 0, 0, 800, 600);
  239. c3->SetLeftMargin(0.20);
  240. h2_Hlt1_flags_B_Mass->SetStats(0);
  241. h2_Hlt1_flags_B_Mass->Draw("COLZ");
  242. c3->Draw();
  243. // CreateRooFitAndDraw(h1_B_M_jpsi, "jpsi");
  244. // CreateRooFitAndDraw(h1_B_M_psi2s, "psi2s");
  245. return 0;
  246. }
  247. void CreateRooFitAndDraw(TH1D *hist, const char *name)
  248. {
  249. auto suffix = [name](const char *text)
  250. {
  251. return TString::Format("%s_%s", text, name);
  252. };
  253. std::cout << " ###### " << suffix("AABB") << std::endl;
  254. RooRealVar roo_var_mass(suffix("var_mass"), "B Mass Variable", MASS_HIST_FIT_MIN, MASS_HIST_FIT_MAX);
  255. roo_var_mass.setRange("fitting_range", MASS_HIST_FIT_MIN, MASS_HIST_FIT_MAX);
  256. TString hist_name = suffix("roohist_B_M");
  257. RooDataHist roohist_B_M(hist_name, "B Mass Histogram", roo_var_mass, RooFit::Import(*hist));
  258. RooPlot *roo_frame_mass = roo_var_mass.frame(RooFit::Title(TITLE));
  259. roohist_B_M.plotOn(roo_frame_mass, RooFit::Binning(nBins), RooFit::Name(hist_name));
  260. roo_frame_mass->GetXaxis()->SetTitle(MASS_LITERAL);
  261. // Crystal Ball for Signal
  262. RooRealVar var_mass_x0(suffix("var_mass_x0"), "#mu", 5278., 5170., 5500.);
  263. RooRealVar var_mass_sigmaLR(suffix("var_mass_sigmaLR"), "#sigma_{LR}", 16., 5., 40.);
  264. // Same Variables for Left and Right Tail
  265. // RooRealVar var_mass_alphaL("var_mass_alphaL", "#alpha_{L}", 2., 0., 4.);
  266. // RooRealVar var_mass_nL("var_mass_nL", "n_{L}", 5., 0., 15.);
  267. // RooRealVar var_mass_alphaR("var_mass_alphaR", "#alpha_{R}", 2., 0., 4.);
  268. // RooRealVar var_mass_nR("var_mass_nR", "n_{R}", 5., 0., 15.);
  269. auto var_mass_alphaL = RooRealConstant::value(1.894);
  270. auto var_mass_nL = RooRealConstant::value(1.120);
  271. auto var_mass_alphaR = RooRealConstant::value(2.446);
  272. auto var_mass_nR = RooRealConstant::value(1.146);
  273. RooCrystalBall sig_cb(suffix("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);
  274. // Exponential for Background
  275. RooRealVar var_mass_bkg_c(suffix("var_mass_bkg_c"), "#lambda", -0.0014, -0.004, -0.000);
  276. RooExponential bkg_exp(suffix("bkg_exp"), "Exp Background", roo_var_mass, var_mass_bkg_c);
  277. RooRealVar var_mass_nsig(suffix("nsig"), "Mass N Signal", 0., hist->GetEntries());
  278. RooRealVar var_mass_nbkg(suffix("nbkg"), "Mass N Background", 0., hist->GetEntries());
  279. TString pdf_name = suffix("sigplusbkg");
  280. RooAddPdf sigplusbkg(pdf_name, "Sig and Bkg PDF", RooArgList(sig_cb, bkg_exp), RooArgList(var_mass_nsig, var_mass_nbkg));
  281. RooFitResult *fitres = sigplusbkg.fitTo(roohist_B_M, RooFit::Save(), RooFit::PrintLevel(1), RooFit::Range("fitting_range"));
  282. // sigplusbkg.plotOn(roo_frame_mass, RooFit::VisualizeError(*fitres, 1), RooFit::FillColor(kOrange + 1), RooFit::FillStyle(3144));
  283. sigplusbkg.plotOn(roo_frame_mass, RooFit::LineColor(kRed), RooFit::LineStyle(kSolid), RooFit::Range("fitting_range"), RooFit::Name(pdf_name));
  284. 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"));
  285. 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"));
  286. RooPlot *roo_frame_pull = roo_var_mass.frame(RooFit::Title("Pull Distribution"));
  287. roo_frame_pull->addPlotable(roo_frame_mass->pullHist(hist_name, pdf_name), "P");
  288. TCanvas *c = new TCanvas(suffix("roofit_c"), "roofit_c", 0, 0, 800, 600);
  289. auto *p2 = new TPad(suffix("p2"), "Pull", 0., 0., 1., 0.3);
  290. p2->Draw();
  291. p2->SetTopMargin(0.001);
  292. p2->SetBottomMargin(0.3);
  293. p2->SetGrid();
  294. auto *p1 = new TPad(suffix("p1"), "Fit", 0., 0.32, 1., 1.);
  295. p1->Draw();
  296. p1->SetBottomMargin(0.001);
  297. p1->cd();
  298. roo_frame_mass->Draw();
  299. TLegend *leg1 = new TLegend(0.50, 0.50, 0.87, 0.89);
  300. leg1->SetFillColor(kWhite);
  301. leg1->SetLineColor(kBlack);
  302. leg1->AddEntry(roo_frame_mass->findObject(pdf_name), "Signal + Background", "LP");
  303. leg1->AddEntry(roo_frame_mass->findObject("sig_cb"), "Signal", "LP");
  304. leg1->AddEntry(roo_frame_mass->findObject("bkg_exp"), "Background", "LP");
  305. leg1->AddEntry((TObject *)0, "", "");
  306. leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", var_mass_x0.getTitle().Data(), var_mass_x0.getVal(), var_mass_x0.getError()).Data(), "");
  307. leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", var_mass_sigmaLR.getTitle().Data(), var_mass_sigmaLR.getVal(), var_mass_sigmaLR.getError()).Data(), "");
  308. 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(), "");
  309. leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", "S", var_mass_nsig.getVal(), var_mass_nsig.getError()).Data(), "");
  310. leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", "B", var_mass_nbkg.getVal(), var_mass_nbkg.getError()).Data(), "");
  311. leg1->AddEntry((TObject *)0, TString::Format("Entries: %.1f", var_mass_nsig.getVal() + var_mass_nbkg.getVal()), "");
  312. leg1->Draw();
  313. p2->cd();
  314. roo_frame_pull->Draw();
  315. c->Draw();
  316. c->SaveAs(TString::Format("images/data/roofit_hist_%s_bmass.pdf", suffix(FILE_NAME).Data()).Data());
  317. }