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.

623 lines
24 KiB

7 months ago
7 months ago
  1. #include <string>
  2. #include <iostream>
  3. #include <cmath>
  4. #include <algorithm>
  5. #include <filesystem>
  6. #include <string_view>
  7. #include "TH1D.h"
  8. #include "TH2D.h"
  9. #include "THStack.h"
  10. #include "TGraph.h"
  11. #include "TTree.h"
  12. #include "TChain.h"
  13. #include "TFile.h"
  14. #include "TCanvas.h"
  15. #include "TROOT.h"
  16. #include "TStyle.h"
  17. #include "TColor.h"
  18. #include "TLorentzVector.h"
  19. #include "TRandom3.h"
  20. #include "TLegend.h"
  21. #include "RooDataHist.h"
  22. #include "RooRealVar.h"
  23. #include "RooPlot.h"
  24. #include "RooGaussian.h"
  25. #include "RooExponential.h"
  26. #include "RooRealConstant.h"
  27. #include "RooAddPdf.h"
  28. #include "RooFitResult.h"
  29. #include "RooProduct.h"
  30. #include "RooCrystalBall.h"
  31. #include "RooBreitWigner.h"
  32. #include "RooArgSet.h"
  33. #include "RooFFTConvPdf.h"
  34. #include "RooNovosibirsk.h"
  35. #include "constants.h"
  36. #include "basic_analysis.h"
  37. #include "hlt1_decision_analysis.h"
  38. #include "bdt_classification.h"
  39. // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% Bu To Hp Mu Mu %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
  40. int new_analysis_bu2hpmumu()
  41. {
  42. const char *analysis_name = "BuToHpMuMu";
  43. const char *data_tree_name = "BuToHpMuMu";
  44. const char *sim_tree_name = "BuToHpMuMu_noPID";
  45. const char *end_state_mass_literal = "m(#pi^{+}_{(#rightarrow K^{+})}#mu^{+}#mu^{-})";
  46. const bool retrain_bdt = false;
  47. TChain *data_chain = new TChain(TString::Format("%s/DecayTree", data_tree_name));
  48. // data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/spruce_magdown_2023_v0_tuple_90000000_v0r0p6288631.root");
  49. data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/spruce_magdown_2023_v0r1_tuple_90000000_2023_v0r0p6288631.root");
  50. FourVect *l14v_data = FourVect::Init(data_chain, "L1");
  51. FourVect *l24v_data = FourVect::Init(data_chain, "L2");
  52. FourVect *hp4v_data = FourVect::Init(data_chain, "Hp");
  53. TChain *sim_chain = new TChain(TString::Format("%s/DecayTree", sim_tree_name));
  54. sim_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/MC/rd_btoxll_simulation_fullstream_v0r0p6671378_BuToKpMuMu_12143001_magdown.root");
  55. FourVect *l14v_sim = FourVect::Init(sim_chain, "L1");
  56. FourVect *l24v_sim = FourVect::Init(sim_chain, "L2");
  57. FourVect *hp4v_sim = FourVect::Init(sim_chain, "Hp");
  58. Int_t B_BKGCAT, L1_TRUEID, L2_TRUEID, Hp_TRUEID;
  59. sim_chain->SetBranchAddress("L1_TRUEID", &L1_TRUEID);
  60. sim_chain->SetBranchAddress("L2_TRUEID", &L2_TRUEID);
  61. sim_chain->SetBranchAddress("Hp_TRUEID", &Hp_TRUEID);
  62. sim_chain->SetBranchAddress("B_BKGCAT", &B_BKGCAT);
  63. TH1D *h1_B_Mass_jpsi = new TH1D("h1_B_Mass_jpsi", TString::Format("B Mass, J/#psi Mode (%s)", analysis_name), B_MASS_HIST_BINS, B_MASS_VAR_MIN, B_MASS_VAR_MAX);
  64. TH1D *h1_B_Mass_psi2s = new TH1D("h1_B_Mass_psi2s", TString::Format("B Mass, #psi(2S) Mode (%s)", analysis_name), B_MASS_HIST_BINS, B_MASS_VAR_MIN, B_MASS_VAR_MAX);
  65. TH1D *h1_B_Mass_sim = new TH1D("h1_B_Mass_sim", TString::Format("B Mass, Simualted (%s_noPID)", analysis_name), B_MASS_HIST_BINS, B_MASS_VAR_MIN, B_MASS_VAR_MAX);
  66. TH2D *h2_Hlt1_flags_B_Mass = new TH2D("h2_Hlt1_flags_B_Mass", "Hlt1 Decision vs B Mass", 50, 5100, 5400, 13, 1., 14.);
  67. TH2D *h2_Hlt1_flags_excl_B_Mass = new TH2D("h2_Hlt1_flags_excl_B_Mass", "Excl Hlt1 Decision vs B Mass", 50, 5100, 5400, 13, 1., 14.);
  68. TH1D *h1_bdt_probs = new TH1D("h1_bdt_probs", "BDT Probabilities", 100, -2, 2);
  69. h1_B_Mass_jpsi->GetXaxis()->SetTitle(end_state_mass_literal);
  70. h1_B_Mass_psi2s->GetXaxis()->SetTitle(end_state_mass_literal);
  71. h2_Hlt1_flags_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  72. h2_Hlt1_flags_excl_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  73. ConnectHlt1Decisions(data_chain, h2_Hlt1_flags_B_Mass, h2_Hlt1_flags_excl_B_Mass);
  74. auto hlt1_decision_histos = CreateHlt1DecisionHistos(analysis_name);
  75. std::map<std::string, int> exclusive_hits{};
  76. std::vector<TV *> vars{
  77. TV::Float("B_PT", "B_PT"),
  78. TV::Float("B_BPVFDCHI2", "B_BPVFDCHI2"),
  79. TV::Float("B_BPVDIRA", "B_BPVDIRA"),
  80. TV::Float("Jpsi_BPVIPCHI2", "Jpsi_BPVIPCHI2"),
  81. // TV::Float("Jpsi_BPVDIRA", "Jpsi_BPVDIRA"),
  82. TV::Float("Jpsi_PT", "Jpsi_PT"),
  83. TV::Float("Hp_BPVIPCHI2", "Hp_BPVIPCHI2"),
  84. TV::Float("Hp_PT", "Hp_PT"),
  85. // TV::Double("Kplus_PID_K", "K_PID_K"),
  86. TV::Double("Hp_PROBNN_K", "Hp_PROBNN_K"),
  87. TV::Float("L1_BPVIPCHI2", "L1_BPVIPCHI2"),
  88. TV::Float("L2_BPVIPCHI2", "L2_BPVIPCHI2"),
  89. };
  90. TTree *sig_tree = new TTree("TreeS", "tree containing signal data");
  91. TTree *bkg_tree = new TTree("TreeB", "tree containing background data");
  92. ConnectVarsToData(vars, data_chain, sim_chain, sig_tree, bkg_tree);
  93. unsigned int data_entries = data_chain->GetEntries();
  94. unsigned int sim_entries = sim_chain->GetEntries();
  95. unsigned int bkg_events = 0;
  96. if (retrain_bdt)
  97. {
  98. std::cout << "# Starting with BDT retrain." << std::endl;
  99. for (unsigned int i = 0; i < data_entries; i++)
  100. {
  101. data_chain->GetEntry(i);
  102. Double_t reconstructed_B_Mass = (hp4v_data->LorentzVector(K_MASS) + l14v_data->LorentzVector() + l24v_data->LorentzVector()).M();
  103. if (std::all_of(vars.begin(), vars.end(), [](TV *v)
  104. { return v->IsDataFinite(); }))
  105. {
  106. if (reconstructed_B_Mass > 5500.)
  107. {
  108. bkg_tree->Fill();
  109. bkg_events++;
  110. }
  111. }
  112. PrintProgress(TString::Format("%s BKG Collection", analysis_name), data_entries, 10000, i);
  113. }
  114. }
  115. else {
  116. std::cout << "# Starting without BDT retrain." << std::endl;
  117. }
  118. for (unsigned int i = 0; i < sim_entries; i++)
  119. {
  120. sim_chain->GetEntry(i);
  121. Double_t reconstructed_B_Mass = (hp4v_sim->LorentzVector(K_MASS) + l14v_sim->LorentzVector() + l24v_sim->LorentzVector()).M();
  122. if (B_BKGCAT == 30 && TMath::Abs(L1_TRUEID) == PID_MUON && L2_TRUEID == -L1_TRUEID && TMath::Abs(Hp_TRUEID) == PID_KAON)
  123. {
  124. if (std::all_of(vars.begin(), vars.end(), [](TV *v)
  125. { return v->IsMCFinite(); }))
  126. {
  127. if (retrain_bdt)
  128. {
  129. sig_tree->Fill();
  130. }
  131. h1_B_Mass_sim->Fill(reconstructed_B_Mass);
  132. }
  133. }
  134. PrintProgress(TString::Format("%s SIG Collection", analysis_name), sim_entries, 10000, i);
  135. }
  136. if (retrain_bdt)
  137. {
  138. TrainBDT(vars, "BuToHpMuMu", sig_tree, bkg_tree);
  139. std::cout << "# Finished BDT retrain." << std::endl;
  140. }
  141. std::cout << "# Starting evaluation of data." << std::endl;
  142. Float_t *train_vars = new Float_t[vars.size()];
  143. auto reader = SetupReader(vars, train_vars, "BuToHpMuMu");
  144. const double mva_cut_value = -0.03;
  145. for (unsigned int i = 0; i < data_entries; i++)
  146. {
  147. data_chain->GetEntry(i);
  148. TLorentzVector dimuon = l14v_data->LorentzVector() + l24v_data->LorentzVector();
  149. Double_t reconstructed_B_Mass = (hp4v_data->LorentzVector(K_MASS) + dimuon).M();
  150. if (std::all_of(vars.begin(), vars.end(), [](TV *v)
  151. { return v->IsDataFinite(); }))
  152. {
  153. for (size_t j = 0; j < vars.size(); j++)
  154. {
  155. if (vars[j]->IsDouble())
  156. {
  157. train_vars[j] = vars[j]->GetDataDouble();
  158. }
  159. else if (vars[j]->IsFloat())
  160. {
  161. train_vars[j] = vars[j]->GetDataFloat();
  162. }
  163. }
  164. }
  165. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  166. {
  167. CheckHlt1Decisioins(h2_Hlt1_flags_B_Mass, h2_Hlt1_flags_excl_B_Mass, exclusive_hits, reconstructed_B_Mass);
  168. FillHlt1DecisionHistos(hlt1_decision_histos, reconstructed_B_Mass);
  169. }
  170. double mva_response = reader->EvaluateMVA("BDT");
  171. h1_bdt_probs->Fill(mva_response);
  172. if (mva_response > mva_cut_value)
  173. {
  174. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  175. {
  176. h1_B_Mass_jpsi->Fill(reconstructed_B_Mass);
  177. }
  178. else if (TMath::Abs(dimuon.M() - PSI2S_MASS) < 100.)
  179. {
  180. h1_B_Mass_psi2s->Fill(reconstructed_B_Mass);
  181. }
  182. }
  183. PrintProgress(TString::Format("%s BDT Evaluation", analysis_name), data_entries, 10000, i);
  184. }
  185. std::cout << "# Exclusive Hits" << std::endl;
  186. for (const auto &[line, hits] : exclusive_hits)
  187. {
  188. std::cout << line << ": " << hits << std::endl;
  189. }
  190. DrawInDefaultCanvas(h2_Hlt1_flags_B_Mass, analysis_name, 0.16, "COLZ");
  191. DrawInDefaultCanvas(h2_Hlt1_flags_excl_B_Mass, analysis_name, 0.16, "COLZ");
  192. DrawInDefaultCanvas(h1_B_Mass_jpsi, analysis_name, 0.1);
  193. DrawInDefaultCanvas(h1_B_Mass_psi2s, analysis_name, 0.1);
  194. auto roofit_hist_sim = CreateRooFitHistogramAndFitCB(h1_B_Mass_sim, false, false, ShapeParamters{});
  195. auto roofit_hist_jpsi_fitsum = CreateRooFitHistogramAndFitCB(h1_B_Mass_jpsi, true, true, roofit_hist_sim.shape_parameters);
  196. auto roofit_hist_psi2s_fitsum = CreateRooFitHistogramAndFitCB(h1_B_Mass_psi2s, true, true, roofit_hist_sim.shape_parameters);
  197. DrawInDefaultCanvas(roofit_hist_jpsi_fitsum, analysis_name);
  198. DrawInDefaultCanvas(roofit_hist_psi2s_fitsum, analysis_name);
  199. DrawInDefaultCanvas(roofit_hist_sim, analysis_name);
  200. DrawHlt1DecisionHistos(analysis_name, hlt1_decision_histos);
  201. DrawBDTProbs(h1_bdt_probs, mva_cut_value, analysis_name);
  202. return 0;
  203. }
  204. int analyze_bu2hpmumu_simulation()
  205. {
  206. const char *analysis_name = "BuToHpMuMu_noPID";
  207. const char *end_state_mass_literal = "m(#pi^{+}_{(#rightarrow K^{+})}#mu^{+}#mu^{-})";
  208. TChain *sim_chain = new TChain(TString::Format("%s/DecayTree", analysis_name));
  209. sim_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/MC/rd_btoxll_simulation_fullstream_v0r0p6671378_BuToKpMuMu_12143001_magdown.root");
  210. FourVect *l14v = FourVect::Init(sim_chain, "L1");
  211. FourVect *l24v = FourVect::Init(sim_chain, "L2");
  212. FourVect *hp4v = FourVect::Init(sim_chain, "Hp");
  213. Int_t B_BKGCAT, L1_TRUEID, L2_TRUEID, Hp_TRUEID;
  214. sim_chain->SetBranchAddress("L1_TRUEID", &L1_TRUEID);
  215. sim_chain->SetBranchAddress("L2_TRUEID", &L2_TRUEID);
  216. sim_chain->SetBranchAddress("Hp_TRUEID", &Hp_TRUEID);
  217. sim_chain->SetBranchAddress("B_BKGCAT", &B_BKGCAT);
  218. TH1D *h1_B_Mass_jpsi = new TH1D("h1_B_Mass_jpsi", TString::Format("B Mass, J/#psi Mode (%s)", analysis_name), B_MASS_HIST_BINS, B_MASS_VAR_MIN, B_MASS_VAR_MAX);
  219. TH1D *h1_B_Mass_psi2s = new TH1D("h1_B_Mass_psi2s", TString::Format("B Mass, #psi(2S) Mode (%s)", analysis_name), B_MASS_HIST_BINS, B_MASS_VAR_MIN, B_MASS_VAR_MAX);
  220. TH2D *h2_Hlt1_flags_B_Mass = new TH2D("h2_Hlt1_flags_B_Mass", "Hlt1 Decision vs B Mass", 50, 5100, 5400, 13, 1., 14.);
  221. TH2D *h2_Hlt1_flags_excl_B_Mass = new TH2D("h2_Hlt1_flags_excl_B_Mass", "Excl Hlt1 Decision vs B Mass", 50, 5100, 5400, 13, 1., 14.);
  222. h1_B_Mass_jpsi->GetXaxis()->SetTitle(end_state_mass_literal);
  223. h1_B_Mass_psi2s->GetXaxis()->SetTitle(end_state_mass_literal);
  224. h2_Hlt1_flags_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  225. h2_Hlt1_flags_excl_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  226. ConnectHlt1Decisions(sim_chain, h2_Hlt1_flags_B_Mass, h2_Hlt1_flags_excl_B_Mass);
  227. auto hlt1_decision_histos = CreateHlt1DecisionHistos(analysis_name);
  228. std::map<std::string, int> exclusive_hits{};
  229. unsigned int entries = sim_chain->GetEntries();
  230. for (unsigned int i = 0; i < entries; i++)
  231. {
  232. sim_chain->GetEntry(i);
  233. if (B_BKGCAT == 30 && TMath::Abs(L1_TRUEID) == PID_MUON && L2_TRUEID == -L1_TRUEID && TMath::Abs(Hp_TRUEID) == PID_KAON)
  234. {
  235. TLorentzVector dimuon = l14v->LorentzVector() + l24v->LorentzVector();
  236. Double_t reconstructed_B_Mass = (hp4v->LorentzVector(K_MASS) + dimuon).M();
  237. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  238. {
  239. CheckHlt1Decisioins(h2_Hlt1_flags_B_Mass, h2_Hlt1_flags_excl_B_Mass, exclusive_hits, reconstructed_B_Mass);
  240. FillHlt1DecisionHistos(hlt1_decision_histos, reconstructed_B_Mass);
  241. }
  242. if (CutHlt1DecisionsOr({"TwoTrackMuonMVA", "DiMuonHighMass"}))
  243. {
  244. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  245. {
  246. h1_B_Mass_jpsi->Fill(reconstructed_B_Mass);
  247. }
  248. else if (TMath::Abs(dimuon.M() - PSI2S_MASS) < 100.)
  249. {
  250. h1_B_Mass_psi2s->Fill(reconstructed_B_Mass);
  251. }
  252. }
  253. }
  254. PrintProgress(analysis_name, entries, 10000, i);
  255. }
  256. std::cout << "# Exclusive Hits" << std::endl;
  257. for (const auto &[line, hits] : exclusive_hits)
  258. {
  259. std::cout << line << ": " << hits << std::endl;
  260. }
  261. DrawInDefaultCanvas(h2_Hlt1_flags_B_Mass, analysis_name, 0.16, "COLZ");
  262. DrawInDefaultCanvas(h2_Hlt1_flags_excl_B_Mass, analysis_name, 0.16, "COLZ");
  263. DrawInDefaultCanvas(h1_B_Mass_jpsi, analysis_name, 0.1);
  264. DrawInDefaultCanvas(h1_B_Mass_psi2s, analysis_name, 0.1);
  265. auto roofit_hist_jpsi = CreateRooFitHistogram(h1_B_Mass_jpsi);
  266. auto roofit_hist_psi2s = CreateRooFitHistogram(h1_B_Mass_psi2s);
  267. DrawInDefaultCanvas(roofit_hist_jpsi, analysis_name, 0.1);
  268. DrawInDefaultCanvas(roofit_hist_psi2s, analysis_name, 0.1);
  269. DrawHlt1DecisionHistos(analysis_name, hlt1_decision_histos);
  270. return 0;
  271. }
  272. // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% B0 To Hp Hm Mu Mu %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
  273. int analyze_b02hphmmumu()
  274. {
  275. const char *analysis_name = "B0ToHpHmMuMu";
  276. const char *end_state_mass_literal = "m(#pi^{+}_{(#rightarrow K^{+})} #pi^{-} #mu^{+}#mu^{-})";
  277. TChain *data_chain = new TChain(TString::Format("%s/DecayTree", analysis_name));
  278. data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/spruce_magdown_2023_v0_tuple_90000000_v0r0p6288631.root");
  279. data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/spruce_magdown_2023_v0r1_tuple_90000000_2023_v0r0p6288631.root");
  280. Double_t Hp_PID_K, Hm_PID_K, Hp_PID_PI, Hm_PID_PI;
  281. data_chain->SetBranchAddress("Hp_PID_K", &Hp_PID_K);
  282. data_chain->SetBranchAddress("Hm_PID_K", &Hm_PID_K);
  283. data_chain->SetBranchAddress("Hp_PID_PI", &Hp_PID_PI);
  284. data_chain->SetBranchAddress("Hm_PID_PI", &Hm_PID_PI);
  285. FourVect *l14v = FourVect::Init(data_chain, "L1");
  286. FourVect *l24v = FourVect::Init(data_chain, "L2");
  287. FourVect *hp4v = FourVect::Init(data_chain, "Hp");
  288. FourVect *hm4v = FourVect::Init(data_chain, "Hm");
  289. TH1D *h1_B_Mass_jpsi = new TH1D("h1_B_Mass_jpsi", TString::Format("B Mass, J/#psi Mode (%s)", analysis_name), B_MASS_HIST_BINS, B_MASS_VAR_MIN, B_MASS_VAR_MAX);
  290. TH1D *h1_B_Mass_psi2s = new TH1D("h1_B_Mass_psi2s", TString::Format("B Mass, #psi(2S) Mode (%s)", analysis_name), B_MASS_HIST_BINS, B_MASS_VAR_MIN, B_MASS_VAR_MAX);
  291. TH2D *h2_Hlt1_flags_B_Mass = new TH2D("h2_Hlt1_flags_B_Mass", "Hlt1 Decision vs B Mass", 50, 5100, 5400, 13, 1., 14.);
  292. TH2D *h2_Hlt1_flags_excl_B_Mass = new TH2D("h2_Hlt1_flags_excl_B_Mass", "Excl Hlt1 Decision vs B Mass", 50, 5100, 5400, 13, 1., 14.);
  293. h1_B_Mass_jpsi->GetXaxis()->SetTitle(end_state_mass_literal);
  294. h1_B_Mass_psi2s->GetXaxis()->SetTitle(end_state_mass_literal);
  295. h2_Hlt1_flags_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  296. h2_Hlt1_flags_excl_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  297. ConnectHlt1Decisions(data_chain, h2_Hlt1_flags_B_Mass, h2_Hlt1_flags_excl_B_Mass);
  298. auto hlt1_decision_histos = CreateHlt1DecisionHistos(analysis_name);
  299. std::map<std::string, int> exclusive_hits{};
  300. unsigned int entries = data_chain->GetEntries();
  301. for (unsigned int i = 0; i < entries; i++)
  302. {
  303. data_chain->GetEntry(i);
  304. TLorentzVector dimuon = l14v->LorentzVector() + l24v->LorentzVector();
  305. TLorentzVector lv_hp{};
  306. TLorentzVector lv_hm{};
  307. if (Hp_PID_K > 0 && Hm_PID_K < 0)
  308. {
  309. continue;
  310. lv_hp = hp4v->LorentzVector(K_MASS);
  311. lv_hm = hm4v->LorentzVector();
  312. }
  313. else if (Hp_PID_K < 0 && Hm_PID_K > 0)
  314. {
  315. lv_hp = hp4v->LorentzVector();
  316. lv_hm = hm4v->LorentzVector(K_MASS);
  317. }
  318. else
  319. {
  320. continue;
  321. }
  322. Double_t reconstructed_Kstar_Mass = (lv_hp + lv_hm).M();
  323. Double_t reconstructed_B_Mass = (lv_hp + lv_hm + dimuon).M();
  324. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  325. {
  326. CheckHlt1Decisioins(h2_Hlt1_flags_B_Mass, h2_Hlt1_flags_excl_B_Mass, exclusive_hits, reconstructed_B_Mass);
  327. FillHlt1DecisionHistos(hlt1_decision_histos, reconstructed_B_Mass);
  328. }
  329. if (CutHlt1DecisionsOr({"TwoTrackMuonMVA", "DiMuonHighMass"}))
  330. {
  331. if (TMath::Abs(reconstructed_Kstar_Mass - KSTAR_MASS) < 50.)
  332. {
  333. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  334. {
  335. h1_B_Mass_jpsi->Fill(reconstructed_B_Mass);
  336. }
  337. else if (TMath::Abs(dimuon.M() - PSI2S_MASS) < 100.)
  338. {
  339. h1_B_Mass_psi2s->Fill(reconstructed_B_Mass);
  340. }
  341. }
  342. }
  343. PrintProgress(analysis_name, entries, 10000, i);
  344. }
  345. std::cout << "# Exclusive Hits" << std::endl;
  346. for (const auto &[line, hits] : exclusive_hits)
  347. {
  348. std::cout << line << ": " << hits << std::endl;
  349. }
  350. DrawInDefaultCanvas(h2_Hlt1_flags_B_Mass, analysis_name, 0.16, "COLZ");
  351. DrawInDefaultCanvas(h2_Hlt1_flags_excl_B_Mass, analysis_name, 0.16, "COLZ");
  352. DrawInDefaultCanvas(h1_B_Mass_jpsi, analysis_name, 0.1);
  353. DrawInDefaultCanvas(h1_B_Mass_psi2s, analysis_name, 0.1);
  354. auto roofit_hist_jpsi = CreateRooFitHistogram(h1_B_Mass_jpsi);
  355. auto roofit_hist_psi2s = CreateRooFitHistogram(h1_B_Mass_psi2s);
  356. DrawInDefaultCanvas(roofit_hist_jpsi, analysis_name, 0.1);
  357. DrawInDefaultCanvas(roofit_hist_psi2s, analysis_name, 0.1);
  358. DrawHlt1DecisionHistos(analysis_name, hlt1_decision_histos);
  359. return 0;
  360. }
  361. // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% Hlt2RD_BuToKpMuMu %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
  362. int analyze_Hlt2RD_BuToKpMuMu()
  363. {
  364. const char *analysis_name = "Hlt2RD_BuToKpMuMu";
  365. const char *end_state_mass_literal = "m(K^{+} #mu^{+}#mu^{-})";
  366. TChain *data_chain = new TChain(TString::Format("%s/DecayTree", analysis_name));
  367. data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/Collision23_Beam6800GeV-VeloClosed-MagDown-Excl-UT_RealData_SprucingPass23r1_94000000_RD.root");
  368. FourVect *l14v = FourVect::Init(data_chain, "muplus");
  369. FourVect *l24v = FourVect::Init(data_chain, "muminus");
  370. FourVect *hp4v = FourVect::Init(data_chain, "Kplus");
  371. TH1D *h1_B_Mass_jpsi = new TH1D("h1_B_Mass_jpsi", TString::Format("B Mass, J/#psi Mode (%s)", analysis_name), B_MASS_HIST_BINS, B_MASS_VAR_MIN, B_MASS_VAR_MAX);
  372. TH1D *h1_B_Mass_psi2s = new TH1D("h1_B_Mass_psi2s", TString::Format("B Mass, #psi(2S) Mode (%s)", analysis_name), B_MASS_HIST_BINS, B_MASS_VAR_MIN, B_MASS_VAR_MAX);
  373. TH2D *h2_Hlt1_flags_B_Mass = new TH2D("h2_Hlt1_flags_B_Mass", "Hlt1 Decision vs B Mass", 50, 5100, 5400, 13, 1., 14.);
  374. TH2D *h2_Hlt1_flags_excl_B_Mass = new TH2D("h2_Hlt1_flags_excl_B_Mass", "Excl Hlt1 Decision vs B Mass", 50, 5100, 5400, 13, 1., 14.);
  375. h1_B_Mass_jpsi->GetXaxis()->SetTitle(end_state_mass_literal);
  376. h1_B_Mass_psi2s->GetXaxis()->SetTitle(end_state_mass_literal);
  377. h2_Hlt1_flags_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  378. h2_Hlt1_flags_excl_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  379. ConnectHlt1Decisions(data_chain, h2_Hlt1_flags_B_Mass, h2_Hlt1_flags_excl_B_Mass);
  380. auto hlt1_decision_histos = CreateHlt1DecisionHistos(analysis_name);
  381. std::map<std::string, int> exclusive_hits{};
  382. unsigned int entries = data_chain->GetEntries();
  383. for (unsigned int i = 0; i < entries; i++)
  384. {
  385. data_chain->GetEntry(i);
  386. TLorentzVector dimuon = l14v->LorentzVector() + l24v->LorentzVector();
  387. TLorentzVector lv_hp = hp4v->LorentzVector();
  388. Double_t reconstructed_B_Mass = (lv_hp + dimuon).M();
  389. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  390. {
  391. CheckHlt1Decisioins(h2_Hlt1_flags_B_Mass, h2_Hlt1_flags_excl_B_Mass, exclusive_hits, reconstructed_B_Mass);
  392. FillHlt1DecisionHistos(hlt1_decision_histos, reconstructed_B_Mass);
  393. }
  394. if (CutHlt1DecisionsOr({"TwoTrackMuonMVA", "DiMuonHighMass"}))
  395. {
  396. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  397. {
  398. h1_B_Mass_jpsi->Fill(reconstructed_B_Mass);
  399. }
  400. else if (TMath::Abs(dimuon.M() - PSI2S_MASS) < 100.)
  401. {
  402. h1_B_Mass_psi2s->Fill(reconstructed_B_Mass);
  403. }
  404. }
  405. PrintProgress("Hlt2RD_BuToKpMuMu", entries, 10000, i);
  406. }
  407. std::cout << "# Exclusive Hits" << std::endl;
  408. for (const auto &[line, hits] : exclusive_hits)
  409. {
  410. std::cout << line << ": " << hits << std::endl;
  411. }
  412. DrawInDefaultCanvas(h2_Hlt1_flags_B_Mass, analysis_name, 0.16, "COLZ");
  413. DrawInDefaultCanvas(h2_Hlt1_flags_excl_B_Mass, analysis_name, 0.16, "COLZ");
  414. DrawInDefaultCanvas(h1_B_Mass_jpsi, analysis_name, 0.1);
  415. DrawInDefaultCanvas(h1_B_Mass_psi2s, analysis_name, 0.1);
  416. auto roofit_hist_jpsi = CreateRooFitHistogram(h1_B_Mass_jpsi);
  417. auto roofit_hist_psi2s = CreateRooFitHistogram(h1_B_Mass_psi2s);
  418. DrawInDefaultCanvas(roofit_hist_jpsi, analysis_name, 0.1);
  419. DrawInDefaultCanvas(roofit_hist_psi2s, analysis_name, 0.1);
  420. DrawHlt1DecisionHistos(analysis_name, hlt1_decision_histos);
  421. return 0;
  422. }
  423. // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% Hlt2RD_B0ToKpPimMuMu %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
  424. int analyze_Hlt2RD_B0ToKpPimMuMu()
  425. {
  426. const char *analysis_name = "Hlt2RD_B0ToKpPimMuMu";
  427. const char *end_state_mass_literal = "m(K^{+} #pi^{-} #mu^{+}#mu^{-})";
  428. TChain *data_chain = new TChain(TString::Format("%s/DecayTree", analysis_name));
  429. data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/Collision23_Beam6800GeV-VeloClosed-MagDown-Excl-UT_RealData_SprucingPass23r1_94000000_RD.root");
  430. FourVect *l14v = FourVect::Init(data_chain, "muplus");
  431. FourVect *l24v = FourVect::Init(data_chain, "muminus");
  432. FourVect *hp4v = FourVect::Init(data_chain, "Kplus");
  433. FourVect *hm4v = FourVect::Init(data_chain, "piminus");
  434. TH1D *h1_B_Mass_jpsi = new TH1D("h1_B_Mass_jpsi", TString::Format("B Mass, J/#psi Mode (%s)", analysis_name), B_MASS_HIST_BINS, B_MASS_VAR_MIN, B_MASS_VAR_MAX);
  435. TH1D *h1_B_Mass_psi2s = new TH1D("h1_B_Mass_psi2s", TString::Format("B Mass, #psi(2S) Mode (%s)", analysis_name), B_MASS_HIST_BINS, B_MASS_VAR_MIN, B_MASS_VAR_MAX);
  436. TH2D *h2_Hlt1_flags_B_Mass = new TH2D("h2_Hlt1_flags_B_Mass", "Hlt1 Decision vs B Mass", 50, 5100, 5400, 13, 1., 14.);
  437. TH2D *h2_Hlt1_flags_excl_B_Mass = new TH2D("h2_Hlt1_flags_excl_B_Mass", "Excl Hlt1 Decision vs B Mass", 50, 5100, 5400, 13, 1., 14.);
  438. h1_B_Mass_jpsi->GetXaxis()->SetTitle(end_state_mass_literal);
  439. h1_B_Mass_psi2s->GetXaxis()->SetTitle(end_state_mass_literal);
  440. h2_Hlt1_flags_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  441. h2_Hlt1_flags_excl_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  442. ConnectHlt1Decisions(data_chain, h2_Hlt1_flags_B_Mass, h2_Hlt1_flags_excl_B_Mass);
  443. auto hlt1_decision_histos = CreateHlt1DecisionHistos(analysis_name);
  444. std::map<std::string, int> exclusive_hits{};
  445. unsigned int entries = data_chain->GetEntries();
  446. for (unsigned int i = 0; i < entries; i++)
  447. {
  448. data_chain->GetEntry(i);
  449. TLorentzVector dimuon = l14v->LorentzVector() + l24v->LorentzVector();
  450. TLorentzVector lv_hp = hp4v->LorentzVector();
  451. TLorentzVector lv_hm = hm4v->LorentzVector();
  452. Double_t reconstructed_Kstar_Mass = (lv_hp + lv_hm).M();
  453. Double_t reconstructed_B_Mass = (lv_hp + lv_hm + dimuon).M();
  454. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  455. {
  456. CheckHlt1Decisioins(h2_Hlt1_flags_B_Mass, h2_Hlt1_flags_excl_B_Mass, exclusive_hits, reconstructed_B_Mass);
  457. FillHlt1DecisionHistos(hlt1_decision_histos, reconstructed_B_Mass);
  458. }
  459. if (CutHlt1DecisionsOr({"TwoTrackMuonMVA", "DiMuonHighMass"}))
  460. {
  461. if (TMath::Abs(reconstructed_Kstar_Mass - KSTAR_MASS) < 50.)
  462. {
  463. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  464. {
  465. h1_B_Mass_jpsi->Fill(reconstructed_B_Mass);
  466. }
  467. else if (TMath::Abs(dimuon.M() - PSI2S_MASS) < 100.)
  468. {
  469. h1_B_Mass_psi2s->Fill(reconstructed_B_Mass);
  470. }
  471. }
  472. }
  473. PrintProgress(analysis_name, entries, 10000, i);
  474. }
  475. std::cout << "# Exclusive Hits" << std::endl;
  476. for (const auto &[line, hits] : exclusive_hits)
  477. {
  478. std::cout << line << ": " << hits << std::endl;
  479. }
  480. DrawInDefaultCanvas(h2_Hlt1_flags_B_Mass, analysis_name, 0.16, "COLZ");
  481. DrawInDefaultCanvas(h2_Hlt1_flags_excl_B_Mass, analysis_name, 0.16, "COLZ");
  482. DrawInDefaultCanvas(h1_B_Mass_jpsi, analysis_name, 0.1);
  483. DrawInDefaultCanvas(h1_B_Mass_psi2s, analysis_name, 0.1);
  484. auto roofit_hist_jpsi = CreateRooFitHistogram(h1_B_Mass_jpsi);
  485. auto roofit_hist_psi2s = CreateRooFitHistogram(h1_B_Mass_psi2s);
  486. DrawInDefaultCanvas(roofit_hist_jpsi, analysis_name, 0.1);
  487. DrawInDefaultCanvas(roofit_hist_psi2s, analysis_name, 0.1);
  488. DrawHlt1DecisionHistos(analysis_name, hlt1_decision_histos);
  489. return 0;
  490. }
  491. // int new_analysis()
  492. // {
  493. // return analyze_bu2hpmumu_data();
  494. // }