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.

199 lines
7.6 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. const int nBins = 70;
  33. const Double_t MASS_HIST_MIN = 5100.;
  34. const Double_t MASS_HIST_MAX = 6000.;
  35. const Double_t MASS_HIST_FIT_MIN = 5100.;
  36. const Double_t MASS_HIST_FIT_MAX = 6000.;
  37. const Double_t K_MASS = 493.677;
  38. void CreateRooFitAndDraw(TH1D *hist, int fitting_entries);
  39. const char* TITLE = "SpruceRD_BuToHpMuMu (#pi^{+} #rightarrow K^{+})";
  40. const char* FILE_NAME = "SpruceRD_BuToHpMuMu_Pip2Kp";
  41. const char* MASS_LITERAL = "m(#pi^{+}_{(#rightarrow K^{+})}#mu^{+}#mu^{-})";
  42. int analysis_fullstream_bu2hpmumu()
  43. {
  44. TChain *data_chain = new TChain("BuToHpMuMu/DecayTree");
  45. data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/spruce_magdown_2023_v0r1_tuple_90000000_2023_v0r0p6288631.root");
  46. data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/spruce_magdown_2023_v0_tuple_90000000_v0r0p6288631.root");
  47. Float_t B_BPVFDCHI2,
  48. B_BPVIPCHI2,
  49. L1_BPVIPCHI2,
  50. L2_BPVIPCHI2,
  51. L1_PT,
  52. L2_PT,
  53. Jpsi_BPVFDCHI2,
  54. Hp_PT,
  55. Hp_BPVIPCHI2,
  56. Hp_P;
  57. Double_t L1_PID_MU,
  58. L2_PID_MU,
  59. B_CHI2VXNDOF,
  60. Jpsi_MAXDOCACHI2,
  61. Jpsi_CHI2DOF,
  62. Hp_PID_K,
  63. Jpsi_M,
  64. B_M;
  65. Bool_t L1_ISMUON,
  66. L2_ISMUON,
  67. Hlt2RD_BuToKpMuMuDecision,
  68. Hlt2_InclDetDiMuon_4BodyDecision,
  69. Hlt2_InclDetDiMuon_3BodyDecision,
  70. Hlt2_InclDetDiMuonDecision,
  71. Hlt1TrackMVADecision,
  72. Hlt1TwoTrackMVADecision;
  73. data_chain->SetBranchAddress("B_M", &B_M);
  74. data_chain->SetBranchAddress("B_BPVFDCHI2", &B_BPVFDCHI2);
  75. data_chain->SetBranchAddress("B_BPVIPCHI2", &B_BPVIPCHI2);
  76. data_chain->SetBranchAddress("L1_BPVIPCHI2", &L1_BPVIPCHI2);
  77. data_chain->SetBranchAddress("L2_BPVIPCHI2", &L2_BPVIPCHI2);
  78. data_chain->SetBranchAddress("L1_PID_MU", &L1_PID_MU);
  79. data_chain->SetBranchAddress("L2_PID_MU", &L2_PID_MU);
  80. data_chain->SetBranchAddress("L1_ISMUON", &L1_ISMUON);
  81. data_chain->SetBranchAddress("L2_ISMUON", &L2_ISMUON);
  82. data_chain->SetBranchAddress("L1_PT", &L1_PT);
  83. data_chain->SetBranchAddress("L2_PT", &L2_PT);
  84. data_chain->SetBranchAddress("B_CHI2VXNDOF", &B_CHI2VXNDOF);
  85. data_chain->SetBranchAddress("Jpsi_MAXDOCACHI2", &Jpsi_MAXDOCACHI2);
  86. data_chain->SetBranchAddress("Jpsi_CHI2DOF", &Jpsi_CHI2DOF);
  87. data_chain->SetBranchAddress("Hp_PT", &Hp_PT);
  88. data_chain->SetBranchAddress("Hp_BPVIPCHI2", &Hp_BPVIPCHI2);
  89. data_chain->SetBranchAddress("Hp_P", &Hp_P);
  90. data_chain->SetBranchAddress("Hp_PID_K", &Hp_PID_K);
  91. data_chain->SetBranchAddress("Jpsi_BPVFDCHI2", &Jpsi_BPVFDCHI2);
  92. data_chain->SetBranchAddress("Jpsi_M", &Jpsi_M);
  93. data_chain->SetBranchAddress("Hlt2RD_BuToKpMuMuDecision", &Hlt2RD_BuToKpMuMuDecision);
  94. data_chain->SetBranchAddress("Hlt2_InclDetDiMuon_4BodyDecision", &Hlt2_InclDetDiMuon_4BodyDecision);
  95. data_chain->SetBranchAddress("Hlt2_InclDetDiMuon_3BodyDecision", &Hlt2_InclDetDiMuon_3BodyDecision);
  96. data_chain->SetBranchAddress("Hlt2_InclDetDiMuonDecision", &Hlt2_InclDetDiMuonDecision);
  97. data_chain->SetBranchAddress("Hlt1TrackMVADecision", &Hlt1TrackMVADecision);
  98. data_chain->SetBranchAddress("Hlt1TwoTrackMVADecision", &Hlt1TwoTrackMVADecision);
  99. Float_t L1_PX, L1_PY, L1_PZ, L1_ENERGY, L2_PX, L2_PY, L2_PZ, L2_ENERGY, Hp_PX, Hp_PY, Hp_PZ, Hp_ENERGY;
  100. data_chain->SetBranchAddress("L1_PX", &L1_PX);
  101. data_chain->SetBranchAddress("L1_PY", &L1_PY);
  102. data_chain->SetBranchAddress("L1_PZ", &L1_PZ);
  103. data_chain->SetBranchAddress("L1_ENERGY", &L1_ENERGY);
  104. data_chain->SetBranchAddress("L2_PX", &L2_PX);
  105. data_chain->SetBranchAddress("L2_PY", &L2_PY);
  106. data_chain->SetBranchAddress("L2_PZ", &L2_PZ);
  107. data_chain->SetBranchAddress("L2_ENERGY", &L2_ENERGY);
  108. data_chain->SetBranchAddress("Hp_PX", &Hp_PX);
  109. data_chain->SetBranchAddress("Hp_PY", &Hp_PY);
  110. data_chain->SetBranchAddress("Hp_PZ", &Hp_PZ);
  111. data_chain->SetBranchAddress("Hp_ENERGY", &Hp_ENERGY);
  112. TH1D *h1_B_M = new TH1D("h1_B_M", "B Mass", nBins, MASS_HIST_MIN, MASS_HIST_MAX);
  113. unsigned int entries = data_chain->GetEntries();
  114. unsigned int fitting_entries = 0;
  115. for (unsigned int i = 0; i < entries; i++)
  116. {
  117. data_chain->GetEntry(i);
  118. TVector3 K_momentum(Hp_PX, Hp_PY, Hp_PZ);
  119. double K_energy = TMath::Sqrt(TMath::Sq(K_MASS) + K_momentum.Mag2());
  120. TLorentzVector K_4v(K_momentum, K_energy);
  121. TLorentzVector l1_4v(L1_PX, L1_PY, L1_PZ, L1_ENERGY);
  122. TLorentzVector l2_4v(L2_PX, L2_PY, L2_PZ, L2_ENERGY);
  123. Double_t reconstructed_B_Mass = (K_4v + l1_4v + l2_4v).M();
  124. if ((((B_BPVFDCHI2 > 36) & (B_CHI2VXNDOF<16) & ( B_BPVIPCHI2<25) & (Jpsi_MAXDOCACHI2 < 36 )&(L1_BPVIPCHI2 > 9) & (L2_BPVIPCHI2 > 9) & (L1_PID_MU > -3) & (L2_PID_MU > -3)& (L1_ISMUON==1) & (L2_ISMUON==1) & (L1_PT>350) & (L2_PT>350) & (Jpsi_CHI2DOF<9 ) & (Jpsi_BPVFDCHI2>16) & (Jpsi_M<5500) & (Hp_PT>400) &( Hp_BPVIPCHI2>6) & (Hp_PT>400) & (Hp_BPVIPCHI2>6) & (Hp_P>2000) & (Hp_PID_K > -4)) & (Hlt2RD_BuToKpMuMuDecision==1) &(((Hlt2_InclDetDiMuon_4BodyDecision==1) | (Hlt2_InclDetDiMuon_3BodyDecision==1) | (Hlt2_InclDetDiMuonDecision==1))))&(TMath::Abs(Jpsi_M - 3096.9) < 100) & ((Hlt1TrackMVADecision==1) | (Hlt1TwoTrackMVADecision==1)))
  125. {
  126. h1_B_M->Fill(reconstructed_B_Mass);
  127. if (MASS_HIST_FIT_MIN <= reconstructed_B_Mass && reconstructed_B_Mass <= MASS_HIST_FIT_MAX) {
  128. fitting_entries++;
  129. }
  130. }
  131. if (i % 10000 == 0)
  132. {
  133. std::cout << "["
  134. << "SpruceRD_BuToHpMuMu"
  135. << "] Processed event: " << i << " (" << TString::Format("%.2f", ((double)i / (double)entries) * 100.) << "%)" << std::endl;
  136. }
  137. }
  138. TCanvas *c1 = new TCanvas("c1", "c1", 0, 0, 800, 600);
  139. h1_B_M->Draw();
  140. c1->Draw();
  141. c1->SaveAs(TString::Format("images/root_hist_%s_bmass.png", FILE_NAME).Data());
  142. CreateRooFitAndDraw(h1_B_M, fitting_entries);
  143. return 0;
  144. }
  145. void CreateRooFitAndDraw(TH1D *hist, int fitting_entries)
  146. {
  147. RooRealVar roo_var_mass("var_mass", "B Mass Variable", MASS_HIST_FIT_MIN, MASS_HIST_FIT_MAX);
  148. roo_var_mass.setRange("fitting_range", MASS_HIST_FIT_MIN, MASS_HIST_FIT_MAX);
  149. RooDataHist roohist_B_M("roohist_B_M", "B Mass Histogram", roo_var_mass, RooFit::Import(*hist));
  150. RooPlot *roo_frame_mass = roo_var_mass.frame(RooFit::Title(TITLE));
  151. roohist_B_M.plotOn(roo_frame_mass, RooFit::Binning(nBins), RooFit::Name("B Mass Distribution"));
  152. roo_frame_mass->GetXaxis()->SetTitle(MASS_LITERAL);
  153. TCanvas *c = new TCanvas("roofit_c", "roofit_c", 0, 0, 800, 600);
  154. roo_frame_mass->Draw();
  155. TLegend *leg1 = new TLegend(0.65, 0.7, 0.87, 0.8);
  156. leg1->SetFillColor(kWhite);
  157. leg1->SetLineColor(kBlack);
  158. leg1->AddEntry((TObject*)0, TString::Format("Entries: %d", fitting_entries), "");
  159. leg1->Draw();
  160. c->Draw();
  161. c->SaveAs(TString::Format("images/roofit_hist_%s_bmass.png", FILE_NAME).Data());
  162. }