CMS 3D CMS Logo

SysShiftMETcorrInputProducer.cc
Go to the documentation of this file.
2 
4 
8 
9 #include <TString.h>
10 
12  : moduleLabel_(cfg.getParameter<std::string>("@module_label")),
13  useNvtx(false),
14  corrPx_(nullptr),
15  corrPy_(nullptr)
16 {
17  token_ = consumes<edm::View<reco::MET> >(cfg.getParameter<edm::InputTag>("src"));
18 
19  edm::ParameterSet cfgCorrParameter = cfg.getParameter<edm::ParameterSet>("parameter");
20  TString corrPxFormula = cfgCorrParameter.getParameter<std::string>("px");
21  TString corrPyFormula = cfgCorrParameter.getParameter<std::string>("py").data();
22  if ( corrPxFormula.Contains("Nvtx") || corrPyFormula.Contains("Nvtx") )
23  {
24  useNvtx = true,
25  verticesToken_ = consumes<reco::VertexCollection>(cfg.getParameter<edm::InputTag>("srcVertices"));
26  }
27 
28  corrPxFormula.ReplaceAll("sumEt", "x");
29  corrPxFormula.ReplaceAll("Nvtx", "y");
30  std::string corrPxName = std::string(moduleLabel_).append("_corrPx");
31  corrPx_ = new TFormula(corrPxName.data(), corrPxFormula.Data());
32 
33  corrPyFormula.ReplaceAll("sumEt", "x");
34  corrPyFormula.ReplaceAll("Nvtx", "y");
35  std::string corrPyName = std::string(moduleLabel_).append("_corrPy");
36  corrPy_ = new TFormula(corrPyName.data(), corrPyFormula.Data());
37 
38  produces<CorrMETData>();
39 }
40 
42 {
43 // nothing to be done yet...
44 }
45 
47 {
48  //std::cout << "<SysShiftMETcorrInputProducer::produce>:" << std::endl;
49 
50  typedef edm::View<reco::MET> METView;
52  evt.getByToken(token_, met);
53  if ( met->size() != 1 )
54  throw cms::Exception("SysShiftMETcorrInputProducer::produce")
55  << "Failed to find unique MET object !!\n";
56 
57  double sumEt = met->front().sumEt();
58  //std::cout << " sumEt = " << sumEt << std::endl;
59 
60  size_t Nvtx = 0;
61  if ( useNvtx )
62  {
64  evt.getByToken(verticesToken_, vertices);
65  Nvtx = vertices->size();
66  }
67  //std::cout << " Nvtx = " << Nvtx << std::endl;
68 
69  std::unique_ptr<CorrMETData> metCorr(new CorrMETData());
70  metCorr->mex = -corrPx_->Eval(sumEt, Nvtx);
71  metCorr->mey = -corrPy_->Eval(sumEt, Nvtx);
72  //std::cout << "--> metCorr: Px = " << metCorr->mex << ", Py = " << metCorr->mey << std::endl;
73 
74  evt.put(std::move(metCorr));
75 }
76 
78 
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
SysShiftMETcorrInputProducer(const edm::ParameterSet &)
#define nullptr
edm::EDGetTokenT< reco::VertexCollection > verticesToken_
edm::EDGetTokenT< edm::View< reco::MET > > token_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void produce(edm::Event &, const edm::EventSetup &) override
met
===> hadronic RAZOR
a MET correction term
Definition: CorrMETData.h:14
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
def move(src, dest)
Definition: eostools.py:511