File indexing completed on 2024-04-06 12:24:04
0001
0002
0003
0004 #include "FWCore/Framework/interface/Frameworkfwd.h"
0005 #include "FWCore/Framework/interface/stream/EDProducer.h"
0006 #include "FWCore/Framework/interface/Event.h"
0007 #include "FWCore/Framework/interface/EventSetup.h"
0008 #include "FWCore/Framework/interface/MakerMacros.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0010 #include "DataFormats/Candidate/interface/Candidate.h"
0011 #include "FWCore/Utilities/interface/InputTag.h"
0012
0013 #include "DataFormats/METReco/interface/CorrMETData.h"
0014
0015 #include <vector>
0016
0017
0018 class CorrMETDataExtractor : public edm::stream::EDProducer<> {
0019 public:
0020 explicit CorrMETDataExtractor(const edm::ParameterSet& cfg) {
0021 std::vector<edm::InputTag> corrInputTags = cfg.getParameter<std::vector<edm::InputTag> >("corrections");
0022 std::vector<edm::EDGetTokenT<CorrMETData> > corrTokens;
0023 for (std::vector<edm::InputTag>::const_iterator inputTag = corrInputTags.begin(); inputTag != corrInputTags.end();
0024 ++inputTag) {
0025 corrTokens_.push_back(consumes<CorrMETData>(*inputTag));
0026 }
0027
0028 produces<float>("corX");
0029 produces<float>("corY");
0030 produces<float>("corSumEt");
0031 }
0032
0033 ~CorrMETDataExtractor() override {}
0034
0035 private:
0036 std::vector<edm::EDGetTokenT<CorrMETData> > corrTokens_;
0037
0038 void produce(edm::Event& evt, const edm::EventSetup& es) override {
0039 CorrMETData sumCor;
0040 edm::Handle<CorrMETData> corr;
0041 for (std::vector<edm::EDGetTokenT<CorrMETData> >::const_iterator corrToken = corrTokens_.begin();
0042 corrToken != corrTokens_.end();
0043 ++corrToken) {
0044 evt.getByToken(*corrToken, corr);
0045 sumCor += (*corr);
0046 }
0047
0048 float cX = (float)sumCor.mex;
0049 float cY = (float)sumCor.mey;
0050 float cSEt = (float)sumCor.sumet;
0051
0052 std::unique_ptr<float> corX(new float(0));
0053 std::unique_ptr<float> corY(new float(0));
0054 std::unique_ptr<float> corSumEt(new float(0));
0055
0056 *corX = cX;
0057 *corY = cY;
0058 *corSumEt = cSEt;
0059
0060 evt.put(std::move(corX), "corX");
0061 evt.put(std::move(corY), "corY");
0062 evt.put(std::move(corSumEt), "corSumEt");
0063 }
0064 };
0065
0066
0067
0068 DEFINE_FWK_MODULE(CorrMETDataExtractor);