CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/DataFormats/PatCandidates/src/MET.cc

Go to the documentation of this file.
00001 //
00002 // $Id: MET.cc,v 1.15 2010/08/09 14:39:02 mbluj Exp $
00003 //
00004 
00005 #include "DataFormats/PatCandidates/interface/MET.h"
00006 
00007 
00008 using namespace pat;
00009 
00010 
00012 MET::MET(): uncorInfo_(0) {
00013 }
00014 
00015 
00017 MET::MET(const reco::MET & aMET) : PATObject<reco::MET>(aMET), uncorInfo_(0) {
00018     const reco::CaloMET * calo = dynamic_cast<const reco::CaloMET *>(&aMET);
00019     if (calo != 0) caloMET_.push_back(calo->getSpecific());
00020     const reco::PFMET * pf = dynamic_cast<const reco::PFMET *>(&aMET);
00021     if (pf != 0) pfMET_.push_back(pf->getSpecific());
00022 }
00023 
00024 
00026 MET::MET(const edm::RefToBase<reco::MET> & aMETRef) : PATObject<reco::MET>(aMETRef), uncorInfo_(0) {
00027     const reco::CaloMET * calo = dynamic_cast<const reco::CaloMET *>(aMETRef.get());
00028     if (calo != 0) caloMET_.push_back(calo->getSpecific());
00029     const reco::PFMET * pf = dynamic_cast<const reco::PFMET *>(aMETRef.get());
00030     if (pf != 0) pfMET_.push_back(pf->getSpecific());
00031 }
00032 
00034 MET::MET(const edm::Ptr<reco::MET> & aMETRef) : PATObject<reco::MET>(aMETRef), uncorInfo_(0) {
00035     const reco::CaloMET * calo = dynamic_cast<const reco::CaloMET *>(aMETRef.get());
00036     if (calo != 0) caloMET_.push_back(calo->getSpecific());
00037     const reco::PFMET * pf = dynamic_cast<const reco::PFMET *>(aMETRef.get());
00038     if (pf != 0) pfMET_.push_back(pf->getSpecific());
00039 }
00040 
00041 
00043 MET::~MET() {
00044 }
00045 
00046 
00048 const reco::GenMET * MET::genMET() const {
00049   return (genMET_.size() > 0 ? &genMET_.front() : 0 );
00050 }
00051 
00053 void MET::setGenMET(const reco::GenMET & gm) {
00054   genMET_.clear();
00055   genMET_.push_back(gm);
00056 }
00057 
00059 unsigned int MET::nCorrections() const { checkUncor_(); return nCorrections_; }
00060 
00061 float MET::corEx(UncorrectionType ix) const { 
00062   if (ix == uncorrNONE) return 0;
00063   checkUncor_(); return uncorInfo_[ix].corEx; 
00064 }
00065 float MET::corEy(UncorrectionType ix) const { 
00066   if (ix == uncorrNONE) return 0;
00067   checkUncor_(); return uncorInfo_[ix].corEy; 
00068 }
00069 float MET::corSumEt(UncorrectionType ix) const { 
00070   if (ix == uncorrNONE) return 0;
00071   checkUncor_(); return uncorInfo_[ix].corSumEt; 
00072 }
00073 float MET::uncorrectedPt(UncorrectionType ix) const { 
00074   if (ix == uncorrNONE) return pt();
00075   checkUncor_(); return uncorInfo_[ix].pt; 
00076 }
00077 float MET::uncorrectedPhi(UncorrectionType ix) const { 
00078   if (ix == uncorrNONE) return phi();
00079   checkUncor_(); return uncorInfo_[ix].phi; 
00080 }
00081 
00082 
00084 void MET::checkUncor_() const {
00085   if (uncorInfo_.size() == uncorrMAXN && oldPt_ == pt() ) return;
00086 
00087   oldPt_ = pt();
00088   std::vector<CorrMETData> corrs(mEtCorr());
00089   nCorrections_ = corrs.size();
00090 
00091   uncorInfo_.resize(uncorrMAXN);
00092   UncorrectionType ix;
00093 
00096   ix = uncorrALL;
00097   uncorInfo_[ix] = UncorInfo();
00098   for (unsigned int iC=0; iC < nCorrections_; ++iC){
00099     uncorInfo_[ix].corEx +=    corrs[iC].mex;
00100     uncorInfo_[ix].corEy +=    corrs[iC].mey;
00101     uncorInfo_[ix].corSumEt += corrs[iC].sumet;
00102   }
00103   setPtPhi_(uncorInfo_[ix]);
00104 
00106   ix = uncorrJES;
00107   uncorInfo_[ix] = UncorInfo();
00108   if (nCorrections_ >=1 ){
00109     unsigned int iC = 0;
00110     uncorInfo_[ix].corEx +=    corrs[iC].mex;
00111     uncorInfo_[ix].corEy +=    corrs[iC].mey;
00112     uncorInfo_[ix].corSumEt += corrs[iC].sumet;
00113   }
00114   setPtPhi_(uncorInfo_[ix]);
00115 
00117   ix = uncorrMUON;
00118   uncorInfo_[ix] = UncorInfo();
00119   if (nCorrections_ >=2 ){
00120     unsigned int iC = 1;
00121     uncorInfo_[ix].corEx +=    corrs[iC].mex;
00122     uncorInfo_[ix].corEy +=    corrs[iC].mey;
00123     uncorInfo_[ix].corSumEt += corrs[iC].sumet;
00124   }
00125   setPtPhi_(uncorInfo_[ix]);
00126 
00128   ix = uncorrTAU;
00129   uncorInfo_[ix] = UncorInfo();
00130   if (nCorrections_ >=3 ){
00131     unsigned int iC = 1;
00132     uncorInfo_[ix].corEx +=    corrs[iC].mex;
00133     uncorInfo_[ix].corEy +=    corrs[iC].mey;
00134     uncorInfo_[ix].corSumEt += corrs[iC].sumet;
00135   }
00136   setPtPhi_(uncorInfo_[ix]);
00137 
00138 }
00139 
00140 void MET::setPtPhi_(UncorInfo& uci) const {
00141   float lpx = px() - uci.corEx;
00142   float lpy = py() - uci.corEy;
00143   uci.pt = sqrt(lpx*lpx + lpy*lpy);
00144   uci.phi = atan2(lpy, lpx);  
00145 }