File indexing completed on 2024-04-06 12:24:49
0001 #include "CommonTools/CandAlgos/interface/ModifyObjectValueBase.h"
0002 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0003 #include "FWCore/Framework/interface/stream/EDProducer.h"
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/Framework/interface/EventSetup.h"
0006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0007 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0010
0011 class LowPtGsfElectronFinalizer : public edm::stream::EDProducer<> {
0012 public:
0013 explicit LowPtGsfElectronFinalizer(const edm::ParameterSet&);
0014
0015 void produce(edm::Event&, const edm::EventSetup&) override;
0016
0017 static void fillDescriptions(edm::ConfigurationDescriptions&);
0018
0019 private:
0020 const edm::EDGetTokenT<reco::GsfElectronCollection> previousGsfElectrons_;
0021 std::unique_ptr<ModifyObjectValueBase> regression_;
0022
0023 const edm::EDPutTokenT<reco::GsfElectronCollection> putToken_;
0024 };
0025
0026 using edm::InputTag;
0027 using reco::GsfElectronCollection;
0028
0029 LowPtGsfElectronFinalizer::LowPtGsfElectronFinalizer(const edm::ParameterSet& cfg)
0030 : previousGsfElectrons_{consumes<GsfElectronCollection>(cfg.getParameter<InputTag>("previousGsfElectronsTag"))},
0031 putToken_{produces<GsfElectronCollection>()} {
0032 auto const& iconf = cfg.getParameterSet("regressionConfig");
0033 auto const& mname = iconf.getParameter<std::string>("modifierName");
0034 auto cc = consumesCollector();
0035 regression_ = ModifyObjectValueFactory::get()->create(mname, iconf, cc);
0036 }
0037
0038 void LowPtGsfElectronFinalizer::produce(edm::Event& event, const edm::EventSetup& setup) {
0039
0040 regression_->setEvent(event);
0041 regression_->setEventContent(setup);
0042
0043
0044 reco::GsfElectronCollection outputElectrons;
0045 for (auto const& electron : event.get(previousGsfElectrons_)) {
0046 outputElectrons.emplace_back(electron);
0047 auto& newElectron = outputElectrons.back();
0048 regression_->modifyObject(newElectron);
0049 }
0050
0051
0052 event.emplace(putToken_, std::move(outputElectrons));
0053 }
0054
0055 void LowPtGsfElectronFinalizer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0056 edm::ParameterSetDescription desc;
0057 desc.add<edm::InputTag>("previousGsfElectronsTag", {});
0058 edm::ParameterSetDescription psd;
0059 psd.setUnknown();
0060 desc.add<edm::ParameterSetDescription>("regressionConfig", psd);
0061 descriptions.addWithDefaultLabel(desc);
0062 }
0063
0064 #include "FWCore/Framework/interface/MakerMacros.h"
0065 DEFINE_FWK_MODULE(LowPtGsfElectronFinalizer);