Study of the CMS Monte Carlo and Run open data for the dimuon system  1.0
A C++ software for the measurement of the Z resonance, the angle of the negative muon in the Collins–Soper frame of the dimuon system and the forward-backward asymmetry of Drell–Yan events in pp collisions at 8 TeV
utilities.h
Go to the documentation of this file.
1 
7 #ifndef UTILITIES_H
8 #define UTILITIES_H
9 
10 #include <ROOT/RDataFrame.hxx>
11 #include <ROOT/RVec.hxx>
12 #include <Math/Vector4Dfwd.h>
13 #include <Math/Vector4D.h>
14 
15 using namespace ROOT::VecOps;
16 
17 ROOT::Math::PtEtaPhiMVector quadrivectot(RVec<float>& pt, RVec<float>& eta, RVec<float>& phi, RVec<float>& mass){
18 
31  ROOT::Math::PtEtaPhiMVector m1(pt[0], eta[0], phi[0], mass[0]);
32  ROOT::Math::PtEtaPhiMVector m2(pt[1], eta[1], phi[1], mass[1]);
33 
34  return m1+m2;
35 
36 }
37 
38 
39 ROOT::Math::PtEtaPhiMVector quadrivec1(RVec<float>& pt, RVec<float>& eta, RVec<float>& phi, RVec<float>& mass){
40 
52  ROOT::Math::PtEtaPhiMVector m1(pt[0], eta[0], phi[0], mass[0]);
53 
54  return m1;
55 
56 }
57 
58 
59 ROOT::Math::PtEtaPhiMVector quadrivec2(RVec<float>& pt, RVec<float>& eta, RVec<float>& phi, RVec<float>& mass){
60 
72  ROOT::Math::PtEtaPhiMVector m2(pt[1], eta[1], phi[1], mass[1]);
73 
74  return m2;
75 
76 }
77 
78 
79 auto allquantities(ROOT::RDataFrame df){
80 
90  auto df_2mu_MC = df.Define("quadrivectot", quadrivectot,{"Muon_pt", "Muon_eta", "Muon_phi", "Muon_mass"});
91 
92  df_2mu_MC = df_2mu_MC.Define("quadrivec1", quadrivec1,{"Muon_pt", "Muon_eta", "Muon_phi", "Muon_mass"});
93 
94  df_2mu_MC = df_2mu_MC.Define("quadrivec2", quadrivec2,{"Muon_pt", "Muon_eta", "Muon_phi", "Muon_mass"});
95 
96  df_2mu_MC = df_2mu_MC.Define("dimuon_mass", "quadrivectot.mass()");
97 
98  df_2mu_MC = df_2mu_MC.Define("y","quadrivectot.Rapidity()");
99 
100  df_2mu_MC = df_2mu_MC.Define("ptll","quadrivectot.Pt()");
101 
102  df_2mu_MC = df_2mu_MC.Define("Pztot","quadrivectot.Pz()");
103 
104  df_2mu_MC = df_2mu_MC.Define("Pz1", "quadrivec1.Pz()");
105 
106  df_2mu_MC = df_2mu_MC.Define("Pz2", "quadrivec2.Pz()");
107 
108  df_2mu_MC = df_2mu_MC.Define("E", "quadrivectot.E()");
109 
110  df_2mu_MC = df_2mu_MC.Define("E1", "quadrivec1.E()");
111 
112  df_2mu_MC = df_2mu_MC.Define("E2", "quadrivec2.E()");
113 
114  df_2mu_MC = df_2mu_MC.Define("Pp1","(E1+Pz1)/pow(2,0.5)");
115 
116  df_2mu_MC = df_2mu_MC.Define("Pp2","(E2+Pz2)/pow(2,0.5)");
117 
118  df_2mu_MC = df_2mu_MC.Define("Pm1","(E1-Pz1)/pow(2,0.5)");
119 
120  df_2mu_MC = df_2mu_MC.Define("Pm2","(E2-Pz2)/pow(2,0.5)");
121 
122  df_2mu_MC = df_2mu_MC.Define("costheta", "2*(Pp1*Pm2-Pm1*Pp2)*Pztot/(pow(pow(dimuon_mass,2)*(pow(dimuon_mass,2)+pow(ptll,2)),0.5)*fabs(Pztot))");
123 
124  return df_2mu_MC;
125 }
126 #endif /* UTILITIES_H */
ROOT::Math::PtEtaPhiMVector quadrivec1(RVec< float > &pt, RVec< float > &eta, RVec< float > &phi, RVec< float > &mass)
Definition: utilities.h:39
ROOT::Math::PtEtaPhiMVector quadrivectot(RVec< float > &pt, RVec< float > &eta, RVec< float > &phi, RVec< float > &mass)
Definition: utilities.h:17
auto allquantities(ROOT::RDataFrame df)
Definition: utilities.h:79
ROOT::Math::PtEtaPhiMVector quadrivec2(RVec< float > &pt, RVec< float > &eta, RVec< float > &phi, RVec< float > &mass)
Definition: utilities.h:59