CMS 3D CMS Logo

GEDGsfElectronValueMapProducer.cc
Go to the documentation of this file.
9 
11 public:
13  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
14  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
15 
16 private:
20 };
21 
23  : electronsToken_(consumes<reco::GsfElectronCollection>(cfg.getParameter<edm::InputTag>("gedGsfElectrons"))),
24  pfCandsToken_(consumes<reco::PFCandidateCollection>(cfg.getParameter<edm::InputTag>("egmPFCandidatesTag"))),
25  putToken_{produces<edm::ValueMap<reco::GsfElectronRef>>()} {}
26 
28  auto electrons = event.getHandle(electronsToken_);
29  auto pfCandidatesHandle = event.getHandle(pfCandsToken_);
30 
31  // ValueMap
34 
35  //Loop over the collection of PFFCandidates
36  std::vector<reco::GsfElectronRef> values;
37 
38  for (auto const& pfCandidate : *pfCandidatesHandle) {
40  // First check that the GsfTrack is non null
41  if (pfCandidate.gsfTrackRef().isNonnull()) {
42  // now look for the corresponding GsfElectron
43  const auto itcheck = std::find_if(electrons->begin(), electrons->end(), [&pfCandidate](const auto& ele) {
44  return (ele.gsfTrack() == pfCandidate.gsfTrackRef());
45  });
46  if (itcheck != electrons->end()) {
47  // Build the Ref from the handle and the index
48  myRef = reco::GsfElectronRef(electrons, itcheck - electrons->begin());
49  }
50  }
51  values.push_back(myRef);
52  }
53  valMapFiller.insert(pfCandidatesHandle, values.begin(), values.end());
54 
55  valMapFiller.fill();
56  event.emplace(putToken_, valMap);
57 }
58 
61  desc.add<edm::InputTag>("gedGsfElectrons", {"gedGsfElectronsTmp"});
62  desc.add<edm::InputTag>("egmPFCandidatesTag", {"particleFlowEGamma"});
63  descriptions.add("gedGsfElectronValueMapsTmp", desc);
64 }
65 
std::vector< l1t::PFCandidate > PFCandidateCollection
Definition: PFCandidate.h:82
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:53
const edm::EDPutTokenT< edm::ValueMap< reco::GsfElectronRef > > putToken_
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
GEDGsfElectronValueMapProducer(const edm::ParameterSet &)
const edm::EDGetTokenT< reco::GsfElectronCollection > electronsToken_
const edm::EDGetTokenT< reco::PFCandidateCollection > pfCandsToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::Ref< GsfElectronCollection > GsfElectronRef
reference to an object in a collection of GsfElectron objects
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
HLT enums.
Definition: event.py:1