Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "JetMETCorrections/Type1MET/plugins/PFchsMETcorrInputProducer.h"
0002 
0003 #include "DataFormats/Common/interface/View.h"
0004 #include "DataFormats/VertexReco/interface/Vertex.h"
0005 
0006 PFchsMETcorrInputProducer::PFchsMETcorrInputProducer(const edm::ParameterSet& cfg)
0007     : moduleLabel_(cfg.getParameter<std::string>("@module_label")),
0008       token_(consumes<reco::VertexCollection>(cfg.getParameter<edm::InputTag>("src"))),
0009       goodVtxNdof_(cfg.getParameter<unsigned int>("goodVtxNdof")),
0010       goodVtxZ_(cfg.getParameter<double>("goodVtxZ")) {
0011   produces<CorrMETData>("type0");
0012 }
0013 
0014 PFchsMETcorrInputProducer::~PFchsMETcorrInputProducer() {}
0015 
0016 void PFchsMETcorrInputProducer::produce(edm::Event& evt, const edm::EventSetup& es) {
0017   edm::Handle<reco::VertexCollection> recVtxs;
0018   evt.getByToken(token_, recVtxs);
0019 
0020   std::unique_ptr<CorrMETData> chsSum(new CorrMETData());
0021 
0022   for (unsigned i = 1; i < recVtxs->size(); ++i) {
0023     const reco::Vertex& v = recVtxs->at(i);
0024     if (v.isFake())
0025       continue;
0026     if (v.ndof() < goodVtxNdof_)
0027       continue;
0028     if (fabs(v.z()) > goodVtxZ_)
0029       continue;
0030 
0031     for (reco::Vertex::trackRef_iterator track = v.tracks_begin(); track != v.tracks_end(); ++track) {
0032       if ((*track)->charge() != 0) {
0033         chsSum->mex += (*track)->px();
0034         chsSum->mey += (*track)->py();
0035         chsSum->sumet += (*track)->pt();
0036       }
0037     }
0038   }
0039 
0040   evt.put(std::move(chsSum), "type0");
0041 }
0042 
0043 #include "FWCore/Framework/interface/MakerMacros.h"
0044 
0045 DEFINE_FWK_MODULE(PFchsMETcorrInputProducer);