CMS 3D CMS Logo

CorrectedPatMETProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 
3 //____________________________________________________________________________||
12 
16 
18 
19 #include <vector>
20 
21 //____________________________________________________________________________||
23 {
24 
25 public:
26 
28  : corrector()
29  //token_(consumes<METCollection>(cfg.getParameter<edm::InputTag>("src")))
30  {
31  isMiniAod = (cfg.exists("isMiniAod") ) ? cfg.getParameter<bool>("isMiniAod"): true;
32  if(isMiniAod)
33  {
34  patToken_=consumes<patMETCollection>(cfg.getParameter<edm::InputTag>("src"));
35  }else{
36  pfToken_=consumes<pfMETCollection>(cfg.getParameter<edm::InputTag>("src"));
37  }
38 
39 
40 
41  std::vector<edm::InputTag> corrInputTags = cfg.getParameter<std::vector<edm::InputTag> >("srcCorrections");
42  std::vector<edm::EDGetTokenT<CorrMETData> > corrTokens;
43  for (std::vector<edm::InputTag>::const_iterator inputTag = corrInputTags.begin(); inputTag != corrInputTags.end(); ++inputTag) {
44  corrTokens.push_back(consumes<CorrMETData>(*inputTag));
45  }
46 
47  corrector.setCorTokens(corrTokens);
48 
49  produces<patMETCollection>("");
50  }
51 
52  ~CorrectedPatMETProducer() override { }
53 
54 private:
55 
56  bool isMiniAod;
57 
59 
60  typedef std::vector<pat::MET> patMETCollection;
61  typedef std::vector<reco::PFMET> pfMETCollection;
62 
65 
66 
67  void produce(edm::Event& evt, const edm::EventSetup& es) override
68  {
69  edm::Handle<patMETCollection> srcPatMETCollection;
70  edm::Handle<pfMETCollection> srcPfMETCollection;
71  if(isMiniAod)
72  {
73  evt.getByToken(patToken_, srcPatMETCollection);
74  }else{
75  evt.getByToken(pfToken_, srcPfMETCollection);
76  }
77 
78 
79 
80  if(isMiniAod){
81  //std::unique_ptr<patMETCollection> product(new patMETCollection);
82  std::unique_ptr<patMETCollection> product(new patMETCollection);
83  const reco::MET& srcMET = (*srcPatMETCollection)[0];
84  pat::MET outMEtPat = corrector.getCorrectedMET(srcMET, evt, es);
85  product->push_back(outMEtPat);
86  evt.put(std::move(product));
87  }else{
88  std::unique_ptr<pfMETCollection> product(new pfMETCollection);
89  const reco::PFMET& srcMET = (*srcPfMETCollection)[0];
90  reco::PFMET outPfMEtReco = corrector.getCorrectedPFMET(srcMET, evt, es);
91  product->push_back(outPfMEtReco);
92  evt.put(std::move(product));
93  }
94 
95  }
96 
97 };
98 
99 //____________________________________________________________________________||
100 
Analysis-level MET class.
Definition: MET.h:40
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
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:517
AddCorrectionsToGenericMET corrector
bool exists(std::string const &parameterName) const
checks if a parameter exists
reco::PFMET getCorrectedPFMET(const reco::PFMET &srcMET, edm::Event &evt, const edm::EventSetup &es)
std::vector< reco::PFMET > pfMETCollection
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Definition: MET.h:42
void produce(edm::Event &evt, const edm::EventSetup &es) override
CorrectedPatMETProducer(const edm::ParameterSet &cfg)
edm::EDGetTokenT< patMETCollection > patToken_
edm::EDGetTokenT< pfMETCollection > pfToken_
std::vector< pat::MET > patMETCollection
def move(src, dest)
Definition: eostools.py:511
void setCorTokens(std::vector< edm::EDGetTokenT< CorrMETData > > const &corrTokens)