CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CorrectedPFMETProducer2.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 
3 //____________________________________________________________________________||
12 
14 
16 
17 #include <vector>
18 
19 //____________________________________________________________________________||
21 {
22 
23 public:
24 
26  : src_(cfg.getParameter<edm::InputTag>("src")),
27  srcCorrections_(cfg.getParameter<std::vector<edm::InputTag> >("srcCorrections"))
28  {
29  produces<METCollection>("");
30  }
31 
33 
34 private:
35 
36  typedef std::vector<reco::PFMET> METCollection;
37 
39  std::vector<edm::InputTag> srcCorrections_;
40 
41  void produce(edm::Event& evt, const edm::EventSetup& es)
42  {
43  edm::Handle<METCollection> srcMETCollection;
44  evt.getByLabel(src_, srcMETCollection);
45 
46  const reco::PFMET& srcMET = (*srcMETCollection)[0];
47 
48  CorrMETData correction = readAndSumCorrections(evt, es);
49 
50  reco::PFMET outMET = applyCorrection(srcMET, correction);
51 
52  std::auto_ptr<METCollection> product(new METCollection);
53  product->push_back(outMET);
54  evt.put(product);
55  }
56 
58  {
60 
62  for (std::vector<edm::InputTag>::const_iterator inputTag = srcCorrections_.begin(); inputTag != srcCorrections_.end(); ++inputTag)
63  {
64  evt.getByLabel(*inputTag, corr);
65  ret += (*corr);
66  }
67  return ret;
68  }
69 
70  reco::PFMET applyCorrection(const reco::PFMET& srcMET, const CorrMETData& correction)
71  {
72  return reco::PFMET(srcMET.getSpecific(), srcMET.sumEt() + correction.sumet, constructP4From(srcMET, correction), srcMET.vertex());
73  }
74 
76  {
77  double px = met.px() + correction.mex;
78  double py = met.py() + correction.mey;
79  double pt = sqrt(px*px + py*py);
80  return reco::Candidate::LorentzVector(px, py, 0., pt);
81  }
82 };
83 
84 //____________________________________________________________________________||
85 
87 
SpecificPFMETData getSpecific() const
Definition: PFMET.h:72
std::vector< reco::PFMET > METCollection
virtual const Point & vertex() const
vertex position
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
reco::PFMET applyCorrection(const reco::PFMET &srcMET, const CorrMETData &correction)
CorrectedPFMETProducer2(const edm::ParameterSet &cfg)
double sumEt() const
Definition: MET.h:48
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
T sqrt(T t)
Definition: SSEVec.h:46
Collection of MET.
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
double sumet
Definition: CorrMETData.h:20
JetCorrectorParameters corr
Definition: classes.h:9
virtual double px() const
x coordinate of momentum vector
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:38
a MET correction term
Definition: CorrMETData.h:14
CorrMETData readAndSumCorrections(edm::Event &evt, const edm::EventSetup &es)
void produce(edm::Event &evt, const edm::EventSetup &es)
double mey
Definition: CorrMETData.h:18
double mex
Definition: CorrMETData.h:17
reco::Candidate::LorentzVector constructP4From(const reco::PFMET &met, const CorrMETData &correction)
virtual double py() const
y coordinate of momentum vector
std::vector< edm::InputTag > srcCorrections_