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.

178 lines
7.4 KiB

  1. #include <string>
  2. #include <iostream>
  3. #include <cmath>
  4. #include <algorithm>
  5. #include <filesystem>
  6. #include <string_view>
  7. #include "TH1D.h"
  8. #include "TH2D.h"
  9. #include "THStack.h"
  10. #include "TGraph.h"
  11. #include "TTree.h"
  12. #include "TChain.h"
  13. #include "TFile.h"
  14. #include "TCanvas.h"
  15. #include "TROOT.h"
  16. #include "TStyle.h"
  17. #include "TColor.h"
  18. #include "TLorentzVector.h"
  19. #include "TRandom3.h"
  20. #include "TLegend.h"
  21. #include "RooDataHist.h"
  22. #include "RooRealVar.h"
  23. #include "RooPlot.h"
  24. #include "RooGaussian.h"
  25. #include "RooExponential.h"
  26. #include "RooRealConstant.h"
  27. #include "RooAddPdf.h"
  28. #include "RooFitResult.h"
  29. #include "RooProduct.h"
  30. #include "RooCrystalBall.h"
  31. #include "RooBreitWigner.h"
  32. #include "RooArgSet.h"
  33. #include "RooFFTConvPdf.h"
  34. #include "RooNovosibirsk.h"
  35. #include "constants.h"
  36. #include "basic_analysis.h"
  37. #include "hlt1_decision_analysis.h"
  38. #include "bdt_classification.h"
  39. // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% Exclusive B2CC Bu To Kp Mu Mu %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
  40. int new_analysis_bu2jpsikplus()
  41. {
  42. const char *analysis_name = "BuToJpsiKplus";
  43. const char *data_tree_name = "Hlt2B2CC_BuToJpsiKplus_JpsiToMuMu_Detached";
  44. const char *sim_tree_name = "BuToKpMuMu_noPID";
  45. const char *end_state_mass_literal = "m(K^{+}#mu^{+}#mu^{-})";
  46. const bool retrain_bdt = false;
  47. TChain *data_chain = new TChain(TString::Format("%s/DecayTree", data_tree_name));
  48. data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/b2cc_turbo_superliaison_collision23_butojpsikplus_turbopass_magdown_v0r0p6205916.root");
  49. FourVect *l14v_data = FourVect::Init(data_chain, "muplus");
  50. FourVect *l24v_data = FourVect::Init(data_chain, "muminus");
  51. FourVect *hp4v_data = FourVect::Init(data_chain, "Kplus");
  52. Bool_t muplus_ISMUON, muminus_ISMUON, Hlt1TrackMVADecision, Hlt1TwoTrackMVADecision;
  53. data_chain->SetBranchAddress("muplus_ISMUON", &muplus_ISMUON);
  54. data_chain->SetBranchAddress("muminus_ISMUON", &muminus_ISMUON);
  55. data_chain->SetBranchAddress("Hlt1TrackMVADecision", &Hlt1TrackMVADecision);
  56. data_chain->SetBranchAddress("Hlt1TwoTrackMVADecision", &Hlt1TwoTrackMVADecision);
  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_v0r0p6657752_BuToKpMuMu_12143001_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. Int_t B_BKGCAT, L1_TRUEID, L2_TRUEID, Hp_TRUEID;
  63. sim_chain->SetBranchAddress("L1_TRUEID", &L1_TRUEID);
  64. sim_chain->SetBranchAddress("L2_TRUEID", &L2_TRUEID);
  65. sim_chain->SetBranchAddress("K_TRUEID", &Hp_TRUEID);
  66. sim_chain->SetBranchAddress("B_BKGCAT", &B_BKGCAT);
  67. Double_t B_Mass_jpsi_var, B_Mass_psi2s_var, B_Mass_sim_var;
  68. TString B_Mass_jpsi_var_name = "B_Mass_jpsi_var";
  69. TString B_Mass_psi2s_var_name = "B_Mass_psi2s_var";
  70. TString B_Mass_sim_var_name = "B_Mass_sim_var";
  71. TTree *tree_B_Mass_jpsi = new TTree("tree_B_Mass_jpsi", TString::Format("B Mass, J/#psi Mode (%s)", analysis_name));
  72. TTree *tree_B_Mass_psi2s = new TTree("tree_B_Mass_psi2s", TString::Format("B Mass, #psi(2S) Mode (%s)", analysis_name));
  73. TTree *tree_B_Mass_sim = new TTree("tree_B_Mass_sim", TString::Format("B Mass, Simualted (%s_noPID)", analysis_name));
  74. tree_B_Mass_jpsi->Branch(B_Mass_jpsi_var_name, &B_Mass_jpsi_var);
  75. tree_B_Mass_psi2s->Branch(B_Mass_psi2s_var_name, &B_Mass_psi2s_var);
  76. tree_B_Mass_sim->Branch(B_Mass_sim_var_name, &B_Mass_sim_var);
  77. TH1D *h1_B_Mass_unf = new TH1D("h1_B_Mass_unf", TString::Format("B Mass (%s), Unfiltered", data_tree_name), B_MASS_HIST_BINS, B_MASS_VAR_MIN, B_MASS_VAR_MAX);
  78. TH1D *h1_B_Mass_sim_unf = new TH1D("h1_B_Mass_sim_unf", TString::Format("B Mass, Simualted (%s), Unfiltered", sim_tree_name), B_MASS_HIST_BINS, B_MASS_VAR_MIN, B_MASS_VAR_MAX);
  79. TH2D *h2_Hlt1_flags_B_Mass = new TH2D("h2_Hlt1_flags_B_Mass", "Hlt1 Decision vs B Mass", 50, 5100, 5400, 13, 1., 14.);
  80. 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.);
  81. h1_B_Mass_unf->GetXaxis()->SetTitle(end_state_mass_literal);
  82. h1_B_Mass_sim_unf->GetXaxis()->SetTitle(end_state_mass_literal);
  83. h2_Hlt1_flags_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  84. h2_Hlt1_flags_excl_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
  85. // ConnectHlt1Decisions(data_chain, h2_Hlt1_flags_B_Mass, h2_Hlt1_flags_excl_B_Mass);
  86. // auto hlt1_decision_histos = CreateHlt1DecisionHistos(analysis_name);
  87. // std::map<std::string, int> exclusive_hits{};
  88. unsigned int data_entries = data_chain->GetEntries();
  89. unsigned int sim_entries = sim_chain->GetEntries();
  90. for (unsigned int i = 0; i < sim_entries; i++)
  91. {
  92. sim_chain->GetEntry(i);
  93. Double_t reconstructed_B_Mass = (hp4v_sim->LorentzVector() + l14v_sim->LorentzVector() + l24v_sim->LorentzVector()).M();
  94. h1_B_Mass_sim_unf->Fill(reconstructed_B_Mass);
  95. if (B_BKGCAT == 0 && TMath::Abs(L1_TRUEID) == PID_MUON && L2_TRUEID == -L1_TRUEID && TMath::Abs(Hp_TRUEID) == PID_KAON)
  96. {
  97. B_Mass_sim_var = reconstructed_B_Mass;
  98. tree_B_Mass_sim->Fill();
  99. }
  100. PrintProgress(TString::Format("%s SIG Collection", analysis_name), sim_entries, 10000, i);
  101. }
  102. std::cout << "# Starting evaluation of data." << std::endl;
  103. for (unsigned int i = 0; i < data_entries; i++)
  104. {
  105. data_chain->GetEntry(i);
  106. TLorentzVector dimuon = l14v_data->LorentzVector() + l24v_data->LorentzVector();
  107. Double_t reconstructed_B_Mass = (hp4v_data->LorentzVector() + dimuon).M();
  108. // if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  109. // {
  110. // CheckHlt1Decisioins(h2_Hlt1_flags_B_Mass, h2_Hlt1_flags_excl_B_Mass, exclusive_hits, reconstructed_B_Mass);
  111. // FillHlt1DecisionHistos(hlt1_decision_histos, reconstructed_B_Mass);
  112. // }
  113. h1_B_Mass_unf->Fill(reconstructed_B_Mass);
  114. if (muplus_ISMUON && muminus_ISMUON && (Hlt1TrackMVADecision || Hlt1TwoTrackMVADecision))
  115. {
  116. if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
  117. {
  118. B_Mass_jpsi_var = reconstructed_B_Mass;
  119. tree_B_Mass_jpsi->Fill();
  120. }
  121. else if (TMath::Abs(dimuon.M() - PSI2S_MASS) < 100.)
  122. {
  123. B_Mass_psi2s_var = reconstructed_B_Mass;
  124. tree_B_Mass_psi2s->Fill();
  125. }
  126. }
  127. PrintProgress(TString::Format("%s BDT Evaluation", analysis_name), data_entries, 10000, i);
  128. }
  129. // std::cout << "# Exclusive Hits" << std::endl;
  130. // for (const auto &[line, hits] : exclusive_hits)
  131. // {
  132. // std::cout << line << ": " << hits << std::endl;
  133. // }
  134. DrawInDefaultCanvas(h1_B_Mass_unf, analysis_name, 0.1);
  135. DrawInDefaultCanvas(h1_B_Mass_sim_unf, analysis_name, 0.1);
  136. auto roofit_hist_sim = CreateRooDataSetAndFitCB(tree_B_Mass_sim, B_Mass_sim_var_name, end_state_mass_literal, false, false, ShapeParamters{});
  137. 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, 5140., 5460.);
  138. 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, 5140., 5460.);
  139. DrawInDefaultCanvas(roofit_hist_jpsi_fitsum, analysis_name);
  140. DrawInDefaultCanvas(roofit_hist_psi2s_fitsum, analysis_name);
  141. DrawInDefaultCanvas(roofit_hist_sim, analysis_name);
  142. // DrawHlt1DecisionHistos(analysis_name, hlt1_decision_histos);
  143. return 0;
  144. }