CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ZeeKinematicTools.cc
Go to the documentation of this file.
1 
3 
20 #include "TFile.h"
21 #include "TH1.h"
22 #include "TH2.h"
23 #include "TF1.h"
24 #include "TRandom.h"
25 
26 
27 #include <iostream>
28 #include <string>
29 #include <stdexcept>
30 #include <vector>
31 
33 
35 
37 
38 
39 //--------------------------------------------
40 
41 
42  float ZeeKinematicTools::cosThetaElectrons_SC( const std::pair<calib::CalibElectron*,calib::CalibElectron*>& aZCandidate, float ele1EnergyCorrection, float ele2EnergyCorrection ){
43 
44  float theta1 = 2. * atan( exp(- aZCandidate.first->getRecoElectron()->superCluster()->eta()) );
45  float phi1 = aZCandidate.first->getRecoElectron()->superCluster()->phi();
46 
47  float x1 = aZCandidate.first->getRecoElectron()->superCluster()->energy() * sin(theta1) * cos ( phi1 );
48  float y1 = aZCandidate.first->getRecoElectron()->superCluster()->energy() * sin(theta1) * sin ( phi1 );
49  float z1 = aZCandidate.first->getRecoElectron()->superCluster()->energy() * cos(theta1);
50  float mod1 = sqrt( x1*x1 + y1*y1 + z1*z1 );
51 
52  float theta2 = 2. * atan( exp(- aZCandidate.second->getRecoElectron()->superCluster()->eta()) );
53  float phi2 = aZCandidate.second->getRecoElectron()->superCluster()->phi();
54 
55  float x2 = aZCandidate.second->getRecoElectron()->superCluster()->energy() * sin(theta2) * cos ( phi2 );
56  float y2 = aZCandidate.second->getRecoElectron()->superCluster()->energy() * sin(theta2) * sin ( phi2 );
57  float z2 = aZCandidate.second->getRecoElectron()->superCluster()->energy() * cos(theta2);
58  float mod2 = sqrt( x2*x2 + y2*y2 + z2*z2 );
59 
60  return (x1*x2 + y1*y2 + z1*z2)/( mod1* mod2 );
61 
62 
63 }
64 
65 //--------------------------------------------
66 
67  float ZeeKinematicTools::cosThetaElectrons_TK( const std::pair<calib::CalibElectron*,calib::CalibElectron*>& aZCandidate, float ele1EnergyCorrection, float ele2EnergyCorrection ){
68 
69  float theta1 = 2. * atan( exp(- aZCandidate.first->getRecoElectron()->eta()) );
70  float phi1 = aZCandidate.first->getRecoElectron()->phi();
71 
72  float x1 = aZCandidate.first->getRecoElectron()->superCluster()->energy() * sin(theta1) * cos ( phi1 );
73  float y1 = aZCandidate.first->getRecoElectron()->superCluster()->energy() * sin(theta1) * sin ( phi1 );
74  float z1 = aZCandidate.first->getRecoElectron()->superCluster()->energy() * cos(theta1);
75  float mod1 = sqrt( x1*x1 + y1*y1 + z1*z1 );
76 
77  float theta2 = 2. * atan( exp(- aZCandidate.second->getRecoElectron()->eta()) );
78  float phi2 = aZCandidate.second->getRecoElectron()->phi();
79 
80  float x2 = aZCandidate.second->getRecoElectron()->superCluster()->energy() * sin(theta2) * cos ( phi2 );
81  float y2 = aZCandidate.second->getRecoElectron()->superCluster()->energy() * sin(theta2) * sin ( phi2 );
82  float z2 = aZCandidate.second->getRecoElectron()->superCluster()->energy() * cos(theta2);
83  float mod2 = sqrt( x2*x2 + y2*y2 + z2*z2 );
84 
85  return (x1*x2 + y1*y2 + z1*z2)/( mod1* mod2 );
86 }
87 
88 //--------------------------------------------
89 
90  float ZeeKinematicTools::calculateZMassWithCorrectedElectrons_noTK(const std::pair<calib::CalibElectron*,calib::CalibElectron*>& aZCandidate, float ele1EnergyCorrection, float ele2EnergyCorrection)
91 {
92  return ZIterativeAlgorithmWithFit::invMassCalc(aZCandidate.first->getRecoElectron()->superCluster()->energy() / ele1EnergyCorrection, aZCandidate.first->getRecoElectron()->superCluster()->eta(), aZCandidate.first->getRecoElectron()->superCluster()->phi(), aZCandidate.second->getRecoElectron()->superCluster()->energy() / ele2EnergyCorrection, aZCandidate.second->getRecoElectron()->superCluster()->eta(), aZCandidate.second->getRecoElectron()->superCluster()->phi());
93 
94 }
95 
96 //--------------------------------------------
97 
98  float ZeeKinematicTools::calculateZMassWithCorrectedElectrons_withTK(const std::pair<calib::CalibElectron*,calib::CalibElectron*>& aZCandidate, float ele1EnergyCorrection, float ele2EnergyCorrection)
99 {
100  return ZIterativeAlgorithmWithFit::invMassCalc(aZCandidate.first->getRecoElectron()->superCluster()->energy() / ele1EnergyCorrection, aZCandidate.first->getRecoElectron()->eta(), aZCandidate.first->getRecoElectron()->phi(), aZCandidate.second->getRecoElectron()->superCluster()->energy() / ele2EnergyCorrection, aZCandidate.second->getRecoElectron()->eta(), aZCandidate.second->getRecoElectron()->phi());
101 
102 }
103 
104 //--------------------------------------------
105 
106  float ZeeKinematicTools::calculateZMass_noTK(const std::pair<calib::CalibElectron*,calib::CalibElectron*>& aZCandidate)
107 {
108 
109  return ZIterativeAlgorithmWithFit::invMassCalc(aZCandidate.first->getRecoElectron()->superCluster()->energy(), aZCandidate.first->getRecoElectron()->superCluster()->eta(), aZCandidate.first->getRecoElectron()->superCluster()->phi(), aZCandidate.second->getRecoElectron()->superCluster()->energy(), aZCandidate.second->getRecoElectron()->superCluster()->eta(), aZCandidate.second->getRecoElectron()->superCluster()->phi());
110 
111 }
112 
113 //--------------------------------------------
114 
115  float ZeeKinematicTools::calculateZMass_withTK(const std::pair<calib::CalibElectron*,calib::CalibElectron*>& aZCandidate)
116 {
117 
118  return ZIterativeAlgorithmWithFit::invMassCalc(aZCandidate.first->getRecoElectron()->superCluster()->energy(), aZCandidate.first->getRecoElectron()->eta(), aZCandidate.first->getRecoElectron()->phi(), aZCandidate.second->getRecoElectron()->superCluster()->energy(), aZCandidate.second->getRecoElectron()->eta(), aZCandidate.second->getRecoElectron()->phi());
119 
120 }
121 
122 //--------------------------------------------
123 
124  float ZeeKinematicTools::calculateZRapidity(const std::pair<calib::CalibElectron*,calib::CalibElectron*>& aZCandidate)
125 {
126 
127  TLorentzVector ele1LV( aZCandidate.first->getRecoElectron()->px(), aZCandidate.first->getRecoElectron()->py(), aZCandidate.first->getRecoElectron()->pz(), aZCandidate.first->getRecoElectron()->superCluster()->energy());
128 
129  TLorentzVector ele2LV( aZCandidate.second->getRecoElectron()->px(), aZCandidate.second->getRecoElectron()->py(), aZCandidate.second->getRecoElectron()->pz(), aZCandidate.second->getRecoElectron()->superCluster()->energy());
130 
131 
132  return (ele1LV + ele2LV).Rapidity();
133 
134 }
135 
136 //--------------------------------------------
137 
138  float ZeeKinematicTools::calculateZEta(const std::pair<calib::CalibElectron*,calib::CalibElectron*>& aZCandidate)
139 {
140 
141  TLorentzVector ele1LV( aZCandidate.first->getRecoElectron()->px(), aZCandidate.first->getRecoElectron()->py(), aZCandidate.first->getRecoElectron()->pz(), aZCandidate.first->getRecoElectron()->superCluster()->energy());
142 
143  TLorentzVector ele2LV( aZCandidate.second->getRecoElectron()->px(), aZCandidate.second->getRecoElectron()->py(), aZCandidate.second->getRecoElectron()->pz(), aZCandidate.second->getRecoElectron()->superCluster()->energy());
144 
145  return (ele1LV + ele2LV).Eta();
146 
147 }
148 
149 //--------------------------------------------
150 
151  float ZeeKinematicTools::calculateZTheta(const std::pair<calib::CalibElectron*,calib::CalibElectron*>& aZCandidate)
152 {
153 
154  TLorentzVector ele1LV( aZCandidate.first->getRecoElectron()->px(), aZCandidate.first->getRecoElectron()->py(), aZCandidate.first->getRecoElectron()->pz(), aZCandidate.first->getRecoElectron()->superCluster()->energy());
155 
156  TLorentzVector ele2LV( aZCandidate.second->getRecoElectron()->px(), aZCandidate.second->getRecoElectron()->py(), aZCandidate.second->getRecoElectron()->pz(), aZCandidate.second->getRecoElectron()->superCluster()->energy());
157 
158  return (ele1LV + ele2LV).Theta();
159 
160 }
161 
162 //--------------------------------------------
163 
164  float ZeeKinematicTools::calculateZPhi(const std::pair<calib::CalibElectron*,calib::CalibElectron*>& aZCandidate)
165 {
166 
167  TLorentzVector ele1LV( aZCandidate.first->getRecoElectron()->px(), aZCandidate.first->getRecoElectron()->py(), aZCandidate.first->getRecoElectron()->pz(), aZCandidate.first->getRecoElectron()->superCluster()->energy());
168 
169  TLorentzVector ele2LV( aZCandidate.second->getRecoElectron()->px(), aZCandidate.second->getRecoElectron()->py(), aZCandidate.second->getRecoElectron()->pz(), aZCandidate.second->getRecoElectron()->superCluster()->energy());
170 
171  return (ele1LV + ele2LV).Phi();
172 
173 }
174 
175 //--------------------------------------------
176 
177  float ZeeKinematicTools::calculateZPt(const std::pair<calib::CalibElectron*,calib::CalibElectron*>& aZCandidate)
178 {
179 
180  TLorentzVector ele1LV( aZCandidate.first->getRecoElectron()->px(), aZCandidate.first->getRecoElectron()->py(), aZCandidate.first->getRecoElectron()->pz(), aZCandidate.first->getRecoElectron()->superCluster()->energy());
181 
182  TLorentzVector ele2LV( aZCandidate.second->getRecoElectron()->px(), aZCandidate.second->getRecoElectron()->py(), aZCandidate.second->getRecoElectron()->pz(), aZCandidate.second->getRecoElectron()->superCluster()->energy());
183 
184  return (ele1LV + ele2LV).Pt();
185 
186 }
187 
188 
189 
static float calculateZMassWithCorrectedElectrons_withTK(const std::pair< calib::CalibElectron *, calib::CalibElectron * > &aZCandidate, float ele1EnergyCorrection, float ele2EnergyCorrection)
static float cosThetaElectrons_TK(const std::pair< calib::CalibElectron *, calib::CalibElectron * > &aZCandidate, float ele1EnergyCorrection, float ele2EnergyCorrection)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
static float calculateZMass_withTK(const std::pair< calib::CalibElectron *, calib::CalibElectron * > &aZCandidate)
static float cosThetaElectrons_SC(const std::pair< calib::CalibElectron *, calib::CalibElectron * > &aZCandidate, float ele1EnergyCorrection, float ele2EnergyCorrection)
static float calculateZMassWithCorrectedElectrons_noTK(const std::pair< calib::CalibElectron *, calib::CalibElectron * > &aZCandidate, float ele1EnergyCorrection, float ele2EnergyCorrection)
static float calculateZEta(const std::pair< calib::CalibElectron *, calib::CalibElectron * > &aZCandidate)
static float calculateZPt(const std::pair< calib::CalibElectron *, calib::CalibElectron * > &aZCandidate)
T sqrt(T t)
Definition: SSEVec.h:48
static float calculateZRapidity(const std::pair< calib::CalibElectron *, calib::CalibElectron * > &aZCandidate)
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
static float calculateZMass_noTK(const std::pair< calib::CalibElectron *, calib::CalibElectron * > &aZCandidate)
static float invMassCalc(float Energy1, float Eta1, float Phi1, float Energy2, float Eta2, float Phi2)
static float calculateZTheta(const std::pair< calib::CalibElectron *, calib::CalibElectron * > &aZCandidate)
static float calculateZPhi(const std::pair< calib::CalibElectron *, calib::CalibElectron * > &aZCandidate)