CMS 3D CMS Logo

CorrectedPATMETProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 
3 //____________________________________________________________________________||
12 
16 
19 
20 #include <vector>
21 
22 //____________________________________________________________________________||
24 {
25 
26 public:
27 
29  : corrector(),
30  token_(consumes<METCollection>(cfg.getParameter<edm::InputTag>("src")))
31  {
32  std::vector<edm::InputTag> corrInputTags = cfg.getParameter<std::vector<edm::InputTag> >("srcCorrections");
33  std::vector<edm::EDGetTokenT<CorrMETData> > corrTokens;
34  for (std::vector<edm::InputTag>::const_iterator inputTag = corrInputTags.begin(); inputTag != corrInputTags.end(); ++inputTag) {
35  corrTokens.push_back(consumes<CorrMETData>(*inputTag));
36  }
37 
38  corrector.setCorTokens(corrTokens);
39 
40  produces<METCollection>("");
41  }
42 
43  ~CorrectedPATMETProducer() override { }
44 
45 private:
46 
48 
49  typedef std::vector<pat::MET> METCollection;
50 
52 
53  void produce(edm::Event& evt, const edm::EventSetup& es) override
54  {
55  edm::Handle<METCollection> srcMETCollection;
56  evt.getByToken(token_, srcMETCollection);
57 
58  const pat::MET& srcMET = (*srcMETCollection)[0];
59 
60  //dispatching to be sure we retrieve all the informations
61  reco::MET corrMET = corrector.getCorrectedMET(srcMET, evt, es);
62  pat::MET outMET(corrMET, srcMET);
63 
64  auto product = std::make_unique<METCollection>();
65 
67  if( !(cov(0,0)==0 && cov(0,1)==0 && cov(1,0)==0 && cov(1,1)==0) ) {
68  outMET.setSignificanceMatrix(cov);
69  double metSig=metsig::METSignificance::getSignificance(cov, outMET);
70  outMET.setMETSignificance(metSig);
71  }
72 
73  product->push_back(outMET);
74  evt.put(std::move(product));
75  }
76 
77 };
78 
79 //____________________________________________________________________________||
80 
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
static double getSignificance(const reco::METCovMatrix &cov, const reco::MET &met)
CorrectedPATMETProducer(const edm::ParameterSet &cfg)
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
void setSignificanceMatrix(const reco::METCovMatrix &matrix)
Definition: MET.cc:157
edm::EDGetTokenT< METCollection > token_
ROOT::Math::SMatrix< double, 2 > METCovMatrix
Definition: MET.h:40
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
AddCorrectionsToGenericMET corrector
Definition: MET.h:42
std::vector< pat::MET > METCollection
Collection of MET.
void produce(edm::Event &evt, const edm::EventSetup &es) override
HLT enums.
void setMETSignificance(const double &metSig)
Definition: MET.cc:124
def move(src, dest)
Definition: eostools.py:511
void setCorTokens(std::vector< edm::EDGetTokenT< CorrMETData > > const &corrTokens)
reco::METCovMatrix getSignificanceMatrix(void) const
Definition: MET.cc:139