CMS 3D CMS Logo

ElectronMatchedCandidateProducer.cc
Go to the documentation of this file.
4 
5 #include "DataFormats/Math/interface/deltaR.h" // reco::deltaR
6 
7 
9 {
10 
11  const edm::InputTag allelectrons("gsfElectrons");
13  consumes<edm::View<reco::GsfElectron> >(params.getUntrackedParameter<edm::InputTag>("ReferenceElectronCollection",
14  allelectrons));
16  consumes<edm::View<reco::Candidate> >(params.getParameter<edm::InputTag>("src"));
17 
18  delRMatchingCut_ = params.getUntrackedParameter<double>("deltaR",
19  0.30);
20 
21  produces< edm::PtrVector<reco::Candidate> >();
22  produces< edm::RefToBaseVector<reco::Candidate> >();
23 }
24 
25 
26 
27 
29 {
30 
31 }
32 
33 
34 //
35 // member functions
36 //
37 
38 
39 // ------------ method called to produce the data ------------
40 
42  const edm::EventSetup &eventSetup)
43 {
44  // Create the output collection
45  auto outColRef = std::make_unique<edm::RefToBaseVector<reco::Candidate>>();
46  auto outColPtr = std::make_unique<edm::PtrVector<reco::Candidate>>();
47 
48 
49  // Read electrons
51  event.getByToken(electronCollectionToken_, electrons);
52 
53 
54 
55  //Read candidates
57  event.getByToken( scCollectionToken_ , recoCandColl);
58 
59 
60  unsigned int counter=0;
61 
62  // Loop over candidates
63  for(edm::View<reco::Candidate>::const_iterator scIt = recoCandColl->begin();
64  scIt != recoCandColl->end(); ++scIt, ++counter){
65  // Now loop over electrons
66  for(edm::View<reco::GsfElectron>::const_iterator elec = electrons->begin();
67  elec != electrons->end(); ++elec) {
68 
69  reco::SuperClusterRef eSC = elec->superCluster();
70 
71  double dRval = reco::deltaR((float)eSC->eta(), (float)eSC->phi(),
72  scIt->eta(), scIt->phi());
73 
74  if( dRval < delRMatchingCut_ ) {
75  //outCol->push_back( *scIt );
76  outColRef->push_back( recoCandColl->refAt(counter) );
77  outColPtr->push_back( recoCandColl->ptrAt(counter) );
78  } // end if loop
79  } // end electron loop
80 
81  } // end candidate loop
82 
83  event.put(std::move(outColRef));
84  event.put(std::move(outColPtr));
85 }
86 
87 
88 
89 
90 // ------ method called once each job just before starting event loop ---
91 
92 
93 
95 
96 
97 
99 
100 
101 
102 //define this as a plug-in
104 
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
virtual void produce(edm::Event &, const edm::EventSetup &) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
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:510
Definition: event.py:1