CMS 3D CMS Logo

CorrectedMETProducerT.h
Go to the documentation of this file.
1 #ifndef JetMETCorrections_Type1MET_CorrectedMETProducer_h
2 #define JetMETCorrections_Type1MET_CorrectedMETProducer_h
3 
29 
33 
35 
36 #include <vector>
37 
39 {
40  template <typename T>
41  reco::Candidate::LorentzVector correctedP4(const T& rawMEt, const CorrMETData& correction)
42  {
43  double correctedMEtPx = rawMEt.px() + correction.mex;
44  double correctedMEtPy = rawMEt.py() + correction.mey;
45  double correctedMEtPt = sqrt(correctedMEtPx*correctedMEtPx + correctedMEtPy*correctedMEtPy);
46  return reco::Candidate::LorentzVector(correctedMEtPx, correctedMEtPy, 0., correctedMEtPt);
47  }
48 
49  template <typename T>
50  double correctedSumEt(const T& rawMEt, const CorrMETData& correction)
51  {
52  return rawMEt.sumEt() + correction.sumet;
53  }
54 
55  template <typename T>
57  {
58  public:
59 
60  T operator()(const T&, const CorrMETData&) const
61  {
62  assert(0); // "place-holder" for template instantiations for concrete T types only, **not** to be called
63  }
64  };
65 }
66 
67 template<typename T>
69 {
70  typedef std::vector<T> METCollection;
71 
72  public:
73 
75  : moduleLabel_(cfg.getParameter<std::string>("@module_label")),
76  algorithm_(0)
77  {
78  token_ = consumes<METCollection>(cfg.getParameter<edm::InputTag>("src"));
79 
80  algorithm_ = new METCorrectionAlgorithm(cfg, consumesCollector());
81 
82  produces<METCollection>("");
83  }
85  {
86  delete algorithm_;
87  }
88 
89  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
91  desc.add<edm::InputTag>("src",edm::InputTag("corrPfMetType1", "type1"));
92  descriptions.add(defaultModuleLabel<CorrectedMETProducerT<T> >(),desc);
93  }
94 
95  private:
96 
97  void produce(edm::Event& evt, const edm::EventSetup& es)
98  {
99  std::unique_ptr<METCollection> correctedMEtCollection(new METCollection);
100 
101  edm::Handle<METCollection> rawMEtCollection;
102  evt.getByToken(token_, rawMEtCollection);
103 
104  for ( typename METCollection::const_iterator rawMEt = rawMEtCollection->begin();
105  rawMEt != rawMEtCollection->end(); ++rawMEt ) {
106  CorrMETData correction = algorithm_->compMETCorrection(evt, es);
107 
108  static const CorrectedMETProducer_namespace::CorrectedMETFactoryT<T> correctedMET_factory {};
109  T correctedMEt = correctedMET_factory(*rawMEt, correction);
110 
111  correctedMEtCollection->push_back(correctedMEt);
112  }
113 
114 //--- add collection of MET objects with Type 1 / Type 1 + 2 corrections applied to the event
115  evt.put(std::move(correctedMEtCollection));
116  }
117 
119 
121 
122  METCorrectionAlgorithm* algorithm_; // algorithm for computing Type 1 / Type 1 + 2 MET corrections
123 };
124 
125 #endif
126 
127 
128 
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
std::string defaultModuleLabel()
edm::EDGetTokenT< METCollection > token_
CorrectedMETProducerT(const edm::ParameterSet &cfg)
T sqrt(T t)
Definition: SSEVec.h:18
ParameterDescriptionBase * add(U const &iLabel, T const &value)
reco::Candidate::LorentzVector correctedP4(const T &rawMEt, const CorrMETData &correction)
double sumet
Definition: CorrMETData.h:20
T operator()(const T &, const CorrMETData &) const
void produce(edm::Event &evt, const edm::EventSetup &es)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
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_
def move(src, dest)
Definition: eostools.py:511
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)