CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AddCorrectionsToGenericMET.cc
Go to the documentation of this file.
2 
3 void
5  corrTokens_=corrTokens;
6 }
7 
10  const CorrMETData& correction) {
11  double px = met.px() + correction.mex;
12  double py = met.py() + correction.mey;
13  double pt = sqrt(px*px + py*py);
14  return reco::Candidate::LorentzVector(px, py, 0., pt);
15 }
16 
17 
20  CorrMETData sumCor;
22  for (std::vector<edm::EDGetTokenT<CorrMETData> >::const_iterator corrToken = corrTokens_.begin(); corrToken != corrTokens_.end(); ++corrToken) {
23  evt.getByToken(*corrToken, corr);
24  sumCor += (*corr);
25  }
26 
27  return sumCor;
28 }
29 
30 
33 
34  CorrMETData corr = getCorrection(srcMET, evt, es);
35  reco::MET outMET(srcMET.sumEt()+corr.sumet, constructP4From(srcMET, corr), srcMET.vertex() );
36 
37  return outMET;
38 }
39 
40 //specific flavors ================================
43 
44  CorrMETData corr = getCorrection(srcMET, evt, es);
45  reco::PFMET outMET(srcMET.getSpecific(),srcMET.sumEt()+corr.sumet, constructP4From(srcMET, corr), srcMET.vertex() );
46 
47  return outMET;
48 }
49 
50 
53 
54  CorrMETData corr = getCorrection(srcMET, evt, es);
55  reco::CaloMET outMET(srcMET.getSpecific(),srcMET.sumEt()+corr.sumet, constructP4From(srcMET, corr), srcMET.vertex() );
56 
57  return outMET;
58 }
SpecificPFMETData getSpecific() const
Definition: PFMET.h:72
reco::MET getCorrectedMET(const reco::MET &srcMET, edm::Event &evt, const edm::EventSetup &es)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
virtual const Point & vertex() const
vertex position (overwritten by PF...)
reco::PFMET getCorrectedPFMET(const reco::PFMET &srcMET, edm::Event &evt, const edm::EventSetup &es)
SpecificCaloMETData getSpecific() const
Definition: CaloMET.h:79
double sumEt() const
Definition: MET.h:56
std::vector< edm::EDGetTokenT< CorrMETData > > corrTokens_
Definition: MET.h:42
T sqrt(T t)
Definition: SSEVec.h:48
reco::Candidate::LorentzVector constructP4From(const reco::MET &met, const CorrMETData &correction)
double sumet
Definition: CorrMETData.h:20
JetCorrectorParameters corr
Definition: classes.h:5
virtual double px() const
x coordinate of momentum vector
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
a MET correction term
Definition: CorrMETData.h:14
CorrMETData getCorrection(const reco::MET &srcMET, edm::Event &evt, const edm::EventSetup &es)
reco::CaloMET getCorrectedCaloMET(const reco::CaloMET &srcMET, edm::Event &evt, const edm::EventSetup &es)
double mey
Definition: CorrMETData.h:18
double mex
Definition: CorrMETData.h:17
virtual double py() const
y coordinate of momentum vector
void setCorTokens(std::vector< edm::EDGetTokenT< CorrMETData > > const &corrTokens)