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