CMS 3D CMS Logo

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 #include <iostream>
27 #include <string>
28 #include <stdexcept>
29 #include <vector>
30 
32 
34 
36 
37 //--------------------------------------------
38 
39 float ZeeKinematicTools::cosThetaElectrons_SC(const std::pair<calib::CalibElectron*, calib::CalibElectron*>& aZCandidate,
40  float ele1EnergyCorrection,
41  float ele2EnergyCorrection) {
42  float theta1 = 2. * atan(exp(-aZCandidate.first->getRecoElectron()->superCluster()->eta()));
43  float phi1 = aZCandidate.first->getRecoElectron()->superCluster()->phi();
44 
45  float x1 = aZCandidate.first->getRecoElectron()->superCluster()->energy() * sin(theta1) * cos(phi1);
46  float y1 = aZCandidate.first->getRecoElectron()->superCluster()->energy() * sin(theta1) * sin(phi1);
47  float z1 = aZCandidate.first->getRecoElectron()->superCluster()->energy() * cos(theta1);
48  float mod1 = sqrt(x1 * x1 + y1 * y1 + z1 * z1);
49 
50  float theta2 = 2. * atan(exp(-aZCandidate.second->getRecoElectron()->superCluster()->eta()));
51  float phi2 = aZCandidate.second->getRecoElectron()->superCluster()->phi();
52 
53  float x2 = aZCandidate.second->getRecoElectron()->superCluster()->energy() * sin(theta2) * cos(phi2);
54  float y2 = aZCandidate.second->getRecoElectron()->superCluster()->energy() * sin(theta2) * sin(phi2);
55  float z2 = aZCandidate.second->getRecoElectron()->superCluster()->energy() * cos(theta2);
56  float mod2 = sqrt(x2 * x2 + y2 * y2 + z2 * z2);
57 
58  return (x1 * x2 + y1 * y2 + z1 * z2) / (mod1 * mod2);
59 }
60 
61 //--------------------------------------------
62 
63 float ZeeKinematicTools::cosThetaElectrons_TK(const std::pair<calib::CalibElectron*, calib::CalibElectron*>& aZCandidate,
64  float ele1EnergyCorrection,
65  float ele2EnergyCorrection) {
66  float theta1 = 2. * atan(exp(-aZCandidate.first->getRecoElectron()->eta()));
67  float phi1 = aZCandidate.first->getRecoElectron()->phi();
68 
69  float x1 = aZCandidate.first->getRecoElectron()->superCluster()->energy() * sin(theta1) * cos(phi1);
70  float y1 = aZCandidate.first->getRecoElectron()->superCluster()->energy() * sin(theta1) * sin(phi1);
71  float z1 = aZCandidate.first->getRecoElectron()->superCluster()->energy() * cos(theta1);
72  float mod1 = sqrt(x1 * x1 + y1 * y1 + z1 * z1);
73 
74  float theta2 = 2. * atan(exp(-aZCandidate.second->getRecoElectron()->eta()));
75  float phi2 = aZCandidate.second->getRecoElectron()->phi();
76 
77  float x2 = aZCandidate.second->getRecoElectron()->superCluster()->energy() * sin(theta2) * cos(phi2);
78  float y2 = aZCandidate.second->getRecoElectron()->superCluster()->energy() * sin(theta2) * sin(phi2);
79  float z2 = aZCandidate.second->getRecoElectron()->superCluster()->energy() * cos(theta2);
80  float mod2 = sqrt(x2 * x2 + y2 * y2 + z2 * z2);
81 
82  return (x1 * x2 + y1 * y2 + z1 * z2) / (mod1 * mod2);
83 }
84 
85 //--------------------------------------------
86 
88  const std::pair<calib::CalibElectron*, calib::CalibElectron*>& aZCandidate,
89  float ele1EnergyCorrection,
90  float ele2EnergyCorrection) {
92  aZCandidate.first->getRecoElectron()->superCluster()->energy() / ele1EnergyCorrection,
93  aZCandidate.first->getRecoElectron()->superCluster()->eta(),
94  aZCandidate.first->getRecoElectron()->superCluster()->phi(),
95  aZCandidate.second->getRecoElectron()->superCluster()->energy() / ele2EnergyCorrection,
96  aZCandidate.second->getRecoElectron()->superCluster()->eta(),
97  aZCandidate.second->getRecoElectron()->superCluster()->phi());
98 }
99 
100 //--------------------------------------------
101 
103  const std::pair<calib::CalibElectron*, calib::CalibElectron*>& aZCandidate,
104  float ele1EnergyCorrection,
105  float ele2EnergyCorrection) {
107  aZCandidate.first->getRecoElectron()->superCluster()->energy() / ele1EnergyCorrection,
108  aZCandidate.first->getRecoElectron()->eta(),
109  aZCandidate.first->getRecoElectron()->phi(),
110  aZCandidate.second->getRecoElectron()->superCluster()->energy() / ele2EnergyCorrection,
111  aZCandidate.second->getRecoElectron()->eta(),
112  aZCandidate.second->getRecoElectron()->phi());
113 }
114 
115 //--------------------------------------------
116 
118  const std::pair<calib::CalibElectron*, calib::CalibElectron*>& aZCandidate) {
119  return ZIterativeAlgorithmWithFit::invMassCalc(aZCandidate.first->getRecoElectron()->superCluster()->energy(),
120  aZCandidate.first->getRecoElectron()->superCluster()->eta(),
121  aZCandidate.first->getRecoElectron()->superCluster()->phi(),
122  aZCandidate.second->getRecoElectron()->superCluster()->energy(),
123  aZCandidate.second->getRecoElectron()->superCluster()->eta(),
124  aZCandidate.second->getRecoElectron()->superCluster()->phi());
125 }
126 
127 //--------------------------------------------
128 
130  const std::pair<calib::CalibElectron*, calib::CalibElectron*>& aZCandidate) {
131  return ZIterativeAlgorithmWithFit::invMassCalc(aZCandidate.first->getRecoElectron()->superCluster()->energy(),
132  aZCandidate.first->getRecoElectron()->eta(),
133  aZCandidate.first->getRecoElectron()->phi(),
134  aZCandidate.second->getRecoElectron()->superCluster()->energy(),
135  aZCandidate.second->getRecoElectron()->eta(),
136  aZCandidate.second->getRecoElectron()->phi());
137 }
138 
139 //--------------------------------------------
140 
141 float ZeeKinematicTools::calculateZRapidity(const std::pair<calib::CalibElectron*, calib::CalibElectron*>& aZCandidate) {
142  TLorentzVector ele1LV(aZCandidate.first->getRecoElectron()->px(),
143  aZCandidate.first->getRecoElectron()->py(),
144  aZCandidate.first->getRecoElectron()->pz(),
145  aZCandidate.first->getRecoElectron()->superCluster()->energy());
146 
147  TLorentzVector ele2LV(aZCandidate.second->getRecoElectron()->px(),
148  aZCandidate.second->getRecoElectron()->py(),
149  aZCandidate.second->getRecoElectron()->pz(),
150  aZCandidate.second->getRecoElectron()->superCluster()->energy());
151 
152  return (ele1LV + ele2LV).Rapidity();
153 }
154 
155 //--------------------------------------------
156 
157 float ZeeKinematicTools::calculateZEta(const std::pair<calib::CalibElectron*, calib::CalibElectron*>& aZCandidate) {
158  TLorentzVector ele1LV(aZCandidate.first->getRecoElectron()->px(),
159  aZCandidate.first->getRecoElectron()->py(),
160  aZCandidate.first->getRecoElectron()->pz(),
161  aZCandidate.first->getRecoElectron()->superCluster()->energy());
162 
163  TLorentzVector ele2LV(aZCandidate.second->getRecoElectron()->px(),
164  aZCandidate.second->getRecoElectron()->py(),
165  aZCandidate.second->getRecoElectron()->pz(),
166  aZCandidate.second->getRecoElectron()->superCluster()->energy());
167 
168  return (ele1LV + ele2LV).Eta();
169 }
170 
171 //--------------------------------------------
172 
173 float ZeeKinematicTools::calculateZTheta(const std::pair<calib::CalibElectron*, calib::CalibElectron*>& aZCandidate) {
174  TLorentzVector ele1LV(aZCandidate.first->getRecoElectron()->px(),
175  aZCandidate.first->getRecoElectron()->py(),
176  aZCandidate.first->getRecoElectron()->pz(),
177  aZCandidate.first->getRecoElectron()->superCluster()->energy());
178 
179  TLorentzVector ele2LV(aZCandidate.second->getRecoElectron()->px(),
180  aZCandidate.second->getRecoElectron()->py(),
181  aZCandidate.second->getRecoElectron()->pz(),
182  aZCandidate.second->getRecoElectron()->superCluster()->energy());
183 
184  return (ele1LV + ele2LV).Theta();
185 }
186 
187 //--------------------------------------------
188 
189 float ZeeKinematicTools::calculateZPhi(const std::pair<calib::CalibElectron*, calib::CalibElectron*>& aZCandidate) {
190  TLorentzVector ele1LV(aZCandidate.first->getRecoElectron()->px(),
191  aZCandidate.first->getRecoElectron()->py(),
192  aZCandidate.first->getRecoElectron()->pz(),
193  aZCandidate.first->getRecoElectron()->superCluster()->energy());
194 
195  TLorentzVector ele2LV(aZCandidate.second->getRecoElectron()->px(),
196  aZCandidate.second->getRecoElectron()->py(),
197  aZCandidate.second->getRecoElectron()->pz(),
198  aZCandidate.second->getRecoElectron()->superCluster()->energy());
199 
200  return (ele1LV + ele2LV).Phi();
201 }
202 
203 //--------------------------------------------
204 
205 float ZeeKinematicTools::calculateZPt(const std::pair<calib::CalibElectron*, calib::CalibElectron*>& aZCandidate) {
206  TLorentzVector ele1LV(aZCandidate.first->getRecoElectron()->px(),
207  aZCandidate.first->getRecoElectron()->py(),
208  aZCandidate.first->getRecoElectron()->pz(),
209  aZCandidate.first->getRecoElectron()->superCluster()->energy());
210 
211  TLorentzVector ele2LV(aZCandidate.second->getRecoElectron()->px(),
212  aZCandidate.second->getRecoElectron()->py(),
213  aZCandidate.second->getRecoElectron()->pz(),
214  aZCandidate.second->getRecoElectron()->superCluster()->energy());
215 
216  return (ele1LV + ele2LV).Pt();
217 }
static float calculateZPt(const std::pair< calib::CalibElectron *, calib::CalibElectron *> &aZCandidate)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T sqrt(T t)
Definition: SSEVec.h:23
static float calculateZMassWithCorrectedElectrons_withTK(const std::pair< calib::CalibElectron *, calib::CalibElectron *> &aZCandidate, float ele1EnergyCorrection, float ele2EnergyCorrection)
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 cosThetaElectrons_TK(const std::pair< calib::CalibElectron *, calib::CalibElectron *> &aZCandidate, float ele1EnergyCorrection, float ele2EnergyCorrection)
static float calculateZMass_withTK(const std::pair< calib::CalibElectron *, calib::CalibElectron *> &aZCandidate)
static float calculateZPhi(const std::pair< calib::CalibElectron *, calib::CalibElectron *> &aZCandidate)
static float calculateZEta(const std::pair< calib::CalibElectron *, calib::CalibElectron *> &aZCandidate)
static float calculateZMassWithCorrectedElectrons_noTK(const std::pair< calib::CalibElectron *, calib::CalibElectron *> &aZCandidate, float ele1EnergyCorrection, float ele2EnergyCorrection)
static float calculateZRapidity(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 cosThetaElectrons_SC(const std::pair< calib::CalibElectron *, calib::CalibElectron *> &aZCandidate, float ele1EnergyCorrection, float ele2EnergyCorrection)