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.

251 lines
11 KiB

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
  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. // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% Exclusive Bu To Kp Mu Mu %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
  42. int new_analysis_bu2kpmumu()
  43. {
  44. const char *analysis_name = "BuToKpMuMu";
  45. const char *data_tree_name = "Hlt2RD_BuToKpMuMu";
  46. const char *sim_tree_name = "BuToKpMuMu_noPID";
  47. const char *end_state_mass_literal = "m(K^{+}#mu^{+}#mu^{-}) [MeV]";
  48. TChain *data_chain = new TChain(TString::Format("%s/DecayTree", data_tree_name));
  49. data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/Collision23_Beam6800GeV-VeloClosed-MagDown-Excl-UT_RealData_SprucingPass23r1_94000000_RD.root");
  50. FourVect *l14v_data = FourVect::Init(data_chain, "muplus");
  51. FourVect *l24v_data = FourVect::Init(data_chain, "muminus");
  52. FourVect *hp4v_data = FourVect::Init(data_chain, "Kplus");
  53. Double_t Kplus_PID_K;
  54. data_chain->SetBranchAddress("Kplus_PID_K", &Kplus_PID_K);
  55. TChain *sim_chain = new TChain(TString::Format("%s/DecayTree", sim_tree_name));
  56. sim_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/MC/rd_btoxll_simulation_turbo_v0r0p6657752_BuToKpMuMu_12143001_magdown.root");
  57. FourVect *l14v_sim = FourVect::Init(sim_chain, "L1");
  58. FourVect *l24v_sim = FourVect::Init(sim_chain, "L2");
  59. FourVect *hp4v_sim = FourVect::Init(sim_chain, "K");
  60. Int_t B_BKGCAT, L1_TRUEID, L2_TRUEID, Hp_TRUEID;
  61. Double_t L1_PID_MU, L2_PID_MU, K_PID_K;
  62. sim_chain->SetBranchAddress("L1_TRUEID", &L1_TRUEID);
  63. sim_chain->SetBranchAddress("L2_TRUEID", &L2_TRUEID);
  64. sim_chain->SetBranchAddress("K_TRUEID", &Hp_TRUEID);
  65. sim_chain->SetBranchAddress("B_BKGCAT", &B_BKGCAT);
  66. sim_chain->SetBranchAddress("L1_PID_MU", &L1_PID_MU);
  67. sim_chain->SetBranchAddress("L2_PID_MU", &L2_PID_MU);
  68. sim_chain->SetBranchAddress("K_PID_K", &K_PID_K);
  69. Double_t B_Mass_jpsi_var, B_Mass_psi2s_var, B_Mass_sim_var;
  70. TString B_Mass_jpsi_var_name = "B_Mass_jpsi_var";
  71. TString B_Mass_psi2s_var_name = "B_Mass_psi2s_var";
  72. TString B_Mass_sim_var_name = "B_Mass_sim_var";
  73. TTree *tree_B_Mass_jpsi = new TTree("tree_B_Mass_jpsi", TString::Format("B^{#pm} Mass, J/#psi Mode (%s)", analysis_name));
  74. TTree *tree_B_Mass_psi2s = new TTree("tree_B_Mass_psi2s", TString::Format("B^{#pm} Mass, #psi(2S) Mode (%s)", analysis_name));
  75. TTree *tree_B_Mass_sim = new TTree("tree_B_Mass_sim", TString::Format("B^{#pm} Mass, Simualted (%s_noPID)", analysis_name));
  76. tree_B_Mass_jpsi->Branch(B_Mass_jpsi_var_name, &B_Mass_jpsi_var);
  77. tree_B_Mass_psi2s->Branch(B_Mass_psi2s_var_name, &B_Mass_psi2s_var);
  78. tree_B_Mass_sim->Branch(B_Mass_sim_var_name, &B_Mass_sim_var);
  79. 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);
  80. TH1D *h1_B_Mass_sim_unf = new TH1D("h1_B_Mass_sim_unf", TString::Format("B^{#pm} Mass, Simualted (%s), Unfiltered", sim_tree_name), B_MASS_HIST_BINS, B_MASS_VAR_MIN, B_MASS_VAR_MAX);
  81. 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.);
  82. 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.);
  83. h1_B_Mass_unf->GetXaxis()->SetTitle(end_state_mass_literal);
  84. h1_B_Mass_sim_unf->GetXaxis()->SetTitle(end_state_mass_literal);
  85. h2_Hlt1_flags_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  86. h2_Hlt1_flags_excl_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  87. ConnectHlt1Decisions(data_chain, h2_Hlt1_flags_B_Mass, h2_Hlt1_flags_excl_B_Mass);
  88. auto hlt1_decision_histos = CreateHlt1DecisionHistos(analysis_name);
  89. std::map<std::string, int> exclusive_hits{};
  90. unsigned int data_entries = data_chain->GetEntries();
  91. unsigned int sim_entries = sim_chain->GetEntries();
  92. for (unsigned int i = 0; i < sim_entries; i++)
  93. {
  94. sim_chain->GetEntry(i);
  95. Double_t reconstructed_B_Mass = (hp4v_sim->LorentzVector() + l14v_sim->LorentzVector() + l24v_sim->LorentzVector()).M();
  96. h1_B_Mass_sim_unf->Fill(reconstructed_B_Mass);
  97. if (B_BKGCAT == 0 && TMath::Abs(L1_TRUEID) == PID_MUON && L2_TRUEID == -L1_TRUEID && TMath::Abs(Hp_TRUEID) == PID_KAON && L1_PID_MU > -4. && L2_PID_MU > -4. && K_PID_K > -4.)
  98. {
  99. B_Mass_sim_var = reconstructed_B_Mass;
  100. tree_B_Mass_sim->Fill();
  101. }
  102. PrintProgress(TString::Format("%s SIG Collection", analysis_name), sim_entries, 10000, i);
  103. }
  104. std::cout << "# Starting evaluation of data." << std::endl;
  105. for (unsigned int i = 0; i < data_entries; i++)
  106. {
  107. data_chain->GetEntry(i);
  108. TLorentzVector dimuon = l14v_data->LorentzVector() + l24v_data->LorentzVector();
  109. Double_t reconstructed_B_Mass = (hp4v_data->LorentzVector() + dimuon).M();
  110. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  111. {
  112. CheckHlt1Decisioins(h2_Hlt1_flags_B_Mass, h2_Hlt1_flags_excl_B_Mass, exclusive_hits, reconstructed_B_Mass);
  113. FillHlt1DecisionHistos(hlt1_decision_histos, reconstructed_B_Mass);
  114. }
  115. h1_B_Mass_unf->Fill(reconstructed_B_Mass);
  116. if (Kplus_PID_K > -3)
  117. {
  118. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  119. {
  120. B_Mass_jpsi_var = reconstructed_B_Mass;
  121. tree_B_Mass_jpsi->Fill();
  122. }
  123. else if (TMath::Abs(dimuon.M() - PSI2S_MASS) < 100.)
  124. {
  125. B_Mass_psi2s_var = reconstructed_B_Mass;
  126. tree_B_Mass_psi2s->Fill();
  127. }
  128. }
  129. PrintProgress(TString::Format("%s BDT Evaluation", analysis_name), data_entries, 10000, i);
  130. }
  131. std::cout << "# Exclusive Hits" << std::endl;
  132. for (const auto &[line, hits] : exclusive_hits)
  133. {
  134. std::cout << line << ": " << hits << std::endl;
  135. }
  136. DrawInDefaultCanvas(h1_B_Mass_unf, analysis_name, 0.1);
  137. DrawInDefaultCanvas(h1_B_Mass_sim_unf, analysis_name, 0.1);
  138. auto roofit_hist_sim = CreateRooDataSetAndFitCB(tree_B_Mass_sim, B_Mass_sim_var_name, end_state_mass_literal, false, false, ShapeParamters{});
  139. 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);
  140. roofit_hist_sim.shape_parameters.sigma_lr = roofit_hist_jpsi_fitsum.shape_parameters.sigma_lr;
  141. 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);
  142. DrawInDefaultCanvas(roofit_hist_jpsi_fitsum, analysis_name);
  143. DrawInDefaultCanvas(roofit_hist_psi2s_fitsum, analysis_name);
  144. DrawInDefaultCanvas(roofit_hist_sim, analysis_name);
  145. // DrawHlt1DecisionHistos(analysis_name, hlt1_decision_histos);
  146. 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);
  147. std::time_t t = std::time(nullptr);
  148. std::tm tm = *std::localtime(&t);
  149. ofstream res_file;
  150. res_file.open(TString::Format("%s_results.txt", analysis_name).Data(), ios::out | ios::trunc);
  151. 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);
  152. 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);
  153. auto mumu_br_frac = DivWithErr(BRF_JPSI_MUMU_VAL, BRF_JPSI_MUMU_ERR, BRF_PSI2S_MUMU_VAL, BRF_PSI2S_MUMU_ERR);
  154. res_file << "#### " << analysis_name << " @ " << std::put_time(&tm, "%c") << " ####" << std::endl;
  155. res_file << "J/Psi Mode: " << ErrToStr(roofit_hist_jpsi_fitsum.signal_yield, 0) << " / " << ErrToStr(roofit_hist_jpsi_fitsum.background_yield, 0) << std::endl;
  156. res_file << " Sig/Bkg: " << ErrToStr(jpsi_sigobkg, 2) << std::endl;
  157. res_file << "Psi(2S) Mode: " << ErrToStr(roofit_hist_psi2s_fitsum.signal_yield, 0) << " / " << ErrToStr(roofit_hist_psi2s_fitsum.background_yield, 0) << std::endl;
  158. res_file << " Sig/Bkg: " << ErrToStr(psi2s_sigobkg, 2) << std::endl;
  159. res_file << "Mode Yield Ratio: " << ErrToStr(signal_ratio, 3) << std::endl;
  160. res_file << "Rel Br Frac MuMu: " << ErrToStr(mumu_br_frac, 3) << std::endl;
  161. auto br_frac = MultWithErr(signal_ratio.first, signal_ratio.second, mumu_br_frac.first, mumu_br_frac.second);
  162. res_file << "Rel Br Frac: " << ErrToStr(br_frac, 3) << std::endl
  163. << std::endl;
  164. res_file << "Params from Sim:" << std::endl
  165. << roofit_hist_sim.shape_parameters.ToString() << std::endl;
  166. for (const auto &par : roofit_hist_sim.fitted_params)
  167. {
  168. res_file << par.ToString(true).c_str() << std::endl;
  169. }
  170. res_file << std::endl
  171. << "Fitted Parameters: J/PSI" << std::endl
  172. << std::endl;
  173. for (const auto &par : roofit_hist_jpsi_fitsum.fitted_params)
  174. {
  175. res_file << par.ToString(true).c_str() << std::endl;
  176. }
  177. res_file << std::endl
  178. << "Fitted Parameters: PSI(2S)" << std::endl
  179. << std::endl;
  180. for (const auto &par : roofit_hist_psi2s_fitsum.fitted_params)
  181. {
  182. res_file << par.ToString(true).c_str() << std::endl;
  183. }
  184. auto print_table = [&res_file](std::string name, std::pair<double, double> sig, std::pair<double, double> bkg)
  185. {
  186. res_file << std::endl;
  187. res_file << "# " << name << std::endl;
  188. res_file << "\\begin{tabular}{c|c}" << std::endl;
  189. res_file << "$N_{Sig}$ & $N_{Bkg}$\\\\\\hline" << std::endl;
  190. res_file << TString::Format("$%.0f \\pm %.0f$ & $%.0f \\pm %.0f$", sig.first, sig.second, bkg.first, bkg.second).Data() << std::endl;
  191. res_file << "\\end{tabular}" << std::endl;
  192. };
  193. print_table("J/psi", roofit_hist_jpsi_fitsum.signal_yield, roofit_hist_jpsi_fitsum.background_yield);
  194. print_table("psi(2S)", roofit_hist_psi2s_fitsum.signal_yield, roofit_hist_psi2s_fitsum.background_yield);
  195. res_file.close();
  196. std::cout << "MC Truth Matched Sim Events:" << tree_B_Mass_sim->GetEntries() << std::endl;
  197. return 0;
  198. }