File indexing completed on 2024-04-06 12:19:26
0001 #include "JetMETCorrections/Type1MET/plugins/SysShiftMETcorrInputProducer.h"
0002
0003 #include "FWCore/Utilities/interface/Exception.h"
0004
0005 #include "DataFormats/METReco/interface/CorrMETData.h"
0006 #include "DataFormats/VertexReco/interface/Vertex.h"
0007 #include "DataFormats/Common/interface/View.h"
0008
0009 #include <TString.h>
0010
0011 SysShiftMETcorrInputProducer::SysShiftMETcorrInputProducer(const edm::ParameterSet& cfg)
0012 : moduleLabel_(cfg.getParameter<std::string>("@module_label")), useNvtx(false), corrPx_(nullptr), corrPy_(nullptr) {
0013 token_ = consumes<edm::View<reco::MET> >(cfg.getParameter<edm::InputTag>("src"));
0014
0015 edm::ParameterSet cfgCorrParameter = cfg.getParameter<edm::ParameterSet>("parameter");
0016 TString corrPxFormula = cfgCorrParameter.getParameter<std::string>("px");
0017 TString corrPyFormula = cfgCorrParameter.getParameter<std::string>("py").data();
0018 if (corrPxFormula.Contains("Nvtx") || corrPyFormula.Contains("Nvtx")) {
0019 useNvtx = true, verticesToken_ = consumes<reco::VertexCollection>(cfg.getParameter<edm::InputTag>("srcVertices"));
0020 }
0021
0022 corrPxFormula.ReplaceAll("sumEt", "x");
0023 corrPxFormula.ReplaceAll("Nvtx", "y");
0024 std::string corrPxName = std::string(moduleLabel_).append("_corrPx");
0025 corrPx_ = new TFormula(corrPxName.data(), corrPxFormula.Data());
0026
0027 corrPyFormula.ReplaceAll("sumEt", "x");
0028 corrPyFormula.ReplaceAll("Nvtx", "y");
0029 std::string corrPyName = std::string(moduleLabel_).append("_corrPy");
0030 corrPy_ = new TFormula(corrPyName.data(), corrPyFormula.Data());
0031
0032 produces<CorrMETData>();
0033 }
0034
0035 SysShiftMETcorrInputProducer::~SysShiftMETcorrInputProducer() {
0036
0037 }
0038
0039 void SysShiftMETcorrInputProducer::produce(edm::Event& evt, const edm::EventSetup& es) {
0040
0041
0042 typedef edm::View<reco::MET> METView;
0043 edm::Handle<METView> met;
0044 evt.getByToken(token_, met);
0045 if (met->size() != 1)
0046 throw cms::Exception("SysShiftMETcorrInputProducer::produce") << "Failed to find unique MET object !!\n";
0047
0048 double sumEt = met->front().sumEt();
0049
0050
0051 size_t Nvtx = 0;
0052 if (useNvtx) {
0053 edm::Handle<reco::VertexCollection> vertices;
0054 evt.getByToken(verticesToken_, vertices);
0055 Nvtx = vertices->size();
0056 }
0057
0058
0059 std::unique_ptr<CorrMETData> metCorr(new CorrMETData());
0060 metCorr->mex = -corrPx_->Eval(sumEt, Nvtx);
0061 metCorr->mey = -corrPy_->Eval(sumEt, Nvtx);
0062
0063
0064 evt.put(std::move(metCorr));
0065 }
0066
0067 #include "FWCore/Framework/interface/MakerMacros.h"
0068
0069 DEFINE_FWK_MODULE(SysShiftMETcorrInputProducer);