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.

356 lines
14 KiB

6 months ago
6 months ago
6 months ago
6 months ago
6 months ago
6 months ago
6 months ago
6 months ago
6 months ago
7 months ago
6 months ago
6 months ago
6 months ago
6 months ago
6 months ago
6 months ago
6 months ago
6 months ago
6 months ago
  1. #include <string>
  2. #include <iostream>
  3. #include <cmath>
  4. #include <algorithm>
  5. #include <filesystem>
  6. #include <string_view>
  7. #include <ctime>
  8. #include <fstream>
  9. #include "TH1D.h"
  10. #include "TH2D.h"
  11. #include "THStack.h"
  12. #include "TGraph.h"
  13. #include "TTree.h"
  14. #include "TChain.h"
  15. #include "TFile.h"
  16. #include "TCanvas.h"
  17. #include "TROOT.h"
  18. #include "TStyle.h"
  19. #include "TColor.h"
  20. #include "TLorentzVector.h"
  21. #include "TRandom3.h"
  22. #include "TLegend.h"
  23. #include "RooDataHist.h"
  24. #include "RooRealVar.h"
  25. #include "RooPlot.h"
  26. #include "RooGaussian.h"
  27. #include "RooExponential.h"
  28. #include "RooRealConstant.h"
  29. #include "RooAddPdf.h"
  30. #include "RooFitResult.h"
  31. #include "RooProduct.h"
  32. #include "RooCrystalBall.h"
  33. #include "RooBreitWigner.h"
  34. #include "RooArgSet.h"
  35. #include "RooFFTConvPdf.h"
  36. #include "RooNovosibirsk.h"
  37. #include "constants.h"
  38. #include "basic_analysis.h"
  39. #include "hlt1_decision_analysis.h"
  40. #include "bdt_classification.h"
  41. // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% Bu To Hp Mu Mu %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
  42. int new_analysis_bu2hpmumu()
  43. {
  44. const char *analysis_name = "BuToHpMuMu";
  45. const char *data_tree_name = "SpruceRD_BuToHpMuMu";
  46. const char *sim_tree_name = "BuToHpMuMu_noPID_mapped";
  47. const char *end_state_mass_literal = "m(#pi^{+}_{(#rightarrow K^{+})}#mu^{+}#mu^{-}) [MeV]";
  48. const bool retrain_bdt = false;
  49. const bool skip_fit = false;
  50. TChain *data_chain = new TChain(TString::Format("%s/DecayTree", data_tree_name));
  51. data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/Collision23_Beam6800GeV-VeloClosed-MagDown-Excl-UT_RealData_Sprucing23r1_90000000_RD.root");
  52. // data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/Comm23_B6800GeV-VeCl-MgDo-Excl-UT_RealD_SprPass23r1_94000000_RD_45.root");
  53. FourVect *l14v_data = FourVect::Init(data_chain, "muminus");
  54. FourVect *l24v_data = FourVect::Init(data_chain, "muplus");
  55. FourVect *hp4v_data = FourVect::Init(data_chain, "Kplus");
  56. Double_t Kplus_PID_K;
  57. data_chain->SetBranchAddress("Kplus_PID_K", &Kplus_PID_K);
  58. TChain *sim_chain = new TChain(TString::Format("%s/DecayTree", sim_tree_name));
  59. sim_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/MC/BuToHpMuMu_mapped_mc.root");
  60. FourVect *l14v_sim = FourVect::Init(sim_chain, "muminus");
  61. FourVect *l24v_sim = FourVect::Init(sim_chain, "muplus");
  62. FourVect *hp4v_sim = FourVect::Init(sim_chain, "Kplus");
  63. Double_t B_Mass_jpsi_var, B_Mass_psi2s_var, B_Mass_sim_var;
  64. TString B_Mass_jpsi_var_name = "B_Mass_jpsi_var";
  65. TString B_Mass_psi2s_var_name = "B_Mass_psi2s_var";
  66. TString B_Mass_sim_var_name = "B_Mass_sim_var";
  67. TTree *tree_B_Mass_jpsi = new TTree("tree_B_Mass_jpsi", TString::Format("B^{#pm} Mass, J/#psi Mode (%s)", analysis_name));
  68. TTree *tree_B_Mass_psi2s = new TTree("tree_B_Mass_psi2s", TString::Format("B^{#pm} Mass, #psi(2S) Mode (%s)", analysis_name));
  69. TTree *tree_B_Mass_sim = new TTree("tree_B_Mass_sim", TString::Format("B^{#pm} Mass, Simualted (%s)", analysis_name));
  70. tree_B_Mass_jpsi->Branch(B_Mass_jpsi_var_name, &B_Mass_jpsi_var);
  71. tree_B_Mass_psi2s->Branch(B_Mass_psi2s_var_name, &B_Mass_psi2s_var);
  72. tree_B_Mass_sim->Branch(B_Mass_sim_var_name, &B_Mass_sim_var);
  73. TH1D *h1_B_Mass_unf = new TH1D("h1_B_Mass_unf", TString::Format("B^{#pm} Mass (%s), Unfiltered", analysis_name), B_MASS_HIST_BINS, B_MASS_VAR_MIN, B_MASS_VAR_MAX);
  74. TH1D *h1_B_Mass_bdtf = new TH1D("h1_B_Mass_bdtf", TString::Format("B^{#pm} Mass (%s), BDT Filter", analysis_name), B_MASS_HIST_BINS, B_MASS_VAR_MIN, B_MASS_VAR_MAX);
  75. TH2D *h2_Hlt1_flags_B_Mass = new TH2D("h2_Hlt1_flags_B_Mass", "Hlt1 Decision vs B^{#pm} Mass", 50, 5100, 5400, 13, 1., 14.);
  76. TH2D *h2_Hlt1_flags_excl_B_Mass = new TH2D("h2_Hlt1_flags_excl_B_Mass", "Excl Hlt1 Decision vs B^{#pm} Mass", 50, 5100, 5400, 13, 1., 14.);
  77. TH1D *h1_bdt_probs = new TH1D("h1_bdt_probs", "BDT Probabilities", 100, -2, 2);
  78. h1_B_Mass_unf->GetXaxis()->SetTitle(end_state_mass_literal);
  79. h1_B_Mass_bdtf->GetXaxis()->SetTitle(end_state_mass_literal);
  80. h2_Hlt1_flags_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  81. h2_Hlt1_flags_excl_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  82. ConnectHlt1Decisions(data_chain, h2_Hlt1_flags_B_Mass, h2_Hlt1_flags_excl_B_Mass);
  83. auto hlt1_decision_histos = CreateHlt1DecisionHistos(analysis_name);
  84. std::map<std::string, int> exclusive_hits{};
  85. // TV* kplus_pid_k_var = TV::Double("Kplus_PID_K", "Kplus_PID_K");
  86. std::vector<TV *> vars{
  87. TV::Float("B_PT", "B_PT"),
  88. // TV::Float("B_BPVFDCHI2", "B_BPVFDCHI2"),
  89. TV::Float("B_BPVDIRA", "B_BPVDIRA"),
  90. TV::Double("B_CHI2", "B_CHI2"),
  91. TV::Float("Jpsi_BPVIPCHI2", "Jpsi_BPVIPCHI2"),
  92. TV::Float("Jpsi_PT", "Jpsi_PT"),
  93. TV::Double("Jpsi_CHI2", "Jpsi_CHI2"),
  94. TV::Float("Kplus_BPVIPCHI2", "Kplus_BPVIPCHI2"),
  95. TV::Float("Kplus_PT", "Kplus_PT"),
  96. // kplus_pid_k_var,
  97. TV::Float("muminus_BPVIPCHI2", "muminus_BPVIPCHI2"),
  98. // TV::Float("muminus_PT", "muminus_PT"),
  99. TV::Float("muplus_BPVIPCHI2", "muplus_BPVIPCHI2"),
  100. // TV::Float("muplus_PT", "muplus_PT"),
  101. };
  102. TTree *sig_tree = new TTree("TreeS", "tree containing signal data");
  103. TTree *bkg_tree = new TTree("TreeB", "tree containing background data");
  104. ConnectVarsToData(vars, data_chain, sim_chain, sig_tree, bkg_tree);
  105. unsigned int data_entries = data_chain->GetEntries();
  106. unsigned int sim_entries = sim_chain->GetEntries();
  107. std::cout << "# Got " << data_entries << " data and " << sim_entries << " simulated events." << std::endl;
  108. unsigned int bkg_events = 0;
  109. unsigned int sig_events = 0;
  110. if (retrain_bdt)
  111. {
  112. std::cout << "# Starting with BDT retrain." << std::endl;
  113. for (unsigned int i = 0; i < data_entries; i++)
  114. {
  115. data_chain->GetEntry(i);
  116. TLorentzVector dimuon = l14v_data->LorentzVector() + l24v_data->LorentzVector();
  117. Double_t reconstructed_B_Mass = (hp4v_data->LorentzVector() + dimuon).M();
  118. if (std::all_of(vars.begin(), vars.end(), [](TV *v)
  119. { return v->IsDataFinite(); }))
  120. {
  121. if (reconstructed_B_Mass > 5500. && ((TMath::Abs(dimuon.M() - JPSI_MASS) < 100.) || (TMath::Abs(dimuon.M() - PSI2S_MASS) < 100.)) && Kplus_PID_K > -3)
  122. {
  123. bkg_tree->Fill();
  124. bkg_events++;
  125. }
  126. }
  127. PrintProgress(TString::Format("%s BKG Collection", analysis_name), data_entries, 10000, i);
  128. }
  129. std::cout << "# Added " << bkg_events << " background events." << std::endl;
  130. }
  131. else
  132. {
  133. std::cout << "# Starting without BDT retrain." << std::endl;
  134. bkg_events = data_entries;
  135. }
  136. for (unsigned int i = 0; i < sim_entries; i++)
  137. {
  138. sim_chain->GetEntry(i);
  139. Double_t reconstructed_B_Mass = (hp4v_sim->LorentzVector() + l14v_sim->LorentzVector() + l24v_sim->LorentzVector()).M();
  140. if (sig_events < bkg_events)
  141. {
  142. if (retrain_bdt && std::all_of(vars.begin(), vars.end(), [](TV *v)
  143. { return v->IsMCFinite(); }))
  144. {
  145. sig_tree->Fill();
  146. sig_events++;
  147. }
  148. }
  149. B_Mass_sim_var = reconstructed_B_Mass;
  150. tree_B_Mass_sim->Fill();
  151. PrintProgress(TString::Format("%s SIG Collection", analysis_name), sim_entries, 10000, i);
  152. }
  153. if (retrain_bdt)
  154. {
  155. std::cout << "# Added " << sig_events << " signal events." << std::endl;
  156. TrainBDT(vars, analysis_name, sig_tree, bkg_tree);
  157. std::cout << "# Finished BDT retrain." << std::endl;
  158. }
  159. if (skip_fit)
  160. {
  161. std::cout << "# Skipping evaluation of data." << std::endl;
  162. return 0;
  163. }
  164. std::cout << "# Starting evaluation of data." << std::endl;
  165. Float_t *train_vars = new Float_t[vars.size()];
  166. auto reader = SetupReader(vars, train_vars, analysis_name);
  167. const double mva_cut_value = 0.0;
  168. for (unsigned int i = 0; i < data_entries; i++)
  169. {
  170. data_chain->GetEntry(i);
  171. TLorentzVector dimuon = l14v_data->LorentzVector() + l24v_data->LorentzVector();
  172. Double_t reconstructed_B_Mass = (hp4v_data->LorentzVector() + dimuon).M();
  173. if (std::all_of(vars.begin(), vars.end(), [](TV *v)
  174. { return v->IsDataFinite(); }))
  175. {
  176. for (size_t j = 0; j < vars.size(); j++)
  177. {
  178. if (vars[j]->IsDouble())
  179. {
  180. train_vars[j] = vars[j]->GetDataDouble();
  181. }
  182. else if (vars[j]->IsFloat())
  183. {
  184. train_vars[j] = vars[j]->GetDataFloat();
  185. }
  186. }
  187. }
  188. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  189. {
  190. CheckHlt1Decisioins(h2_Hlt1_flags_B_Mass, h2_Hlt1_flags_excl_B_Mass, exclusive_hits, reconstructed_B_Mass);
  191. FillHlt1DecisionHistos(hlt1_decision_histos, reconstructed_B_Mass);
  192. }
  193. if (Kplus_PID_K > -3 && ((TMath::Abs(dimuon.M() - JPSI_MASS) < 100.) || (TMath::Abs(dimuon.M() - PSI2S_MASS) < 100.)))
  194. {
  195. double mva_response = reader->EvaluateMVA("BDT");
  196. h1_bdt_probs->Fill(mva_response);
  197. h1_B_Mass_unf->Fill(reconstructed_B_Mass);
  198. if (mva_response > mva_cut_value)
  199. {
  200. h1_B_Mass_bdtf->Fill(reconstructed_B_Mass);
  201. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  202. {
  203. B_Mass_jpsi_var = reconstructed_B_Mass;
  204. tree_B_Mass_jpsi->Fill();
  205. }
  206. else if (TMath::Abs(dimuon.M() - PSI2S_MASS) < 100.)
  207. {
  208. B_Mass_psi2s_var = reconstructed_B_Mass;
  209. tree_B_Mass_psi2s->Fill();
  210. }
  211. }
  212. }
  213. PrintProgress(TString::Format("%s BDT Evaluation", analysis_name), data_entries, 10000, i);
  214. }
  215. std::cout << "# Exclusive Hits" << std::endl;
  216. for (const auto &[line, hits] : exclusive_hits)
  217. {
  218. std::cout << line << ": " << hits << std::endl;
  219. }
  220. DrawInDefaultCanvas(h2_Hlt1_flags_B_Mass, analysis_name, 0.16, "COLZ");
  221. DrawInDefaultCanvas(h2_Hlt1_flags_excl_B_Mass, analysis_name, 0.16, "COLZ");
  222. DrawInDefaultCanvas(h1_B_Mass_unf, analysis_name, 0.1);
  223. DrawInDefaultCanvas(h1_B_Mass_bdtf, analysis_name, 0.1);
  224. DrawInDefaultCanvasStacked({h1_B_Mass_unf, h1_B_Mass_bdtf}, {kRed, kBlue}, {0, 3003}, analysis_name, TString::Format("B^{#pm} Mass (%s) Before/After BDT Selection", analysis_name));
  225. auto roofit_hist_sim = CreateRooDataSetAndFitCB(tree_B_Mass_sim, B_Mass_sim_var_name, end_state_mass_literal, false, false, ShapeParamters{});
  226. auto roofit_hist_jpsi_fitsum = CreateRooDataSetAndFitCB(tree_B_Mass_jpsi, B_Mass_jpsi_var_name, end_state_mass_literal, true, true, roofit_hist_sim.shape_parameters);
  227. roofit_hist_sim.shape_parameters.sigma_lr = roofit_hist_jpsi_fitsum.shape_parameters.sigma_lr;
  228. auto roofit_hist_psi2s_fitsum = CreateRooDataSetAndFitCB(tree_B_Mass_psi2s, B_Mass_psi2s_var_name, end_state_mass_literal, true, true, roofit_hist_sim.shape_parameters, true);
  229. DrawInDefaultCanvas(roofit_hist_jpsi_fitsum, analysis_name);
  230. DrawInDefaultCanvas(roofit_hist_psi2s_fitsum, analysis_name);
  231. DrawInDefaultCanvas(roofit_hist_sim, analysis_name);
  232. // DrawHlt1DecisionHistos(analysis_name, hlt1_decision_histos);
  233. DrawBDTProbs(h1_bdt_probs, mva_cut_value, analysis_name);
  234. auto signal_ratio = DivWithErr(roofit_hist_psi2s_fitsum.signal_yield.first, roofit_hist_psi2s_fitsum.signal_yield.second, roofit_hist_jpsi_fitsum.signal_yield.first, roofit_hist_jpsi_fitsum.signal_yield.second);
  235. std::time_t t = std::time(nullptr);
  236. std::tm tm = *std::localtime(&t);
  237. ofstream res_file;
  238. res_file.open(TString::Format("%s_results.txt", analysis_name).Data(), ios::out | ios::trunc);
  239. auto jpsi_sigobkg = DivWithErr(roofit_hist_jpsi_fitsum.signal_yield.first, roofit_hist_jpsi_fitsum.signal_yield.second, roofit_hist_jpsi_fitsum.background_yield.first, roofit_hist_jpsi_fitsum.background_yield.second);
  240. auto psi2s_sigobkg = DivWithErr(roofit_hist_psi2s_fitsum.signal_yield.first, roofit_hist_psi2s_fitsum.signal_yield.second, roofit_hist_psi2s_fitsum.background_yield.first, roofit_hist_psi2s_fitsum.background_yield.second);
  241. auto mumu_br_frac = DivWithErr(BRF_JPSI_MUMU_VAL, BRF_JPSI_MUMU_ERR, BRF_PSI2S_MUMU_VAL, BRF_PSI2S_MUMU_ERR);
  242. res_file << "#### " << analysis_name << " @ " << std::put_time(&tm, "%c") << " ####" << std::endl;
  243. res_file << "J/Psi Mode: " << ErrToStr(roofit_hist_jpsi_fitsum.signal_yield, 0) << " / " << ErrToStr(roofit_hist_jpsi_fitsum.background_yield, 0) << std::endl;
  244. res_file << " Sig/Bkg: " << ErrToStr(jpsi_sigobkg, 2) << std::endl;
  245. res_file << "Psi(2S) Mode: " << ErrToStr(roofit_hist_psi2s_fitsum.signal_yield, 0) << " / " << ErrToStr(roofit_hist_psi2s_fitsum.background_yield, 0) << std::endl;
  246. res_file << " Sig/Bkg: " << ErrToStr(psi2s_sigobkg, 2) << std::endl;
  247. res_file << "Mode Yield Ratio: " << ErrToStr(signal_ratio, 3) << std::endl;
  248. res_file << "Rel Br Frac MuMu: " << ErrToStr(mumu_br_frac, 3) << std::endl;
  249. auto br_frac = MultWithErr(signal_ratio.first, signal_ratio.second, mumu_br_frac.first, mumu_br_frac.second);
  250. res_file << "Rel Br Frac: " << ErrToStr(br_frac, 3) << std::endl
  251. << std::endl;
  252. res_file << "Params from Sim:" << std::endl
  253. << roofit_hist_sim.shape_parameters.ToString() << std::endl;
  254. for (const auto &par : roofit_hist_sim.fitted_params)
  255. {
  256. res_file << par.ToString(true).c_str() << std::endl;
  257. }
  258. res_file << std::endl
  259. << "Fitted Parameters: J/PSI" << std::endl
  260. << std::endl;
  261. for (const auto &par : roofit_hist_jpsi_fitsum.fitted_params)
  262. {
  263. res_file << par.ToString(true).c_str() << std::endl;
  264. }
  265. res_file << std::endl
  266. << "Fitted Parameters: PSI(2S)" << std::endl
  267. << std::endl;
  268. for (const auto &par : roofit_hist_psi2s_fitsum.fitted_params)
  269. {
  270. res_file << par.ToString(true).c_str() << std::endl;
  271. }
  272. auto print_table = [&res_file](std::string name, std::pair<double, double> sig, std::pair<double, double> bkg)
  273. {
  274. res_file << std::endl;
  275. res_file << "# " << name << std::endl;
  276. res_file << "\\begin{tabular}{c|c}" << std::endl;
  277. res_file << "$N_{Sig}$ & $N_{Bkg}$\\\\\\hline" << std::endl;
  278. res_file << TString::Format("$%.0f \\pm %.0f$ & $%.0f \\pm %.0f$", sig.first, sig.second, bkg.first, bkg.second).Data() << std::endl;
  279. res_file << "\\end{tabular}" << std::endl;
  280. };
  281. print_table("J/psi", roofit_hist_jpsi_fitsum.signal_yield, roofit_hist_jpsi_fitsum.background_yield);
  282. print_table("psi(2S)", roofit_hist_psi2s_fitsum.signal_yield, roofit_hist_psi2s_fitsum.background_yield);
  283. res_file.close();
  284. return 0;
  285. }