File indexing completed on 2024-04-06 12:23:55
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0013 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0015 #include "FWCore/Framework/interface/global/EDProducer.h"
0016 #include "FWCore/Framework/interface/Event.h"
0017 #include "FWCore/Framework/interface/ConsumesCollector.h"
0018 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0019 #include "FWCore/Utilities/interface/InputTag.h"
0020 #include "DataFormats/PatCandidates/interface/MET.h"
0021 #include "DataFormats/METReco/interface/MET.h"
0022
0023 #include <memory>
0024
0025 namespace pat {
0026
0027 class RecoMETExtractor : public edm::global::EDProducer<> {
0028 public:
0029 explicit RecoMETExtractor(const edm::ParameterSet& iConfig);
0030 ~RecoMETExtractor() override;
0031
0032 void produce(edm::StreamID streamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
0033
0034 private:
0035 edm::EDGetTokenT<std::vector<pat::MET> > metSrcToken_;
0036
0037 pat::MET::METCorrectionLevel corLevel_;
0038 };
0039
0040 }
0041
0042 using namespace pat;
0043
0044 RecoMETExtractor::RecoMETExtractor(const edm::ParameterSet& iConfig) {
0045 edm::InputTag metIT = iConfig.getParameter<edm::InputTag>("metSource");
0046 metSrcToken_ = consumes<pat::METCollection>(metIT);
0047
0048 std::string corLevel = iConfig.getParameter<std::string>("correctionLevel");
0049
0050
0051 if (corLevel == "raw") {
0052 corLevel_ = pat::MET::Raw;
0053 } else if (corLevel == "type1") {
0054 corLevel_ = pat::MET::Type1;
0055 } else if (corLevel == "type01") {
0056 corLevel_ = pat::MET::Type01;
0057 } else if (corLevel == "typeXY") {
0058 corLevel_ = pat::MET::TypeXY;
0059 } else if (corLevel == "type1XY") {
0060 corLevel_ = pat::MET::Type1XY;
0061 } else if (corLevel == "type01XY") {
0062 corLevel_ = pat::MET::Type01XY;
0063 } else if (corLevel == "type1Smear") {
0064 corLevel_ = pat::MET::Type1Smear;
0065 } else if (corLevel == "type01Smear") {
0066 corLevel_ = pat::MET::Type01Smear;
0067 } else if (corLevel == "type1SmearXY") {
0068 corLevel_ = pat::MET::Type1SmearXY;
0069 } else if (corLevel == "type01SmearXY") {
0070 corLevel_ = pat::MET::Type01SmearXY;
0071 } else if (corLevel == "rawCalo") {
0072 corLevel_ = pat::MET::RawCalo;
0073 } else if (corLevel == "rawDeepResponseTune") {
0074 corLevel_ = pat::MET::RawDeepResponseTune;
0075 } else if (corLevel == "rawDeepResolutionTune") {
0076 corLevel_ = pat::MET::RawDeepResolutionTune;
0077 } else {
0078
0079 }
0080
0081
0082 produces<std::vector<reco::MET> >();
0083 }
0084
0085 RecoMETExtractor::~RecoMETExtractor() {}
0086
0087 void RecoMETExtractor::produce(edm::StreamID streamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0088 edm::Handle<std::vector<pat::MET> > src;
0089 iEvent.getByToken(metSrcToken_, src);
0090 if (src->empty())
0091 edm::LogError("RecoMETExtractor::produce") << "input reco MET collection is empty";
0092
0093 std::vector<reco::MET>* metCol = new std::vector<reco::MET>();
0094
0095 reco::MET met(src->front().corSumEt(corLevel_), src->front().corP4(corLevel_), src->front().vertex());
0096
0097 metCol->push_back(met);
0098
0099 std::unique_ptr<std::vector<reco::MET> > recoMETs(metCol);
0100 iEvent.put(std::move(recoMETs));
0101 }
0102
0103 #include "FWCore/Framework/interface/MakerMacros.h"
0104 DEFINE_FWK_MODULE(RecoMETExtractor);