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.

255 lines
11 KiB

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