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.

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