File indexing completed on 2024-04-06 12:24:50
0001 #include "DataFormats/Common/interface/Handle.h"
0002 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0003 #include "DataFormats/EgammaReco/interface/ElectronSeed.h"
0004 #include "DataFormats/EgammaReco/interface/ElectronSeedFwd.h"
0005 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
0006 #include "DataFormats/GsfTrackReco/interface/GsfTrackExtra.h"
0007 #include "DataFormats/ParticleFlowReco/interface/PreId.h"
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0011 #include "FWCore/Utilities/interface/InputTag.h"
0012 #include "DataFormats/Common/interface/ValueMap.h"
0013 #include "DataFormats/GsfTrackReco/interface/GsfTrackFwd.h"
0014 #include "DataFormats/ParticleFlowReco/interface/PreIdFwd.h"
0015 #include "FWCore/Framework/interface/ESHandle.h"
0016 #include "FWCore/Framework/interface/stream/EDProducer.h"
0017 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0018 #include "FWCore/Utilities/interface/transform.h"
0019
0020 #include <vector>
0021 #include <string>
0022
0023 class LowPtGsfElectronSeedValueMapsProducer : public edm::stream::EDProducer<> {
0024 public:
0025 explicit LowPtGsfElectronSeedValueMapsProducer(const edm::ParameterSet&);
0026
0027 void produce(edm::Event&, const edm::EventSetup&) override;
0028
0029 static void fillDescriptions(edm::ConfigurationDescriptions&);
0030
0031 private:
0032 edm::EDGetTokenT<reco::GsfTrackCollection> gsfTracks_;
0033 edm::EDGetTokenT<edm::ValueMap<reco::PreIdRef> > preIdsValueMap_;
0034 std::vector<std::string> names_;
0035 const bool rekey_;
0036 edm::EDGetTokenT<reco::GsfElectronCollection> gsfElectrons_;
0037 std::vector<edm::EDGetTokenT<edm::ValueMap<float> > > floatValueMaps_;
0038 };
0039
0040
0041
0042 LowPtGsfElectronSeedValueMapsProducer::LowPtGsfElectronSeedValueMapsProducer(const edm::ParameterSet& conf)
0043 : gsfTracks_(),
0044 preIdsValueMap_(),
0045 names_(),
0046 rekey_(conf.getParameter<bool>("rekey")),
0047 gsfElectrons_(),
0048 floatValueMaps_() {
0049 if (rekey_) {
0050 gsfElectrons_ = consumes<reco::GsfElectronCollection>(conf.getParameter<edm::InputTag>("gsfElectrons"));
0051 std::vector<edm::InputTag> tags = conf.getParameter<std::vector<edm::InputTag> >("floatValueMaps");
0052 for (const auto& tag : tags) {
0053 floatValueMaps_ = edm::vector_transform(
0054 tags, [this](edm::InputTag const& tag) { return consumes<edm::ValueMap<float> >(tag); });
0055 names_.push_back(tag.instance());
0056 produces<edm::ValueMap<float> >(tag.instance());
0057 }
0058 } else {
0059 gsfTracks_ = consumes<reco::GsfTrackCollection>(conf.getParameter<edm::InputTag>("gsfTracks"));
0060 preIdsValueMap_ = consumes<edm::ValueMap<reco::PreIdRef> >(conf.getParameter<edm::InputTag>("preIdsValueMap"));
0061 names_ = conf.getParameter<std::vector<std::string> >("ModelNames");
0062 for (const auto& name : names_) {
0063 produces<edm::ValueMap<float> >(name);
0064 }
0065 }
0066 }
0067
0068
0069
0070 void LowPtGsfElectronSeedValueMapsProducer::produce(edm::Event& event, const edm::EventSetup&) {
0071 if (rekey_ == false) {
0072
0073
0074
0075
0076 auto gsfTracks = event.getHandle(gsfTracks_);
0077
0078
0079 auto preIdsValueMap = event.getHandle(preIdsValueMap_);
0080
0081
0082 std::vector<std::vector<float> > output;
0083 for (unsigned int iname = 0; iname < names_.size(); ++iname) {
0084 output.push_back(std::vector<float>(gsfTracks->size(), -999.));
0085 }
0086 auto const& gsfTracksV = *gsfTracks;
0087 for (unsigned int igsf = 0; igsf < gsfTracksV.size(); igsf++) {
0088 const reco::GsfTrack& gsf = gsfTracksV[igsf];
0089 if (gsf.extra().isNonnull() && gsf.extra()->seedRef().isNonnull()) {
0090 reco::ElectronSeedRef seed = gsf.extra()->seedRef().castTo<reco::ElectronSeedRef>();
0091 if (seed.isNonnull() && seed->ctfTrack().isNonnull()) {
0092 const reco::PreIdRef preid = (*preIdsValueMap)[seed->ctfTrack()];
0093 if (preid.isNonnull()) {
0094 for (unsigned int iname = 0; iname < names_.size(); ++iname) {
0095 output[iname][igsf] = preid->mva(iname);
0096 }
0097 }
0098 }
0099 }
0100 }
0101
0102
0103 for (unsigned int iname = 0; iname < names_.size(); ++iname) {
0104 auto ptr = std::make_unique<edm::ValueMap<float> >(edm::ValueMap<float>());
0105 edm::ValueMap<float>::Filler filler(*ptr);
0106 filler.insert(gsfTracks, output[iname].begin(), output[iname].end());
0107 filler.fill();
0108 event.put(std::move(ptr), names_[iname]);
0109 }
0110
0111 } else {
0112
0113
0114
0115
0116 auto gsfElectrons = event.getHandle(gsfElectrons_);
0117
0118
0119 for (unsigned int idx = 0; idx < names_.size(); ++idx) {
0120
0121 auto const& floatValueMap = event.get(floatValueMaps_[idx]);
0122
0123
0124 std::vector<float> output(gsfElectrons->size(), -99.);
0125 auto const& gsfElectronsV = *gsfElectrons;
0126 for (unsigned int iele = 0; iele < gsfElectronsV.size(); iele++) {
0127 const reco::GsfElectron& ele = gsfElectronsV[iele];
0128 reco::GsfTrackRef gsf = ele.gsfTrack();
0129 output[iele] = floatValueMap[gsf];
0130 }
0131
0132 auto ptr = std::make_unique<edm::ValueMap<float> >(edm::ValueMap<float>());
0133 edm::ValueMap<float>::Filler filler(*ptr);
0134 filler.insert(gsfElectrons, output.begin(), output.end());
0135 filler.fill();
0136 event.put(std::move(ptr), names_[idx]);
0137 }
0138 }
0139 }
0140
0141
0142
0143 void LowPtGsfElectronSeedValueMapsProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0144 edm::ParameterSetDescription desc;
0145 desc.add<edm::InputTag>("gsfTracks", edm::InputTag("lowPtGsfEleGsfTracks"));
0146 desc.add<edm::InputTag>("preIdsValueMap", edm::InputTag("lowPtGsfElectronSeeds"));
0147 desc.add<std::vector<std::string> >("ModelNames", {"unbiased", "ptbiased"});
0148 desc.add<bool>("rekey", false);
0149 desc.add<edm::InputTag>("gsfElectrons", edm::InputTag());
0150 desc.add<std::vector<edm::InputTag> >("floatValueMaps", std::vector<edm::InputTag>());
0151 descriptions.add("lowPtGsfElectronSeedValueMaps", desc);
0152 }
0153
0154
0155
0156 #include "FWCore/Framework/interface/MakerMacros.h"
0157 DEFINE_FWK_MODULE(LowPtGsfElectronSeedValueMapsProducer);