CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CorrectedMETProducerT.h
Go to the documentation of this file.
1 #ifndef JetMETCorrections_Type1MET_CorrectedMETProducer_h
2 #define JetMETCorrections_Type1MET_CorrectedMETProducer_h
3 
28 
32 
33 #include <vector>
34 
35 namespace CorrectedMETProducer_namespace
36 {
37  template <typename T>
38  reco::Candidate::LorentzVector correctedP4(const T& rawMEt, const CorrMETData& correction)
39  {
40  double correctedMEtPx = rawMEt.px() + correction.mex;
41  double correctedMEtPy = rawMEt.py() + correction.mey;
42  double correctedMEtPt = sqrt(correctedMEtPx*correctedMEtPx + correctedMEtPy*correctedMEtPy);
43  return reco::Candidate::LorentzVector(correctedMEtPx, correctedMEtPy, 0., correctedMEtPt);
44  }
45 
46  template <typename T>
47  double correctedSumEt(const T& rawMEt, const CorrMETData& correction)
48  {
49  return rawMEt.sumEt() + correction.sumet;
50  }
51 
52  template <typename T>
54  {
55  public:
56 
57  T operator()(const T&, const CorrMETData&) const
58  {
59  assert(0); // "place-holder" for template instantiations for concrete T types only, **not** to be called
60  }
61  };
62 }
63 
64 template<typename T>
66 {
67  typedef std::vector<T> METCollection;
68 
69  public:
70 
72  : moduleLabel_(cfg.getParameter<std::string>("@module_label")),
73  algorithm_(0)
74  {
75  src_ = cfg.getParameter<edm::InputTag>("src");
76 
78 
79  produces<METCollection>("");
80  }
82  {
83  delete algorithm_;
84  }
85 
86  private:
87 
88  void produce(edm::Event& evt, const edm::EventSetup& es)
89  {
90  std::auto_ptr<METCollection> correctedMEtCollection(new METCollection);
91 
92  edm::Handle<METCollection> rawMEtCollection;
93  evt.getByLabel(src_, rawMEtCollection);
94 
95  for ( typename METCollection::const_iterator rawMEt = rawMEtCollection->begin();
96  rawMEt != rawMEtCollection->end(); ++rawMEt ) {
97  CorrMETData correction = algorithm_->compMETCorrection(evt, es);
98 
100  T correctedMEt = correctedMET_factory(*rawMEt, correction);
101 
102  correctedMEtCollection->push_back(correctedMEt);
103  }
104 
105 //--- add collection of MET objects with Type 1 / Type 1 + 2 corrections applied to the event
106  evt.put(correctedMEtCollection);
107  }
108 
109  std::string moduleLabel_;
110 
111  edm::InputTag src_; // input collection
112 
113  METCorrectionAlgorithm* algorithm_; // algorithm for computing Type 1 / Type 1 + 2 MET corrections
114 };
115 
116 #endif
117 
118 
119 
T getParameter(std::string const &) const
CorrMETData compMETCorrection(edm::Event &, const edm::EventSetup &)
CorrectedMETProducerT(const edm::ParameterSet &cfg)
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.
reco::Candidate::LorentzVector correctedP4(const T &rawMEt, const CorrMETData &correction)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
double sumet
Definition: CorrMETData.h:20
T operator()(const T &, const CorrMETData &) const
void produce(edm::Event &evt, const edm::EventSetup &es)
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:38
std::vector< T > METCollection
a MET correction term
Definition: CorrMETData.h:14
double correctedSumEt(const T &rawMEt, const CorrMETData &correction)
double mey
Definition: CorrMETData.h:18
double mex
Definition: CorrMETData.h:17
long double T
METCorrectionAlgorithm * algorithm_