Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:19:27

0001 #include "JetMETCorrections/Type1MET/interface/AddCorrectionsToGenericMET.h"
0002 
0003 void AddCorrectionsToGenericMET::setCorTokens(std::vector<edm::EDGetTokenT<CorrMETData> > const& corrTokens) {
0004   corrTokens_ = corrTokens;
0005 }
0006 
0007 reco::Candidate::LorentzVector AddCorrectionsToGenericMET::constructP4From(const reco::MET& met,
0008                                                                            const CorrMETData& correction) {
0009   double px = met.px() + correction.mex;
0010   double py = met.py() + correction.mey;
0011   double pt = sqrt(px * px + py * py);
0012   return reco::Candidate::LorentzVector(px, py, 0., pt);
0013 }
0014 
0015 CorrMETData AddCorrectionsToGenericMET::getCorrection(const reco::MET& srcMET, edm::Event& evt) {
0016   CorrMETData sumCor;
0017   edm::Handle<CorrMETData> corr;
0018   for (std::vector<edm::EDGetTokenT<CorrMETData> >::const_iterator corrToken = corrTokens_.begin();
0019        corrToken != corrTokens_.end();
0020        ++corrToken) {
0021     evt.getByToken(*corrToken, corr);
0022     sumCor += (*corr);
0023   }
0024 
0025   return sumCor;
0026 }
0027 
0028 reco::MET AddCorrectionsToGenericMET::getCorrectedMET(const reco::MET& srcMET, edm::Event& evt) {
0029   CorrMETData corr = getCorrection(srcMET, evt);
0030   reco::MET outMET(srcMET.sumEt() + corr.sumet, constructP4From(srcMET, corr), srcMET.vertex(), srcMET.isWeighted());
0031 
0032   return outMET;
0033 }
0034 
0035 //specific flavors ================================
0036 reco::PFMET AddCorrectionsToGenericMET::getCorrectedPFMET(const reco::PFMET& srcMET, edm::Event& evt) {
0037   CorrMETData corr = getCorrection(srcMET, evt);
0038   reco::PFMET outMET(srcMET.getSpecific(),
0039                      srcMET.sumEt() + corr.sumet,
0040                      constructP4From(srcMET, corr),
0041                      srcMET.vertex(),
0042                      srcMET.isWeighted());
0043 
0044   return outMET;
0045 }
0046 
0047 reco::CaloMET AddCorrectionsToGenericMET::getCorrectedCaloMET(const reco::CaloMET& srcMET, edm::Event& evt) {
0048   CorrMETData corr = getCorrection(srcMET, evt);
0049   reco::CaloMET outMET(
0050       srcMET.getSpecific(), srcMET.sumEt() + corr.sumet, constructP4From(srcMET, corr), srcMET.vertex());
0051 
0052   return outMET;
0053 }