CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MET.cc
Go to the documentation of this file.
1 //
2 // $Id: MET.cc,v 1.16 2013/02/19 16:18:45 vadler Exp $
3 //
4 
6 
7 
8 using namespace pat;
9 
10 
12 MET::MET(): uncorInfo_(0) {
13 }
14 
15 
17 MET::MET(const reco::MET & aMET) : PATObject<reco::MET>(aMET), uncorInfo_(0) {
18  const reco::CaloMET * calo = dynamic_cast<const reco::CaloMET *>(&aMET);
19  if (calo != 0) caloMET_.push_back(calo->getSpecific());
20  const reco::PFMET * pf = dynamic_cast<const reco::PFMET *>(&aMET);
21  if (pf != 0) pfMET_.push_back(pf->getSpecific());
22 }
23 
24 
26 MET::MET(const edm::RefToBase<reco::MET> & aMETRef) : PATObject<reco::MET>(aMETRef), uncorInfo_(0) {
27  const reco::CaloMET * calo = dynamic_cast<const reco::CaloMET *>(aMETRef.get());
28  if (calo != 0) caloMET_.push_back(calo->getSpecific());
29  const reco::PFMET * pf = dynamic_cast<const reco::PFMET *>(aMETRef.get());
30  if (pf != 0) pfMET_.push_back(pf->getSpecific());
31 }
32 
34 MET::MET(const edm::Ptr<reco::MET> & aMETRef) : PATObject<reco::MET>(aMETRef), uncorInfo_(0) {
35  const reco::CaloMET * calo = dynamic_cast<const reco::CaloMET *>(aMETRef.get());
36  if (calo != 0) caloMET_.push_back(calo->getSpecific());
37  const reco::PFMET * pf = dynamic_cast<const reco::PFMET *>(aMETRef.get());
38  if (pf != 0) pfMET_.push_back(pf->getSpecific());
39 }
40 
41 
44 }
45 
46 
48 const reco::GenMET * MET::genMET() const {
49  return (genMET_.size() > 0 ? &genMET_.front() : 0 );
50 }
51 
53 void MET::setGenMET(const reco::GenMET & gm) {
54  genMET_.clear();
55  genMET_.push_back(gm);
56 }
57 
59 unsigned int MET::nCorrections() const { checkUncor_(); return nCorrections_; }
60 
61 float MET::corEx(UncorrectionType ix) const {
62  if (ix == uncorrNONE) return 0;
63  checkUncor_(); return uncorInfo_[ix].corEx;
64 }
65 float MET::corEy(UncorrectionType ix) const {
66  if (ix == uncorrNONE) return 0;
67  checkUncor_(); return uncorInfo_[ix].corEy;
68 }
70  if (ix == uncorrNONE) return 0;
71  checkUncor_(); return uncorInfo_[ix].corSumEt;
72 }
74  if (ix == uncorrNONE) return pt();
75  checkUncor_(); return uncorInfo_[ix].pt;
76 }
78  if (ix == uncorrNONE) return phi();
79  checkUncor_(); return uncorInfo_[ix].phi;
80 }
81 
82 
84 void MET::checkUncor_() const {
85  if (uncorInfo_.size() == uncorrMAXN && oldPt_ == pt() ) return;
86 
87  oldPt_ = pt();
88  std::vector<CorrMETData> corrs(mEtCorr());
89  nCorrections_ = corrs.size();
90 
91  uncorInfo_.resize(uncorrMAXN);
93 
96  ix = uncorrALL;
97  uncorInfo_[ix] = UncorInfo();
98  for (unsigned int iC=0; iC < nCorrections_; ++iC){
99  uncorInfo_[ix].corEx += corrs[iC].mex;
100  uncorInfo_[ix].corEy += corrs[iC].mey;
101  uncorInfo_[ix].corSumEt += corrs[iC].sumet;
102  }
103  setPtPhi_(uncorInfo_[ix]);
104 
106  ix = uncorrJES;
107  uncorInfo_[ix] = UncorInfo();
108  if (nCorrections_ >=1 ){
109  unsigned int iC = 0;
110  uncorInfo_[ix].corEx += corrs[iC].mex;
111  uncorInfo_[ix].corEy += corrs[iC].mey;
112  uncorInfo_[ix].corSumEt += corrs[iC].sumet;
113  }
114  setPtPhi_(uncorInfo_[ix]);
115 
117  ix = uncorrMUON;
118  uncorInfo_[ix] = UncorInfo();
119  if (nCorrections_ >=2 ){
120  unsigned int iC = 1;
121  uncorInfo_[ix].corEx += corrs[iC].mex;
122  uncorInfo_[ix].corEy += corrs[iC].mey;
123  uncorInfo_[ix].corSumEt += corrs[iC].sumet;
124  }
125  setPtPhi_(uncorInfo_[ix]);
126 
128  ix = uncorrTAU;
129  uncorInfo_[ix] = UncorInfo();
130  if (nCorrections_ >=3 ){
131  unsigned int iC = 2;
132  uncorInfo_[ix].corEx += corrs[iC].mex;
133  uncorInfo_[ix].corEy += corrs[iC].mey;
134  uncorInfo_[ix].corSumEt += corrs[iC].sumet;
135  }
136  setPtPhi_(uncorInfo_[ix]);
137 
138 }
139 
140 void MET::setPtPhi_(UncorInfo& uci) const {
141  float lpx = px() - uci.corEx;
142  float lpy = py() - uci.corEy;
143  uci.pt = sqrt(lpx*lpx + lpy*lpy);
144  uci.phi = atan2(lpy, lpx);
145 }
Analysis-level MET class.
Definition: MET.h:42
float uncorrectedPhi(UncorrectionType ix=uncorrALL) const
Definition: MET.cc:77
MET()
default constructor
Definition: MET.cc:12
SpecificPFMETData getSpecific() const
Definition: PFMET.h:72
std::vector< UncorInfo > uncorInfo_
Definition: MET.h:173
uncorrect for JES only
Definition: MET.h:79
void checkUncor_() const
check and set transients
Definition: MET.cc:84
float corEx(UncorrectionType ix=uncorrALL) const
Definition: MET.cc:61
float oldPt_
Definition: MET.h:175
virtual double py() const GCC11_FINAL
y coordinate of momentum vector
unsigned int nCorrections_
Definition: MET.h:174
float corSumEt(UncorrectionType ix=uncorrALL) const
Definition: MET.cc:69
uncorrect to bare bones
Definition: MET.h:78
SpecificCaloMETData getSpecific() const
Definition: CaloMET.h:79
virtual float phi() const GCC11_FINAL
momentum azimuthal angle
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:143
Definition: MET.h:32
T sqrt(T t)
Definition: SSEVec.h:48
virtual double px() const GCC11_FINAL
x coordinate of momentum vector
std::vector< reco::GenMET > genMET_
Definition: MET.h:157
uncorrect for MUON only
Definition: MET.h:80
float uncorrectedPt(UncorrectionType ix=uncorrALL) const
Definition: MET.cc:73
std::vector< CorrMETData > mEtCorr() const
Definition: MET.h:63
float corEy(UncorrectionType ix=uncorrALL) const
Definition: MET.cc:65
const reco::GenMET * genMET() const
return the associated GenMET
Definition: MET.cc:48
void setPtPhi_(UncorInfo &uci) const
Definition: MET.cc:140
std::vector< SpecificCaloMETData > caloMET_
Definition: MET.h:159
Templated PAT object container.
Definition: PATObject.h:43
void setGenMET(const reco::GenMET &gm)
set the associated GenMET
Definition: MET.cc:53
std::vector< SpecificPFMETData > pfMET_
Definition: MET.h:161
virtual ~MET()
destructor
Definition: MET.cc:43
uncorrect for TAU only
Definition: MET.h:81
virtual float pt() const GCC11_FINAL
transverse momentum
UncorrectionType
Definition: MET.h:75
value_type const * get() const
Definition: RefToBase.h:212
do nothing
Definition: MET.h:77
unsigned int nCorrections() const
return uncorrrection related stuff
Definition: MET.cc:59