CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MuonMETAlgo.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: METAlgorithms
4 // Class: MuonMETAlgo
5 //
6 // Original Authors: Michael Schmitt, Richard Cavanaugh The University of Florida
7 // Created: August 30, 2007
8 //
9 
10 //____________________________________________________________________________||
11 #include <math.h>
12 #include <vector>
18 
19 #include "TMath.h"
20 
21 
22 using namespace std;
23 using namespace reco;
24 
27 
28 //____________________________________________________________________________||
30  double fSumEt,
31  const std::vector<CorrMETData>& fCorrections,
32  const MET::LorentzVector& fP4) {
33  return CaloMET (fMet.getSpecific (), fSumEt, fCorrections, fP4, fMet.vertex ());
34 }
35 
36 
37 //____________________________________________________________________________||
38 MET MuonMETAlgo::makeMET (const MET& fMet,
39  double fSumEt,
40  const std::vector<CorrMETData>& fCorrections,
41  const MET::LorentzVector& fP4) {
42  return MET (fSumEt, fCorrections, fP4, fMet.vertex ());
43 }
44 
45 
46 
47 //____________________________________________________________________________||
48 template <class T> void MuonMETAlgo::MuonMETAlgo_run(const edm::View<reco::Muon>& inputMuons,
49  const edm::ValueMap<reco::MuonMETCorrectionData>& vm_muCorrData,
50  const edm::View<T>& v_uncorMET,
51  std::vector<T>* v_corMET)
52 {
53  T uncorMETObj = v_uncorMET.front();
54 
55  double corMETX = uncorMETObj.px();
56  double corMETY = uncorMETObj.py();
57 
58  double sumMuPx = 0.;
59  double sumMuPy = 0.;
60  double sumMuPt = 0.;
61  double sumMuDepEx = 0.;
62  double sumMuDepEy = 0.;
63  double sumMuDepEt = 0.;
64 
65  for(unsigned int i = 0; i < inputMuons.size(); ++i)
66  {
67  reco::MuonMETCorrectionData muCorrData = vm_muCorrData[inputMuons.refAt(i)];
68 
69  if(muCorrData.type() == 0) continue;
70 
71  float deltax = muCorrData.corrX();
72  float deltay = muCorrData.corrY();
73 
74  const reco::Muon *mu = &inputMuons[i];
75  LorentzVector mup4 = mu->p4();
76 
77  sumMuPx += mup4.px();
78  sumMuPy += mup4.py();
79  sumMuPt += mup4.pt();
80  sumMuDepEx += deltax;
81  sumMuDepEy += deltay;
82  sumMuDepEt += sqrt(deltax*deltax + deltay*deltay);
83  corMETX = corMETX - mup4.px() + deltax;
84  corMETY = corMETY - mup4.py() + deltay;
85 
86  }
87 
89  delta.mex = sumMuDepEx - sumMuPx;
90  delta.mey = sumMuDepEy - sumMuPy;
91  delta.sumet = sumMuPt - sumMuDepEt;
92  MET::LorentzVector correctedMET4vector(corMETX, corMETY, 0., sqrt(corMETX*corMETX + corMETY*corMETY));
93  std::vector<CorrMETData> corrections = uncorMETObj.mEtCorr();
94  corrections.push_back(delta);
95 
96  T result = makeMET(uncorMETObj, uncorMETObj.sumEt() + delta.sumet, corrections, correctedMET4vector);
97  v_corMET->push_back(result);
98 }
99 
100 //____________________________________________________________________________||
103  bool useTrackAssociatorPositions,
104  bool useRecHits,
105  bool useHO,
106  double towerEtThreshold,
107  double& deltax, double& deltay,
108  double Bfield) {
109 
110 
111  bool useAverage = false;
112  //decide whether or not we want to correct on average based
113  //on isolation information from the muon
114  double sumPt = inputMuon->isIsolationValid()? inputMuon->isolationR03().sumPt : 0.0;
115  double sumEtEcal = inputMuon->isIsolationValid() ? inputMuon->isolationR03().emEt : 0.0;
116  double sumEtHcal = inputMuon->isIsolationValid() ? inputMuon->isolationR03().hadEt : 0.0;
117 
118  if(sumPt > 3 || sumEtEcal + sumEtHcal > 5) useAverage = true;
119 
120  //get the energy using TrackAssociator if
121  //the muon turns out to be isolated
122  MuonMETInfo muMETInfo;
123  muMETInfo.useAverage = useAverage;
124  muMETInfo.useTkAssociatorPositions = useTrackAssociatorPositions;
125  muMETInfo.useHO = useHO;
126 
127 
128  TrackRef mu_track;
129  if(inputMuon->isGlobalMuon()) {
130  mu_track = inputMuon->globalTrack();
131  } else if(inputMuon->isTrackerMuon()) {
132  mu_track = inputMuon->innerTrack();
133  } else
134  mu_track = inputMuon->outerTrack();
135 
136  if(useTrackAssociatorPositions) {
137  muMETInfo.ecalPos = info.trkGlobPosAtEcal;
138  muMETInfo.hcalPos = info.trkGlobPosAtHcal;
139  muMETInfo.hoPos = info.trkGlobPosAtHO;
140  }
141 
142  if(!useAverage) {
143 
144  if(useRecHits) {
145  muMETInfo.ecalE = inputMuon->calEnergy().emS9;
146  muMETInfo.hcalE = inputMuon->calEnergy().hadS9;
147  if(useHO) //muMETInfo.hoE is 0 by default
148  muMETInfo.hoE = inputMuon->calEnergy().hoS9;
149  } else {// use Towers (this is the default)
150  //only include towers whose Et > 0.5 since
151  //by default the MET only includes towers with Et > 0.5
152  std::vector<const CaloTower*> towers = info.crossedTowers;
153  for(vector<const CaloTower*>::const_iterator it = towers.begin();
154  it != towers.end(); it++) {
155  if((*it)->et() < towerEtThreshold) continue;
156  muMETInfo.ecalE += (*it)->emEnergy();
157  muMETInfo.hcalE += (*it)->hadEnergy();
158  if(useHO)
159  muMETInfo.hoE +=(*it)->outerEnergy();
160  }
161  }//use Towers
162  }
163 
164 
165 
166  //This needs to be fixed!!!!!
167  //The tracker has better resolution for pt < 200 GeV
169  if(inputMuon->isGlobalMuon()) {
170  if(inputMuon->globalTrack()->pt() < 200) {
171  mup4 = LorentzVector(inputMuon->innerTrack()->px(), inputMuon->innerTrack()->py(),
172  inputMuon->innerTrack()->pz(), inputMuon->innerTrack()->p());
173  } else {
174  mup4 = LorentzVector(inputMuon->globalTrack()->px(), inputMuon->globalTrack()->py(),
175  inputMuon->globalTrack()->pz(), inputMuon->globalTrack()->p());
176  }
177  } else if(inputMuon->isTrackerMuon()) {
178  mup4 = LorentzVector(inputMuon->innerTrack()->px(), inputMuon->innerTrack()->py(),
179  inputMuon->innerTrack()->pz(), inputMuon->innerTrack()->p());
180  } else
181  mup4 = LorentzVector(inputMuon->outerTrack()->px(), inputMuon->outerTrack()->py(),
182  inputMuon->outerTrack()->pz(), inputMuon->outerTrack()->p());
183 
184 
185  //call function that does the work
186  correctMETforMuon(deltax, deltay, Bfield, inputMuon->charge(),
187  mup4, inputMuon->vertex(),
188  muMETInfo);
189 }
190 
191 
192 //____________________________________________________________________________||
193 void MuonMETAlgo::correctMETforMuon(double& deltax, double& deltay, double bfield, int muonCharge,
194  const math::XYZTLorentzVector& muonP4,const math::XYZPoint& muonVertex,
195  MuonMETInfo& muonMETInfo) {
196 
197  double mu_p = muonP4.P();
198  double mu_pt = muonP4.Pt();
199  double mu_phi = muonP4.Phi();
200  double mu_eta = muonP4.Eta();
201  double mu_vz = muonVertex.z()/100.;
202  double mu_pz = muonP4.Pz();
203 
204  double ecalPhi, ecalTheta;
205  double hcalPhi, hcalTheta;
206  double hoPhi, hoTheta;
207 
208 
209  //should always be false for FWLite
210  //unless you want to supply co-ordinates at
211  //the calorimeter sub-detectors yourself
212  if(muonMETInfo.useTkAssociatorPositions) {
213  ecalPhi = muonMETInfo.ecalPos.Phi();
214  ecalTheta = muonMETInfo.ecalPos.Theta();
215  hcalPhi = muonMETInfo.hcalPos.Phi();
216  hcalTheta = muonMETInfo.hcalPos.Theta();
217  hoPhi = muonMETInfo.hoPos.Phi();
218  hoTheta = muonMETInfo.hoPos.Theta();
219  } else {
220 
221  /*
222  use the analytical solution for the
223  intersection of a helix with a cylinder
224  to find the positions of the muon
225  at the various calo surfaces
226  */
227 
228  //radii of subdetectors in meters
229  double rEcal = 1.290;
230  double rHcal = 1.9;
231  double rHo = 3.82;
232  if(abs(mu_eta) > 0.3) rHo = 4.07;
233  //distance from the center of detector to face of Ecal
234  double zFaceEcal = 3.209;
235  if(mu_eta < 0 ) zFaceEcal = -1*zFaceEcal;
236  //distance from the center of detector to face of Hcal
237  double zFaceHcal = 3.88;
238  if(mu_eta < 0 ) zFaceHcal = -1*zFaceHcal;
239 
240  //now we have to get Phi
241  //bending radius of the muon (units are meters)
242  double bendr = mu_pt*1000/(300*bfield);
243 
244  double tb_ecal = TMath::ACos(1-rEcal*rEcal/(2*bendr*bendr)); //helix time interval parameter
245  double tb_hcal = TMath::ACos(1-rHcal*rHcal/(2*bendr*bendr)); //helix time interval parameter
246  double tb_ho = TMath::ACos(1-rHo*rHo/(2*bendr*bendr)); //helix time interval parameter
247  double xEcal,yEcal,zEcal;
248  double xHcal,yHcal,zHcal;
249  double xHo, yHo,zHo;
250  //Ecal
251  //in the barrel and if not a looper
252  if(fabs(mu_pz*bendr*tb_ecal/mu_pt+mu_vz) < fabs(zFaceEcal) && rEcal < 2*bendr) {
253  xEcal = bendr*(TMath::Sin(tb_ecal+mu_phi)-TMath::Sin(mu_phi));
254  yEcal = bendr*(-TMath::Cos(tb_ecal+mu_phi)+TMath::Cos(mu_phi));
255  zEcal = bendr*tb_ecal*mu_pz/mu_pt + mu_vz;
256  } else { //endcap
257  if(mu_pz > 0) {
258  double te_ecal = (fabs(zFaceEcal) - mu_vz)*mu_pt/(bendr*mu_pz);
259  xEcal = bendr*(TMath::Sin(te_ecal+mu_phi) - TMath::Sin(mu_phi));
260  yEcal = bendr*(-TMath::Cos(te_ecal+mu_phi) + TMath::Cos(mu_phi));
261  zEcal = fabs(zFaceEcal);
262  } else {
263  double te_ecal = -(fabs(zFaceEcal) + mu_vz)*mu_pt/(bendr*mu_pz);
264  xEcal = bendr*(TMath::Sin(te_ecal+mu_phi) - TMath::Sin(mu_phi));
265  yEcal = bendr*(-TMath::Cos(te_ecal+mu_phi) + TMath::Cos(mu_phi));
266  zEcal = -fabs(zFaceEcal);
267  }
268  }
269 
270  //Hcal
271  if(fabs(mu_pz*bendr*tb_hcal/mu_pt+mu_vz) < fabs(zFaceHcal) && rEcal < 2*bendr) { //in the barrel
272  xHcal = bendr*(TMath::Sin(tb_hcal+mu_phi)-TMath::Sin(mu_phi));
273  yHcal = bendr*(-TMath::Cos(tb_hcal+mu_phi)+TMath::Cos(mu_phi));
274  zHcal = bendr*tb_hcal*mu_pz/mu_pt + mu_vz;
275  } else { //endcap
276  if(mu_pz > 0) {
277  double te_hcal = (fabs(zFaceHcal) - mu_vz)*mu_pt/(bendr*mu_pz);
278  xHcal = bendr*(TMath::Sin(te_hcal+mu_phi) - TMath::Sin(mu_phi));
279  yHcal = bendr*(-TMath::Cos(te_hcal+mu_phi) + TMath::Cos(mu_phi));
280  zHcal = fabs(zFaceHcal);
281  } else {
282  double te_hcal = -(fabs(zFaceHcal) + mu_vz)*mu_pt/(bendr*mu_pz);
283  xHcal = bendr*(TMath::Sin(te_hcal+mu_phi) - TMath::Sin(mu_phi));
284  yHcal = bendr*(-TMath::Cos(te_hcal+mu_phi) + TMath::Cos(mu_phi));
285  zHcal = -fabs(zFaceHcal);
286  }
287  }
288 
289  //Ho - just the barrel
290  xHo = bendr*(TMath::Sin(tb_ho+mu_phi)-TMath::Sin(mu_phi));
291  yHo = bendr*(-TMath::Cos(tb_ho+mu_phi)+TMath::Cos(mu_phi));
292  zHo = bendr*tb_ho*mu_pz/mu_pt + mu_vz;
293 
294  ecalTheta = TMath::ACos(zEcal/sqrt(pow(xEcal,2) + pow(yEcal,2)+pow(zEcal,2)));
295  ecalPhi = atan2(yEcal,xEcal);
296  hcalTheta = TMath::ACos(zHcal/sqrt(pow(xHcal,2) + pow(yHcal,2)+pow(zHcal,2)));
297  hcalPhi = atan2(yHcal,xHcal);
298  hoTheta = TMath::ACos(zHo/sqrt(pow(xHo,2) + pow(yHo,2)+pow(zHo,2)));
299  hoPhi = atan2(yHo,xHo);
300 
301  //2d radius in x-y plane
302  double r2dEcal = sqrt(pow(xEcal,2)+pow(yEcal,2));
303  double r2dHcal = sqrt(pow(xHcal,2)+pow(yHcal,2));
304  double r2dHo = sqrt(pow(xHo,2) +pow(yHo,2));
305 
306  /*
307  the above prescription is for right handed helicies only
308  Positively charged muons trace a left handed helix
309  so we correct for that
310  */
311  if(muonCharge > 0) {
312 
313  //Ecal
314  double dphi = mu_phi - ecalPhi;
315  if(fabs(dphi) > TMath::Pi())
316  dphi = 2*TMath::Pi() - fabs(dphi);
317  ecalPhi = mu_phi - fabs(dphi);
318  if(fabs(ecalPhi) > TMath::Pi()) {
319  double temp = 2*TMath::Pi() - fabs(ecalPhi);
320  ecalPhi = -1*temp*ecalPhi/fabs(ecalPhi);
321  }
322  xEcal = r2dEcal*TMath::Cos(ecalPhi);
323  yEcal = r2dEcal*TMath::Sin(ecalPhi);
324 
325  //Hcal
326  dphi = mu_phi - hcalPhi;
327  if(fabs(dphi) > TMath::Pi())
328  dphi = 2*TMath::Pi() - fabs(dphi);
329  hcalPhi = mu_phi - fabs(dphi);
330  if(fabs(hcalPhi) > TMath::Pi()) {
331  double temp = 2*TMath::Pi() - fabs(hcalPhi);
332  hcalPhi = -1*temp*hcalPhi/fabs(hcalPhi);
333  }
334  xHcal = r2dHcal*TMath::Cos(hcalPhi);
335  yHcal = r2dHcal*TMath::Sin(hcalPhi);
336 
337 
338  //Ho
339  dphi = mu_phi - hoPhi;
340  if(fabs(dphi) > TMath::Pi())
341  dphi = 2*TMath::Pi() - fabs(dphi);
342  hoPhi = mu_phi - fabs(dphi);
343  if(fabs(hoPhi) > TMath::Pi()) {
344  double temp = 2*TMath::Pi() - fabs(hoPhi);
345  hoPhi = -1*temp*hoPhi/fabs(hoPhi);
346  }
347  xHo = r2dHo*TMath::Cos(hoPhi);
348  yHo = r2dHo*TMath::Sin(hoPhi);
349 
350  }
351  }
352 
353  //for isolated muons
354  if(!muonMETInfo.useAverage) {
355 
356  double mu_Ex = muonMETInfo.ecalE*sin(ecalTheta)*cos(ecalPhi)
357  + muonMETInfo.hcalE*sin(hcalTheta)*cos(hcalPhi)
358  + muonMETInfo.hoE*sin(hoTheta)*cos(hoPhi);
359  double mu_Ey = muonMETInfo.ecalE*sin(ecalTheta)*sin(ecalPhi)
360  + muonMETInfo.hcalE*sin(hcalTheta)*sin(hcalPhi)
361  + muonMETInfo.hoE*sin(hoTheta)*sin(hoPhi);
362 
363  deltax += mu_Ex;
364  deltay += mu_Ey;
365 
366  } else { //non-isolated muons - derive the correction
367 
368  //dE/dx in matter for iron:
369  //-(11.4 + 0.96*fabs(log(p0*2.8)) + 0.033*p0*(1.0 - pow(p0, -0.33)) )*1e-3
370  //from http://cmslxr.fnal.gov/lxr/source/TrackPropagation/SteppingHelixPropagator/src/SteppingHelixPropagator.ccyes,
371  //line ~1100
372  //normalisation is at 50 GeV
373  double dEdx_normalization = -(11.4 + 0.96*fabs(log(50*2.8)) + 0.033*50*(1.0 - pow(50, -0.33)) )*1e-3;
374  double dEdx_numerator = -(11.4 + 0.96*fabs(log(mu_p*2.8)) + 0.033*mu_p*(1.0 - pow(mu_p, -0.33)) )*1e-3;
375 
376  double temp = 0.0;
377 
378  if(muonMETInfo.useHO) {
379  //for the Towers, with HO
380  if(fabs(mu_eta) < 0.2)
381  temp = 2.75*(1-0.00003*mu_p);
382  if(fabs(mu_eta) > 0.2 && fabs(mu_eta) < 1.0)
383  temp = (2.38+0.0144*fabs(mu_eta))*(1-0.0003*mu_p);
384  if(fabs(mu_eta) > 1.0 && fabs(mu_eta) < 1.3)
385  temp = 7.413-5.12*fabs(mu_eta);
386  if(fabs(mu_eta) > 1.3)
387  temp = 2.084-0.743*fabs(mu_eta);
388  } else {
389  if(fabs(mu_eta) < 1.0)
390  temp = 2.33*(1-0.0004*mu_p);
391  if(fabs(mu_eta) > 1.0 && fabs(mu_eta) < 1.3)
392  temp = (7.413-5.12*fabs(mu_eta))*(1-0.0003*mu_p);
393  if(fabs(mu_eta) > 1.3)
394  temp = 2.084-0.743*fabs(mu_eta);
395  }
396 
397  double dep = temp*dEdx_normalization/dEdx_numerator;
398  if(dep < 0.5) dep = 0;
399  //use the average phi of the 3 subdetectors
400  if(fabs(mu_eta) < 1.3) {
401  deltax += dep*cos((ecalPhi+hcalPhi+hoPhi)/3);
402  deltay += dep*sin((ecalPhi+hcalPhi+hoPhi)/3);
403  } else {
404  deltax += dep*cos( (ecalPhi+hcalPhi)/2);
405  deltay += dep*cos( (ecalPhi+hcalPhi)/2);
406  }
407  }
408 
409 
410 }
411 
412 //____________________________________________________________________________||
413 void MuonMETAlgo::run(const edm::View<reco::Muon>& inputMuons,
414  const edm::ValueMap<reco::MuonMETCorrectionData>& vm_muCorrData,
415  const edm::View<reco::MET>& uncorMET,
416  METCollection *corMET) {
417 
418  MuonMETAlgo_run(inputMuons, vm_muCorrData, uncorMET, corMET);
419 }
420 
421 //____________________________________________________________________________||
422 void MuonMETAlgo::run(const edm::View<reco::Muon>& inputMuons,
423  const edm::ValueMap<reco::MuonMETCorrectionData>& vm_muCorrData,
424  const edm::View<reco::CaloMET>& uncorMET,
425  CaloMETCollection *corMET) {
426 
427  MuonMETAlgo_run(inputMuons, vm_muCorrData, uncorMET, corMET);
428 
429 }
430 
431 //____________________________________________________________________________||
dbl * delta
Definition: mlp_gen.cc:36
const double Pi
float hadEt
hcal sum-Et
Definition: MuonIsolation.h:9
int i
Definition: DBlmapReader.cc:9
math::XYZPoint trkGlobPosAtHO
static const TGPicture * info(bool iBackgroundIsBlack)
static std::vector< std::string > checklist log
math::XYZPoint hcalPos
Definition: MuonMETInfo.h:12
float sumPt
sum-pt of tracks
Definition: MuonIsolation.h:7
std::vector< const CaloTower * > crossedTowers
virtual TrackRef innerTrack() const
Definition: Muon.h:48
virtual const Point & vertex() const
vertex position (overwritten by PF...)
bool isTrackerMuon() const
Definition: Muon.h:219
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
bool isGlobalMuon() const
Definition: Muon.h:218
float emS9
energy deposited in 3x3 ECAL crystal shape around central crystal
Definition: MuonEnergy.h:18
RefToBase< value_type > refAt(size_type i) const
float hcalE
Definition: MuonMETInfo.h:9
std::pair< double, double > Point
Definition: CaloEllipse.h:18
math::XYZPoint trkGlobPosAtHcal
Collection of Calo MET.
void MuonMETAlgo_run(const edm::View< reco::Muon > &inputMuons, const edm::ValueMap< reco::MuonMETCorrectionData > &vm_muCorrData, const edm::View< T > &v_uncorMET, std::vector< T > *v_corMET)
Definition: MuonMETAlgo.cc:48
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
void GetMuDepDeltas(const reco::Muon *inputMuon, TrackDetMatchInfo &info, bool useTrackAssociatorPositions, bool useRecHits, bool useHO, double towerEtThreshold, double &deltax, double &deltay, double Bfield)
Definition: MuonMETAlgo.cc:101
SpecificCaloMETData getSpecific() const
Definition: CaloMET.h:79
Definition: MET.h:39
T sqrt(T t)
Definition: SSEVec.h:48
float emEt
ecal sum-Et
Definition: MuonIsolation.h:8
tuple result
Definition: query.py:137
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
virtual int charge() const
electric charge
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float hoS9
energy deposited in 3x3 HO tower shape around central tower
Definition: MuonEnergy.h:34
const int mu
Definition: Constants.h:22
Collection of MET.
math::XYZPoint hoPos
Definition: MuonMETInfo.h:13
bool useTkAssociatorPositions
Definition: MuonMETInfo.h:16
float ecalE
Definition: MuonMETInfo.h:8
virtual TrackRef outerTrack() const
reference to Track reconstructed in the muon detector only
Definition: Muon.h:51
double sumet
Definition: CorrMETData.h:20
MuonEnergy calEnergy() const
get energy deposition information
Definition: Muon.h:111
bool useAverage
Definition: MuonMETInfo.h:14
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
size_type size() const
a MET correction term
Definition: CorrMETData.h:14
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< float > > LorentzVector
Definition: analysisEnums.h:9
static void correctMETforMuon(double &deltax, double &deltay, double bfield, int muonCharge, const math::XYZTLorentzVector &muonP4, const math::XYZPoint &muonVertex, MuonMETInfo &)
Definition: MuonMETAlgo.cc:193
math::XYZPoint ecalPos
Definition: MuonMETInfo.h:11
math::XYZPoint trkGlobPosAtEcal
Track position at different parts of the calorimeter.
virtual void run(const edm::View< reco::Muon > &inputMuons, const edm::ValueMap< reco::MuonMETCorrectionData > &vm_muCorrData, const edm::View< reco::MET > &uncorMET, reco::METCollection *corMET)
double mey
Definition: CorrMETData.h:18
double mex
Definition: CorrMETData.h:17
long double T
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
reco::CaloMET makeMET(const reco::CaloMET &fMet, double fSumEt, const std::vector< CorrMETData > &fCorrections, const reco::MET::LorentzVector &)
bool isIsolationValid() const
Definition: Muon.h:170
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
float hadS9
energy deposited in 3x3 HCAL tower shape around central tower
Definition: MuonEnergy.h:28
const MuonIsolation & isolationR03() const
Definition: Muon.h:158
math::PtEtaPhiELorentzVectorF LorentzVector
virtual TrackRef globalTrack() const
reference to Track reconstructed in both tracked and muon detector
Definition: Muon.h:54