CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DuplicatedElectronCleaner.cc
Go to the documentation of this file.
1 //
2 //
3 
24 //#include "DataFormats/Common/interface/RefVector.h"
26 //#include "DataFormats/Common/interface/PtrVector.h"
27 
31 
32 namespace pat{
34  public:
35  explicit DuplicatedElectronCleaner(const edm::ParameterSet & iConfig);
37 
38  virtual void produce(edm::StreamID, edm::Event & iEvent, const edm::EventSetup& iSetup) const override final;
39 
40  private:
43  mutable std::atomic<uint64_t> try_, pass_;
44  };
45 } // namespace
46 
48  electronSrcToken_(consumes<edm::View<reco::GsfElectron> >(iConfig.getParameter<edm::InputTag>("electronSource"))),
49  duplicateRemover_(),
50  try_(0), pass_(0)
51 {
52  //produces<edm::RefVector<reco::GsfElectronCollection> >();
53  produces<edm::RefToBaseVector<reco::GsfElectron> >();
54  //produces<edm::PtrVector<reco::GsfElectron> >();
55 }
56 
58 {
59 }
60 
61 void
63 {
64  using namespace edm;
66  iEvent.getByToken(electronSrcToken_, electrons);
67  try_ += electrons->size();
68 
69  //std::auto_ptr<RefVector<reco::GsfElectronCollection> > result(new RefVector<reco::GsfElectronCollection>());
70  std::auto_ptr<RefToBaseVector<reco::GsfElectron> > result(new RefToBaseVector<reco::GsfElectron>());
71  //std::auto_ptr<PtrVector<reco::GsfElectron> > result(new PtrVector<reco::GsfElectron>());
72  std::auto_ptr< std::vector<size_t> > duplicates = duplicateRemover_.duplicatesToRemove(*electrons);
73 
74  std::vector<size_t>::const_iterator itdup = duplicates->begin(), enddup = duplicates->end();
75  for (size_t i = 0, n = electrons->size(); i < n; ++i) {
76  while ((itdup != enddup) && (*itdup < i)) { ++itdup; }
77  if ((itdup != enddup) && (*itdup == i)) continue;
78  //result->push_back(electrons->refAt(i).castTo<edm::Ref<reco::GsfElectronCollection> >());
79  result->push_back(electrons->refAt(i));
80  //result->push_back(electrons->ptrAt(i));
81  }
82  pass_ += result->size();
83  iEvent.put(result);
84 }
85 
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Remove duplicates from the list of electrons.
virtual void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const overridefinal
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
tuple result
Definition: query.py:137
DuplicatedElectronCleaner(const edm::ParameterSet &iConfig)
const edm::EDGetTokenT< edm::View< reco::GsfElectron > > electronSrcToken_
const pat::DuplicatedElectronRemover duplicateRemover_