Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:23:47

0001 //
0002 //
0003 
0004 /**
0005   \class    pat::DuplicatedElectronCleaner DuplicatedElectronCleaner.h "PhysicsTools/PatAlgos/interface/DuplicatedElectronCleaner.h"
0006   \brief    Remove duplicates from the list of electrons
0007 
0008    The DuplicatedElectronCleaner removes duplicates from the input collection.
0009    Two electrons are considered duplicate if they share the same gsfTrack or the same superCluster.
0010    Among the two, the one with |E/P| nearest to 1 is kept.
0011    This is performed by the DuplicatedElectronRemover in PhysicsTools/PatUtils
0012 
0013    The output is an edm::RefVector<reco:::GsfElectron>,
0014    which can be read through edm::View<reco::GsfElectron>
0015 
0016   \author   Giovanni Petrucciani
0017   \version  $Id: DuplicatedElectronCleaner.cc,v 1.4 2010/02/20 21:00:16 wmtan Exp $
0018 */
0019 
0020 #include "FWCore/Framework/interface/Event.h"
0021 #include "FWCore/Framework/interface/global/EDProducer.h"
0022 #include "FWCore/Utilities/interface/InputTag.h"
0023 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0024 //#include "DataFormats/Common/interface/RefVector.h"
0025 #include "DataFormats/Common/interface/RefToBaseVector.h"
0026 //#include "DataFormats/Common/interface/PtrVector.h"
0027 
0028 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0029 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
0030 #include "PhysicsTools/PatUtils/interface/DuplicatedElectronRemover.h"
0031 
0032 namespace pat {
0033   class DuplicatedElectronCleaner : public edm::global::EDProducer<> {
0034   public:
0035     explicit DuplicatedElectronCleaner(const edm::ParameterSet& iConfig);
0036     ~DuplicatedElectronCleaner() override;
0037 
0038     void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const final;
0039 
0040   private:
0041     const edm::EDGetTokenT<edm::View<reco::GsfElectron>> electronSrcToken_;
0042     const pat::DuplicatedElectronRemover duplicateRemover_;
0043     mutable std::atomic<uint64_t> try_, pass_;
0044   };
0045 }  // namespace pat
0046 
0047 pat::DuplicatedElectronCleaner::DuplicatedElectronCleaner(const edm::ParameterSet& iConfig)
0048     : electronSrcToken_(consumes<edm::View<reco::GsfElectron>>(iConfig.getParameter<edm::InputTag>("electronSource"))),
0049       duplicateRemover_(),
0050       try_(0),
0051       pass_(0) {
0052   //produces<edm::RefVector<reco::GsfElectronCollection> >();
0053   produces<edm::RefToBaseVector<reco::GsfElectron>>();
0054   //produces<edm::PtrVector<reco::GsfElectron> >();
0055 }
0056 
0057 pat::DuplicatedElectronCleaner::~DuplicatedElectronCleaner() {}
0058 
0059 void pat::DuplicatedElectronCleaner::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0060   using namespace edm;
0061   Handle<View<reco::GsfElectron>> electrons;
0062   iEvent.getByToken(electronSrcToken_, electrons);
0063   try_ += electrons->size();
0064 
0065   //auto result = std::make_unique<RefVector<reco::GsfElectronCollection>>();
0066   auto result = std::make_unique<RefToBaseVector<reco::GsfElectron>>();
0067   //auto result = std::make_unique<PtrVector<reco::GsfElectron>>();
0068   std::unique_ptr<std::vector<size_t>> duplicates = duplicateRemover_.duplicatesToRemove(*electrons);
0069 
0070   std::vector<size_t>::const_iterator itdup = duplicates->begin(), enddup = duplicates->end();
0071   for (size_t i = 0, n = electrons->size(); i < n; ++i) {
0072     while ((itdup != enddup) && (*itdup < i)) {
0073       ++itdup;
0074     }
0075     if ((itdup != enddup) && (*itdup == i))
0076       continue;
0077     //result->push_back(electrons->refAt(i).castTo<edm::Ref<reco::GsfElectronCollection> >());
0078     result->push_back(electrons->refAt(i));
0079     //result->push_back(electrons->ptrAt(i));
0080   }
0081   pass_ += result->size();
0082   iEvent.put(std::move(result));
0083 }
0084 
0085 #include "FWCore/Framework/interface/MakerMacros.h"
0086 using pat::DuplicatedElectronCleaner;
0087 DEFINE_FWK_MODULE(DuplicatedElectronCleaner);