CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ElectronMatchedCandidateProducer.cc
Go to the documentation of this file.
4 
5 #include "DataFormats/Math/interface/deltaR.h" // reco::deltaR
6 
8  const edm::InputTag allelectrons("gsfElectrons");
9  electronCollectionToken_ = consumes<edm::View<reco::GsfElectron>>(
10  params.getUntrackedParameter<edm::InputTag>("ReferenceElectronCollection", allelectrons));
11  scCollectionToken_ = consumes<edm::View<reco::Candidate>>(params.getParameter<edm::InputTag>("src"));
12 
13  delRMatchingCut_ = params.getUntrackedParameter<double>("deltaR", 0.30);
14 
15  produces<edm::PtrVector<reco::Candidate>>();
16  produces<edm::RefToBaseVector<reco::Candidate>>();
17 }
18 
20 
21 //
22 // member functions
23 //
24 
25 // ------------ method called to produce the data ------------
26 
28  // Create the output collection
29  auto outColRef = std::make_unique<edm::RefToBaseVector<reco::Candidate>>();
30  auto outColPtr = std::make_unique<edm::PtrVector<reco::Candidate>>();
31 
32  // Read electrons
34  event.getByToken(electronCollectionToken_, electrons);
35 
36  //Read candidates
38  event.getByToken(scCollectionToken_, recoCandColl);
39 
40  unsigned int counter = 0;
41 
42  // Loop over candidates
43  for (edm::View<reco::Candidate>::const_iterator scIt = recoCandColl->begin(); scIt != recoCandColl->end();
44  ++scIt, ++counter) {
45  // Now loop over electrons
46  for (edm::View<reco::GsfElectron>::const_iterator elec = electrons->begin(); elec != electrons->end(); ++elec) {
47  reco::SuperClusterRef eSC = elec->superCluster();
48 
49  double dRval = reco::deltaR((float)eSC->eta(), (float)eSC->phi(), scIt->eta(), scIt->phi());
50 
51  if (dRval < delRMatchingCut_) {
52  //outCol->push_back( *scIt );
53  outColRef->push_back(recoCandColl->refAt(counter));
54  outColPtr->push_back(recoCandColl->ptrAt(counter));
55  } // end if loop
56  } // end electron loop
57 
58  } // end candidate loop
59 
60  event.put(std::move(outColRef));
61  event.put(std::move(outColPtr));
62 }
63 
64 // ------ method called once each job just before starting event loop ---
65 
67 
69 
70 //define this as a plug-in
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
void produce(edm::Event &, const edm::EventSetup &) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
edm::EDGetTokenT< edm::View< reco::Candidate > > scCollectionToken_
ElectronMatchedCandidateProducer(const edm::ParameterSet &)
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
static std::atomic< unsigned int > counter
edm::EDGetTokenT< edm::View< reco::GsfElectron > > electronCollectionToken_
def move(src, dest)
Definition: eostools.py:511
Definition: event.py:1