00001 // 00002 // $Id: DuplicatedElectronCleaner.cc,v 1.4 2010/02/20 21:00:16 wmtan Exp $ 00003 // 00004 00021 #include "FWCore/Framework/interface/Event.h" 00022 #include "FWCore/Framework/interface/EDProducer.h" 00023 #include "FWCore/Utilities/interface/InputTag.h" 00024 #include "FWCore/ParameterSet/interface/ParameterSet.h" 00025 //#include "DataFormats/Common/interface/RefVector.h" 00026 #include "DataFormats/Common/interface/RefToBaseVector.h" 00027 //#include "DataFormats/Common/interface/PtrVector.h" 00028 00029 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h" 00030 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h" 00031 #include "PhysicsTools/PatUtils/interface/DuplicatedElectronRemover.h" 00032 00033 namespace pat{ 00034 class DuplicatedElectronCleaner : public edm::EDProducer{ 00035 public: 00036 explicit DuplicatedElectronCleaner(const edm::ParameterSet & iConfig); 00037 ~DuplicatedElectronCleaner(); 00038 00039 virtual void produce(edm::Event & iEvent, const edm::EventSetup & iSetup); 00040 virtual void endJob(); 00041 00042 private: 00043 edm::InputTag electronSrc_; 00044 pat::DuplicatedElectronRemover duplicateRemover_; 00045 uint64_t try_, pass_; 00046 }; 00047 } // namespace 00048 00049 pat::DuplicatedElectronCleaner::DuplicatedElectronCleaner(const edm::ParameterSet & iConfig): 00050 electronSrc_(iConfig.getParameter<edm::InputTag>("electronSource")), 00051 try_(0), pass_(0) 00052 { 00053 //produces<edm::RefVector<reco::GsfElectronCollection> >(); 00054 produces<edm::RefToBaseVector<reco::GsfElectron> >(); 00055 //produces<edm::PtrVector<reco::GsfElectron> >(); 00056 } 00057 00058 pat::DuplicatedElectronCleaner::~DuplicatedElectronCleaner() 00059 { 00060 } 00061 00062 void 00063 pat::DuplicatedElectronCleaner::produce(edm::Event & iEvent, const edm::EventSetup & iSetup) 00064 { 00065 using namespace edm; 00066 Handle<View<reco::GsfElectron> > electrons; 00067 iEvent.getByLabel(electronSrc_, electrons); 00068 try_ += electrons->size(); 00069 00070 //std::auto_ptr<RefVector<reco::GsfElectronCollection> > result(new RefVector<reco::GsfElectronCollection>()); 00071 std::auto_ptr<RefToBaseVector<reco::GsfElectron> > result(new RefToBaseVector<reco::GsfElectron>()); 00072 //std::auto_ptr<PtrVector<reco::GsfElectron> > result(new PtrVector<reco::GsfElectron>()); 00073 std::auto_ptr< std::vector<size_t> > duplicates = duplicateRemover_.duplicatesToRemove(*electrons); 00074 00075 std::vector<size_t>::const_iterator itdup = duplicates->begin(), enddup = duplicates->end(); 00076 for (size_t i = 0, n = electrons->size(); i < n; ++i) { 00077 while ((itdup != enddup) && (*itdup < i)) { ++itdup; } 00078 if ((itdup != enddup) && (*itdup == i)) continue; 00079 //result->push_back(electrons->refAt(i).castTo<edm::Ref<reco::GsfElectronCollection> >()); 00080 result->push_back(electrons->refAt(i)); 00081 //result->push_back(electrons->ptrAt(i)); 00082 } 00083 pass_ += result->size(); 00084 iEvent.put(result); 00085 } 00086 00087 void 00088 pat::DuplicatedElectronCleaner::endJob() 00089 { 00090 } 00091 00092 #include "FWCore/Framework/interface/MakerMacros.h" 00093 using pat::DuplicatedElectronCleaner; 00094 DEFINE_FWK_MODULE(DuplicatedElectronCleaner);