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.

731 lines
23 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. #include <algorithm>
  33. #include <filesystem>
  34. const Double_t B_MASS_VAR_MIN = 5100.;
  35. const Double_t B_MASS_VAR_MAX = 6000.;
  36. const Int_t B_MASS_HIST_BINS = 70;
  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. class FourVect
  42. {
  43. private:
  44. Float_t px, py, pz, energy;
  45. Float_t *PtrPX()
  46. {
  47. return &px;
  48. }
  49. Float_t *PtrPY()
  50. {
  51. return &py;
  52. }
  53. Float_t *PtrPZ()
  54. {
  55. return &pz;
  56. }
  57. Float_t *PtrEnergy()
  58. {
  59. return &energy;
  60. }
  61. public:
  62. static FourVect *Init(TTree *tree, const char *particle)
  63. {
  64. FourVect *v = new FourVect{};
  65. tree->SetBranchAddress(TString::Format("%s_PX", particle).Data(), v->PtrPX());
  66. tree->SetBranchAddress(TString::Format("%s_PY", particle).Data(), v->PtrPY());
  67. tree->SetBranchAddress(TString::Format("%s_PZ", particle).Data(), v->PtrPZ());
  68. tree->SetBranchAddress(TString::Format("%s_ENERGY", particle).Data(), v->PtrEnergy());
  69. return v;
  70. }
  71. TLorentzVector LorentzVector()
  72. {
  73. return TLorentzVector(px, py, pz, energy);
  74. }
  75. TLorentzVector LorentzVector(double sub_mass)
  76. {
  77. TVector3 momentum(px, py, pz);
  78. float energy = TMath::Sqrt(TMath::Sq(sub_mass) + momentum.Mag2());
  79. return TLorentzVector(momentum, energy);
  80. }
  81. std::string ToString()
  82. {
  83. return TString::Format("(%f, %f, %f, %f)", px, py, pz, energy).Data();
  84. }
  85. };
  86. struct Hlt1Decision
  87. {
  88. std::string name;
  89. int index;
  90. Bool_t value;
  91. std::string GetName() const
  92. {
  93. return TString::Format("Hlt1%sDecision", name.c_str()).Data();
  94. }
  95. Bool_t *GetValuePointer()
  96. {
  97. return &value;
  98. }
  99. };
  100. std::vector<Hlt1Decision> Hlt1Decisions{
  101. Hlt1Decision{"DiMuonHighMass", 1},
  102. Hlt1Decision{"DiMuonLowMass", 2},
  103. Hlt1Decision{"DiMuonNoIP", 3},
  104. Hlt1Decision{"DiMuonSoft", 4},
  105. Hlt1Decision{"DisplacedLeptons", 5},
  106. Hlt1Decision{"LowPtDiMuon", 6},
  107. Hlt1Decision{"LowPtMuon", 7},
  108. Hlt1Decision{"OneMuonTrackLine", 8},
  109. Hlt1Decision{"SingleHighEt", 9},
  110. Hlt1Decision{"SingleHighPtMuon", 10},
  111. Hlt1Decision{"TrackMVA", 11},
  112. Hlt1Decision{"TrackMuonMVA", 12},
  113. Hlt1Decision{"TwoTrackMVA", 13},
  114. };
  115. bool CutHlt1DecisionsAnd(const std::set<std::string> &decision_list)
  116. {
  117. bool okay = true;
  118. for (const auto &var : Hlt1Decisions)
  119. {
  120. if (decision_list.find(var.name) != decision_list.end())
  121. {
  122. okay = okay && var.value;
  123. }
  124. }
  125. return okay;
  126. }
  127. bool CutHlt1DecisionsOr(const std::set<std::string> &decision_list)
  128. {
  129. bool okay = false;
  130. for (const auto &var : Hlt1Decisions)
  131. {
  132. if (decision_list.find(var.name) != decision_list.end())
  133. {
  134. okay = okay || var.value;
  135. }
  136. }
  137. return okay;
  138. }
  139. void ConnectHlt1Decisions(TChain *chain)
  140. {
  141. for (auto &var : Hlt1Decisions)
  142. {
  143. if (chain->FindBranch(var.GetName().c_str()))
  144. chain->SetBranchAddress(var.GetName().c_str(), var.GetValuePointer());
  145. }
  146. }
  147. void DrawInDefaultCanvas(TH1 *histogram, const char *folder, double margin_left = 0, Option_t *option = "")
  148. {
  149. std::filesystem::create_directory(TString::Format("output_files/analysis/%s", folder).Data());
  150. TString name = TString::Format("%s_canvas", histogram->GetName());
  151. TCanvas *c = new TCanvas(name, histogram->GetName(), 0, 0, 800, 600);
  152. c->SetLeftMargin(margin_left);
  153. histogram->SetStats(0);
  154. histogram->Draw(option);
  155. c->Draw();
  156. c->SaveAs(TString::Format("output_files/analysis/%s/%s.pdf", folder, name.Data()).Data());
  157. }
  158. void DrawInDefaultCanvas(RooPlot *rooPlot, const char *folder, double margin_left = 0, Option_t *option = "")
  159. {
  160. std::filesystem::create_directory(TString::Format("output_files/analysis/%s", folder).Data());
  161. TString name = TString::Format("%s_canvas", rooPlot->GetName());
  162. TCanvas *c = new TCanvas(name, rooPlot->GetName(), 0, 0, 800, 600);
  163. c->SetLeftMargin(margin_left);
  164. rooPlot->Draw(option);
  165. c->Draw();
  166. c->SaveAs(TString::Format("output_files/analysis/%s/%s.pdf", folder, name.Data()).Data());
  167. }
  168. void PrintProgress(const char *title, unsigned int total, unsigned int every, unsigned int current)
  169. {
  170. if ((current + 1) % every == 0 || current + 1 == total)
  171. {
  172. std::cout << "["
  173. << title
  174. << "] Processed event: " << current + 1 << " (" << TString::Format("%.2f", ((double)(current + 1) / (double)total) * 100.) << "%)" << std::endl;
  175. }
  176. }
  177. RooPlot* CreateRooFitHistogram(TH1D* hist) {
  178. RooRealVar roo_var_mass("var_mass", "B Mass Variable", B_MASS_VAR_MIN, B_MASS_VAR_MAX);
  179. roo_var_mass.setRange("fitting_range", B_MASS_VAR_MIN, B_MASS_VAR_MAX);
  180. RooDataHist roohist_B_M("roohist_B_M", "B Mass Histogram", roo_var_mass, RooFit::Import(*hist));
  181. RooPlot *roo_frame_mass = roo_var_mass.frame(RooFit::Title(hist->GetTitle()));
  182. roohist_B_M.plotOn(roo_frame_mass, RooFit::Binning(B_MASS_HIST_BINS), RooFit::Name("B Mass Distribution"));
  183. roo_frame_mass->GetXaxis()->SetTitle(hist->GetXaxis()->GetTitle());
  184. return roo_frame_mass;
  185. }
  186. // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% Bu To Hp Mu Mu %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
  187. int analyze_bu2hpmumu()
  188. {
  189. const char *analysis_name = "BuToHpMuMu";
  190. const char *end_state_mass_literal = "m(#pi^{+}_{(#rightarrow K^{+})}#mu^{+}#mu^{-})";
  191. TChain *data_chain = new TChain(TString::Format("%s/DecayTree", analysis_name));
  192. data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/spruce_magdown_2023_v0_tuple_90000000_v0r0p6288631.root");
  193. data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/spruce_magdown_2023_v0r1_tuple_90000000_2023_v0r0p6288631.root");
  194. FourVect *l14v = FourVect::Init(data_chain, "L1");
  195. FourVect *l24v = FourVect::Init(data_chain, "L2");
  196. FourVect *hp4v = FourVect::Init(data_chain, "Hp");
  197. 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);
  198. 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);
  199. TH2D *h2_Hlt1_flags_B_Mass = new TH2D("h2_Hlt1_flags_B_Mass", "Hlt1 Decision vs B Mass", 50, 5100, 5400, 13, 1., 14.);
  200. 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.);
  201. h1_B_Mass_jpsi->GetXaxis()->SetTitle(end_state_mass_literal);
  202. h1_B_Mass_psi2s->GetXaxis()->SetTitle(end_state_mass_literal);
  203. h2_Hlt1_flags_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  204. h2_Hlt1_flags_excl_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  205. ConnectHlt1Decisions(data_chain);
  206. for (const auto &var : Hlt1Decisions)
  207. {
  208. if (data_chain->FindBranch(var.GetName().c_str()))
  209. {
  210. h2_Hlt1_flags_B_Mass->Fill(0., var.name.c_str(), 0.);
  211. h2_Hlt1_flags_excl_B_Mass->Fill(0., var.name.c_str(), 0.);
  212. }
  213. }
  214. std::map<std::string, int> exclusive_hits{};
  215. unsigned int entries = data_chain->GetEntries();
  216. for (unsigned int i = 0; i < entries; i++)
  217. {
  218. data_chain->GetEntry(i);
  219. TLorentzVector dimuon = l14v->LorentzVector() + l24v->LorentzVector();
  220. Double_t reconstructed_B_Mass = (hp4v->LorentzVector(K_MASS) + dimuon).M();
  221. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  222. {
  223. for (const auto &var : Hlt1Decisions)
  224. {
  225. if (var.value)
  226. {
  227. h2_Hlt1_flags_B_Mass->Fill(reconstructed_B_Mass, var.name.c_str(), 1);
  228. }
  229. }
  230. bool exclusive = true;
  231. std::string line{};
  232. for (const auto &var : Hlt1Decisions)
  233. {
  234. if (var.value)
  235. {
  236. if (!line.empty())
  237. {
  238. exclusive = false;
  239. }
  240. line = var.name;
  241. }
  242. }
  243. if (!line.empty() && exclusive)
  244. {
  245. int &hits = exclusive_hits[line];
  246. if (hits)
  247. {
  248. hits += 1;
  249. }
  250. else
  251. {
  252. hits = 1;
  253. }
  254. h2_Hlt1_flags_excl_B_Mass->Fill(reconstructed_B_Mass, line.c_str(), 1);
  255. }
  256. }
  257. if (CutHlt1DecisionsOr({"TwoTrackMuonMVA", "TrackMuonMVA"}))
  258. {
  259. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  260. {
  261. h1_B_Mass_jpsi->Fill(reconstructed_B_Mass);
  262. }
  263. else if (TMath::Abs(dimuon.M() - PSI2S_MASS) < 100.)
  264. {
  265. h1_B_Mass_psi2s->Fill(reconstructed_B_Mass);
  266. }
  267. }
  268. PrintProgress(analysis_name, entries, 10000, i);
  269. }
  270. std::cout << "# Exclusive Hits" << std::endl;
  271. for (const auto &[line, hits] : exclusive_hits)
  272. {
  273. std::cout << line << ": " << hits << std::endl;
  274. }
  275. DrawInDefaultCanvas(h2_Hlt1_flags_B_Mass, analysis_name, 0.16, "COLZ");
  276. DrawInDefaultCanvas(h2_Hlt1_flags_excl_B_Mass, analysis_name, 0.16, "COLZ");
  277. DrawInDefaultCanvas(h1_B_Mass_jpsi, analysis_name, 0.1);
  278. DrawInDefaultCanvas(h1_B_Mass_psi2s, analysis_name, 0.1);
  279. auto roofit_hist_jpsi = CreateRooFitHistogram(h1_B_Mass_jpsi);
  280. auto roofit_hist_psi2s = CreateRooFitHistogram(h1_B_Mass_psi2s);
  281. DrawInDefaultCanvas(roofit_hist_jpsi, analysis_name, 0.1);
  282. DrawInDefaultCanvas(roofit_hist_psi2s, analysis_name, 0.1);
  283. return 0;
  284. }
  285. // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% B0 To Hp Hm Mu Mu %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
  286. int analyze_b02hphmmumu()
  287. {
  288. const char *analysis_name = "B0ToHpHmMuMu";
  289. const char *end_state_mass_literal = "m(#pi^{+}_{(#rightarrow K^{+})} #pi^{-} #mu^{+}#mu^{-})";
  290. TChain *data_chain = new TChain(TString::Format("%s/DecayTree", analysis_name));
  291. data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/spruce_magdown_2023_v0_tuple_90000000_v0r0p6288631.root");
  292. data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/spruce_magdown_2023_v0r1_tuple_90000000_2023_v0r0p6288631.root");
  293. Double_t Hp_PID_K, Hm_PID_K, Hp_PID_PI, Hm_PID_PI;
  294. data_chain->SetBranchAddress("Hp_PID_K", &Hp_PID_K);
  295. data_chain->SetBranchAddress("Hm_PID_K", &Hm_PID_K);
  296. data_chain->SetBranchAddress("Hp_PID_PI", &Hp_PID_PI);
  297. data_chain->SetBranchAddress("Hm_PID_PI", &Hm_PID_PI);
  298. FourVect *l14v = FourVect::Init(data_chain, "L1");
  299. FourVect *l24v = FourVect::Init(data_chain, "L2");
  300. FourVect *hp4v = FourVect::Init(data_chain, "Hp");
  301. FourVect *hm4v = FourVect::Init(data_chain, "Hm");
  302. 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);
  303. 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);
  304. TH2D *h2_Hlt1_flags_B_Mass = new TH2D("h2_Hlt1_flags_B_Mass", "Hlt1 Decision vs B Mass", 50, 5100, 5400, 13, 1., 14.);
  305. 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.);
  306. h1_B_Mass_jpsi->GetXaxis()->SetTitle(end_state_mass_literal);
  307. h1_B_Mass_psi2s->GetXaxis()->SetTitle(end_state_mass_literal);
  308. h2_Hlt1_flags_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  309. h2_Hlt1_flags_excl_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  310. ConnectHlt1Decisions(data_chain);
  311. for (const auto &var : Hlt1Decisions)
  312. {
  313. h2_Hlt1_flags_B_Mass->Fill(0., var.name.c_str(), 0.);
  314. h2_Hlt1_flags_excl_B_Mass->Fill(0., var.name.c_str(), 0.);
  315. }
  316. std::map<std::string, int> exclusive_hits{};
  317. unsigned int entries = data_chain->GetEntries();
  318. for (unsigned int i = 0; i < entries; i++)
  319. {
  320. data_chain->GetEntry(i);
  321. TLorentzVector dimuon = l14v->LorentzVector() + l24v->LorentzVector();
  322. TLorentzVector lv_hp{};
  323. TLorentzVector lv_hm{};
  324. if (Hp_PID_K > 0 && Hm_PID_K < 0)
  325. {
  326. lv_hp = hp4v->LorentzVector(K_MASS);
  327. lv_hm = hm4v->LorentzVector();
  328. }
  329. else if (Hp_PID_K < 0 && Hm_PID_K > 0)
  330. {
  331. continue;
  332. lv_hp = hp4v->LorentzVector();
  333. lv_hm = hm4v->LorentzVector(K_MASS);
  334. }
  335. else
  336. {
  337. continue;
  338. }
  339. Double_t reconstructed_Kstar_Mass = (lv_hp + lv_hm).M();
  340. Double_t reconstructed_B_Mass = (lv_hp + lv_hm + dimuon).M();
  341. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  342. {
  343. for (const auto &var : Hlt1Decisions)
  344. {
  345. if (var.value)
  346. {
  347. h2_Hlt1_flags_B_Mass->Fill(reconstructed_B_Mass, var.name.c_str(), 1);
  348. }
  349. }
  350. bool exclusive = true;
  351. std::string line{};
  352. for (const auto &var : Hlt1Decisions)
  353. {
  354. if (var.value)
  355. {
  356. if (!line.empty())
  357. {
  358. exclusive = false;
  359. }
  360. line = var.name;
  361. }
  362. }
  363. if (!line.empty() && exclusive)
  364. {
  365. int &hits = exclusive_hits[line];
  366. if (hits)
  367. {
  368. hits += 1;
  369. }
  370. else
  371. {
  372. hits = 1;
  373. }
  374. h2_Hlt1_flags_excl_B_Mass->Fill(reconstructed_B_Mass, line.c_str(), 1);
  375. }
  376. }
  377. if (CutHlt1DecisionsOr({"TwoTrackMuonMVA", "TrackMuonMVA"}))
  378. {
  379. if (TMath::Abs(reconstructed_Kstar_Mass - KSTAR_MASS) < 50.)
  380. {
  381. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  382. {
  383. h1_B_Mass_jpsi->Fill(reconstructed_B_Mass);
  384. }
  385. else if (TMath::Abs(dimuon.M() - PSI2S_MASS) < 100.)
  386. {
  387. h1_B_Mass_psi2s->Fill(reconstructed_B_Mass);
  388. }
  389. }
  390. }
  391. PrintProgress(analysis_name, entries, 10000, i);
  392. }
  393. std::cout << "# Exclusive Hits" << std::endl;
  394. for (const auto &[line, hits] : exclusive_hits)
  395. {
  396. std::cout << line << ": " << hits << std::endl;
  397. }
  398. DrawInDefaultCanvas(h2_Hlt1_flags_B_Mass, analysis_name, 0.16, "COLZ");
  399. DrawInDefaultCanvas(h2_Hlt1_flags_excl_B_Mass, analysis_name, 0.16, "COLZ");
  400. DrawInDefaultCanvas(h1_B_Mass_jpsi, analysis_name, 0.1);
  401. DrawInDefaultCanvas(h1_B_Mass_psi2s, analysis_name, 0.1);
  402. auto roofit_hist_jpsi = CreateRooFitHistogram(h1_B_Mass_jpsi);
  403. auto roofit_hist_psi2s = CreateRooFitHistogram(h1_B_Mass_psi2s);
  404. DrawInDefaultCanvas(roofit_hist_jpsi, analysis_name, 0.1);
  405. DrawInDefaultCanvas(roofit_hist_psi2s, analysis_name, 0.1);
  406. return 0;
  407. }
  408. // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% Hlt2RD_BuToKpMuMu %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
  409. int analyze_Hlt2RD_BuToKpMuMu()
  410. {
  411. const char *analysis_name = "Hlt2RD_BuToKpMuMu";
  412. const char *end_state_mass_literal = "m(K^{+} #mu^{+}#mu^{-})";
  413. TChain *data_chain = new TChain(TString::Format("%s/DecayTree", analysis_name));
  414. data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/Collision23_Beam6800GeV-VeloClosed-MagDown-Excl-UT_RealData_SprucingPass23r1_94000000_RD.root");
  415. FourVect *l14v = FourVect::Init(data_chain, "muplus");
  416. FourVect *l24v = FourVect::Init(data_chain, "muminus");
  417. FourVect *hp4v = FourVect::Init(data_chain, "Kplus");
  418. 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);
  419. 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);
  420. TH2D *h2_Hlt1_flags_B_Mass = new TH2D("h2_Hlt1_flags_B_Mass", "Hlt1 Decision vs B Mass", 50, 5100, 5400, 13, 1., 14.);
  421. 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.);
  422. h1_B_Mass_jpsi->GetXaxis()->SetTitle(end_state_mass_literal);
  423. h1_B_Mass_psi2s->GetXaxis()->SetTitle(end_state_mass_literal);
  424. h2_Hlt1_flags_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  425. h2_Hlt1_flags_excl_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  426. ConnectHlt1Decisions(data_chain);
  427. for (const auto &var : Hlt1Decisions)
  428. {
  429. h2_Hlt1_flags_B_Mass->Fill(0., var.name.c_str(), 0.);
  430. h2_Hlt1_flags_excl_B_Mass->Fill(0., var.name.c_str(), 0.);
  431. }
  432. std::map<std::string, int> exclusive_hits{};
  433. unsigned int entries = data_chain->GetEntries();
  434. for (unsigned int i = 0; i < entries; i++)
  435. {
  436. data_chain->GetEntry(i);
  437. TLorentzVector dimuon = l14v->LorentzVector() + l24v->LorentzVector();
  438. TLorentzVector lv_hp = hp4v->LorentzVector();
  439. Double_t reconstructed_B_Mass = (lv_hp + dimuon).M();
  440. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  441. {
  442. for (const auto &var : Hlt1Decisions)
  443. {
  444. if (var.value)
  445. {
  446. h2_Hlt1_flags_B_Mass->Fill(reconstructed_B_Mass, var.name.c_str(), 1);
  447. }
  448. }
  449. bool exclusive = true;
  450. std::string line{};
  451. for (const auto &var : Hlt1Decisions)
  452. {
  453. if (var.value)
  454. {
  455. if (!line.empty())
  456. {
  457. exclusive = false;
  458. }
  459. line = var.name;
  460. }
  461. }
  462. if (!line.empty() && exclusive)
  463. {
  464. int &hits = exclusive_hits[line];
  465. if (hits)
  466. {
  467. hits += 1;
  468. }
  469. else
  470. {
  471. hits = 1;
  472. }
  473. h2_Hlt1_flags_excl_B_Mass->Fill(reconstructed_B_Mass, line.c_str(), 1);
  474. }
  475. }
  476. if (CutHlt1DecisionsOr({"TwoTrackMuonMVA", "TrackMuonMVA"}))
  477. {
  478. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  479. {
  480. h1_B_Mass_jpsi->Fill(reconstructed_B_Mass);
  481. }
  482. else if (TMath::Abs(dimuon.M() - PSI2S_MASS) < 100.)
  483. {
  484. h1_B_Mass_psi2s->Fill(reconstructed_B_Mass);
  485. }
  486. }
  487. PrintProgress("Hlt2RD_BuToKpMuMu", entries, 10000, i);
  488. }
  489. std::cout << "# Exclusive Hits" << std::endl;
  490. for (const auto &[line, hits] : exclusive_hits)
  491. {
  492. std::cout << line << ": " << hits << std::endl;
  493. }
  494. DrawInDefaultCanvas(h2_Hlt1_flags_B_Mass, analysis_name, 0.16, "COLZ");
  495. DrawInDefaultCanvas(h2_Hlt1_flags_excl_B_Mass, analysis_name, 0.16, "COLZ");
  496. DrawInDefaultCanvas(h1_B_Mass_jpsi, analysis_name, 0.1);
  497. DrawInDefaultCanvas(h1_B_Mass_psi2s, analysis_name, 0.1);
  498. auto roofit_hist_jpsi = CreateRooFitHistogram(h1_B_Mass_jpsi);
  499. auto roofit_hist_psi2s = CreateRooFitHistogram(h1_B_Mass_psi2s);
  500. DrawInDefaultCanvas(roofit_hist_jpsi, analysis_name, 0.1);
  501. DrawInDefaultCanvas(roofit_hist_psi2s, analysis_name, 0.1);
  502. return 0;
  503. }
  504. // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% Hlt2RD_B0ToKpPimMuMu %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
  505. int analyze_Hlt2RD_B0ToKpPimMuMu()
  506. {
  507. const char *analysis_name = "Hlt2RD_B0ToKpPimMuMu";
  508. const char *end_state_mass_literal = "m(K^{+} #pi^{-} #mu^{+}#mu^{-})";
  509. TChain *data_chain = new TChain(TString::Format("%s/DecayTree", analysis_name));
  510. data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/Collision23_Beam6800GeV-VeloClosed-MagDown-Excl-UT_RealData_SprucingPass23r1_94000000_RD.root");
  511. FourVect *l14v = FourVect::Init(data_chain, "muplus");
  512. FourVect *l24v = FourVect::Init(data_chain, "muminus");
  513. FourVect *hp4v = FourVect::Init(data_chain, "Kplus");
  514. FourVect *hm4v = FourVect::Init(data_chain, "piminus");
  515. 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);
  516. 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);
  517. TH2D *h2_Hlt1_flags_B_Mass = new TH2D("h2_Hlt1_flags_B_Mass", "Hlt1 Decision vs B Mass", 50, 5100, 5400, 13, 1., 14.);
  518. 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.);
  519. h1_B_Mass_jpsi->GetXaxis()->SetTitle(end_state_mass_literal);
  520. h1_B_Mass_psi2s->GetXaxis()->SetTitle(end_state_mass_literal);
  521. h2_Hlt1_flags_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  522. h2_Hlt1_flags_excl_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  523. ConnectHlt1Decisions(data_chain);
  524. for (const auto &var : Hlt1Decisions)
  525. {
  526. h2_Hlt1_flags_B_Mass->Fill(0., var.name.c_str(), 0.);
  527. h2_Hlt1_flags_excl_B_Mass->Fill(0., var.name.c_str(), 0.);
  528. }
  529. std::map<std::string, int> exclusive_hits{};
  530. unsigned int entries = data_chain->GetEntries();
  531. for (unsigned int i = 0; i < entries; i++)
  532. {
  533. data_chain->GetEntry(i);
  534. TLorentzVector dimuon = l14v->LorentzVector() + l24v->LorentzVector();
  535. TLorentzVector lv_hp = hp4v->LorentzVector();
  536. TLorentzVector lv_hm = hm4v->LorentzVector();
  537. Double_t reconstructed_Kstar_Mass = (lv_hp + lv_hm).M();
  538. Double_t reconstructed_B_Mass = (lv_hp + lv_hm + dimuon).M();
  539. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  540. {
  541. for (const auto &var : Hlt1Decisions)
  542. {
  543. if (var.value)
  544. {
  545. h2_Hlt1_flags_B_Mass->Fill(reconstructed_B_Mass, var.name.c_str(), 1);
  546. }
  547. }
  548. bool exclusive = true;
  549. std::string line{};
  550. for (const auto &var : Hlt1Decisions)
  551. {
  552. if (var.value)
  553. {
  554. if (!line.empty())
  555. {
  556. exclusive = false;
  557. }
  558. line = var.name;
  559. }
  560. }
  561. if (!line.empty() && exclusive)
  562. {
  563. int &hits = exclusive_hits[line];
  564. if (hits)
  565. {
  566. hits += 1;
  567. }
  568. else
  569. {
  570. hits = 1;
  571. }
  572. h2_Hlt1_flags_excl_B_Mass->Fill(reconstructed_B_Mass, line.c_str(), 1);
  573. }
  574. }
  575. if (CutHlt1DecisionsOr({"TwoTrackMuonMVA", "TrackMuonMVA"}))
  576. {
  577. if (TMath::Abs(reconstructed_Kstar_Mass - KSTAR_MASS) < 50.)
  578. {
  579. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  580. {
  581. h1_B_Mass_jpsi->Fill(reconstructed_B_Mass);
  582. }
  583. else if (TMath::Abs(dimuon.M() - PSI2S_MASS) < 100.)
  584. {
  585. h1_B_Mass_psi2s->Fill(reconstructed_B_Mass);
  586. }
  587. }
  588. }
  589. PrintProgress(analysis_name, entries, 10000, i);
  590. }
  591. std::cout << "# Exclusive Hits" << std::endl;
  592. for (const auto &[line, hits] : exclusive_hits)
  593. {
  594. std::cout << line << ": " << hits << std::endl;
  595. }
  596. DrawInDefaultCanvas(h2_Hlt1_flags_B_Mass, analysis_name, 0.16, "COLZ");
  597. DrawInDefaultCanvas(h2_Hlt1_flags_excl_B_Mass, analysis_name, 0.16, "COLZ");
  598. DrawInDefaultCanvas(h1_B_Mass_jpsi, analysis_name, 0.1);
  599. DrawInDefaultCanvas(h1_B_Mass_psi2s, analysis_name, 0.1);
  600. auto roofit_hist_jpsi = CreateRooFitHistogram(h1_B_Mass_jpsi);
  601. auto roofit_hist_psi2s = CreateRooFitHistogram(h1_B_Mass_psi2s);
  602. DrawInDefaultCanvas(roofit_hist_jpsi, analysis_name, 0.1);
  603. DrawInDefaultCanvas(roofit_hist_psi2s, analysis_name, 0.1);
  604. return 0;
  605. }
  606. int new_analysis()
  607. {
  608. return analyze_Hlt2RD_B0ToKpPimMuMu();
  609. }