CMS 3D CMS Logo

UniqueElectrons.h
Go to the documentation of this file.
2 
3 #include <vector>
4 
5 std::vector<reco::GsfElectronRef> uniqueElectronFinder(edm::Handle<reco::GsfElectronCollection>& pElectrons) {
6  const reco::GsfElectronCollection* electrons = pElectrons.product();
7  //Remove duplicate electrons which share a supercluster
8  std::vector<reco::GsfElectronRef> UniqueElectrons;
9  int index = 0;
10  for (reco::GsfElectronCollection::const_iterator elec = electrons->begin(); elec != electrons->end(); ++elec) {
11  const reco::GsfElectronRef electronRef(pElectrons, index);
12  reco::GsfElectronCollection::const_iterator BestDuplicate = elec;
13  for (reco::GsfElectronCollection::const_iterator elec2 = electrons->begin(); elec2 != electrons->end(); ++elec2) {
14  if (elec != elec2) {
15  if (elec->superCluster() == elec2->superCluster()) {
16  edm::LogDebug_("", "MySelection.cc", 122) << "e/p Best duplicate = " << BestDuplicate->eSuperClusterOverP()
17  << "\telec2 = " << elec2->eSuperClusterOverP();
18  if (fabs(BestDuplicate->eSuperClusterOverP() - 1.) >= fabs(elec2->eSuperClusterOverP() - 1.)) {
19  BestDuplicate = elec2;
20  edm::LogDebug_("", "MySelection.cc", 122) << "elec2 is now best duplicate";
21  } else
22  edm::LogDebug_("", "MySelection.cc", 122) << "BestDuplicate remains best duplicate";
23  }
24  }
25  }
26  if (BestDuplicate == elec)
27  UniqueElectrons.push_back(electronRef);
28  ++index;
29  }
30  return UniqueElectrons;
31 }
edm::Handle::product
T const * product() const
Definition: Handle.h:70
edm::LogDebug_
Definition: MessageLogger.h:510
reco::GsfElectronCollection
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
Definition: GsfElectronFwd.h:14
edm::Handle< reco::GsfElectronCollection >
edm::Ref
Definition: AssociativeIterator.h:58
uniqueElectronFinder
std::vector< reco::GsfElectronRef > uniqueElectronFinder(edm::Handle< reco::GsfElectronCollection > &pElectrons)
Definition: UniqueElectrons.h:5
GsfElectronFwd.h
pwdgSkimBPark_cfi.electrons
electrons
Definition: pwdgSkimBPark_cfi.py:6
AlignmentPI::index
index
Definition: AlignmentPayloadInspectorHelper.h:46