Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:04

0001 // -*- C++ -*-
0002 //
0003 // Package:    EgammaTools/HGCalElectronFilter
0004 // Class:      HGCalElectronFilter
0005 //
0006 /**\class HGCalElectronFilter HGCalElectronFilter.cc EgammaTools/HGCalElectronFilter/plugins/HGCalElectronFilter.cc
0007 
0008  Description: filtering of duplicate HGCAL electrons
0009 
0010  Implementation:
0011      [Notes on implementation]
0012 */
0013 //
0014 // Original Author:  Florian beaudette
0015 //         Created:  Mon, 06 Nov 2017 21:49:54 GMT
0016 //
0017 //
0018 
0019 // system include files
0020 #include <memory>
0021 
0022 // user include files
0023 #include "FWCore/Framework/interface/Frameworkfwd.h"
0024 #include "FWCore/Framework/interface/stream/EDProducer.h"
0025 
0026 #include "FWCore/Framework/interface/Event.h"
0027 #include "FWCore/Framework/interface/MakerMacros.h"
0028 
0029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0030 #include "FWCore/Utilities/interface/StreamID.h"
0031 
0032 #include "DataFormats/Math/interface/deltaR.h"
0033 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0034 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
0035 
0036 //
0037 // class declaration
0038 //
0039 
0040 class HGCalElectronFilter : public edm::stream::EDProducer<> {
0041 public:
0042   explicit HGCalElectronFilter(const edm::ParameterSet&);
0043   ~HGCalElectronFilter() override;
0044 
0045   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0046 
0047 private:
0048   void beginStream(edm::StreamID) override;
0049   void produce(edm::Event&, const edm::EventSetup&) override;
0050   void endStream() override;
0051 
0052   // ----------member data ---------------------------
0053   edm::EDGetTokenT<edm::View<reco::GsfElectron>> electronsToken_;
0054   std::string outputCollection_;
0055   bool cleanBarrel_;
0056 };
0057 
0058 HGCalElectronFilter::HGCalElectronFilter(const edm::ParameterSet& iConfig)
0059     : electronsToken_(consumes(iConfig.getParameter<edm::InputTag>("inputGsfElectrons"))),
0060       outputCollection_(iConfig.getParameter<std::string>("outputCollection")),
0061       cleanBarrel_(iConfig.getParameter<bool>("cleanBarrel")) {
0062   produces<reco::GsfElectronCollection>(outputCollection_);
0063 }
0064 
0065 HGCalElectronFilter::~HGCalElectronFilter() {}
0066 
0067 //
0068 // member functions
0069 //
0070 
0071 // ------------ method called to produce the data  ------------
0072 void HGCalElectronFilter::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0073   auto gsfElectrons_p = std::make_unique<reco::GsfElectronCollection>();
0074 
0075   auto electronsH = iEvent.getHandle(electronsToken_);
0076 
0077   for (const auto& electron1 : *electronsH) {
0078     bool isBest = true;
0079     if (!cleanBarrel_ && electron1.isEB()) {  // keep all barrel electrons
0080       gsfElectrons_p->push_back(electron1);
0081       continue;
0082     } else {
0083       for (const auto& electron2 : *electronsH) {
0084         if (&electron1 == &electron2)
0085           continue;
0086         if (electron1.superCluster() != electron2.superCluster())
0087           continue;
0088         if (electron1.electronCluster() != electron2.electronCluster()) {
0089           if (electron1.electronCluster()->energy() < electron2.electronCluster()->energy()) {
0090             isBest = false;
0091             break;
0092           }
0093         } else {
0094           if (fabs(electron1.eEleClusterOverPout() - 1.) > fabs(electron2.eEleClusterOverPout() - 1.)) {
0095             isBest = false;
0096             break;
0097           }
0098         }
0099       }
0100       if (isBest)
0101         gsfElectrons_p->push_back(electron1);
0102     }
0103   }
0104 
0105   iEvent.put(std::move(gsfElectrons_p), outputCollection_);
0106 }
0107 
0108 // ------------ method called once each stream before processing any runs, lumis or events  ------------
0109 void HGCalElectronFilter::beginStream(edm::StreamID) {}
0110 
0111 // ------------ method called once each stream after processing all runs, lumis and events  ------------
0112 void HGCalElectronFilter::endStream() {}
0113 
0114 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0115 void HGCalElectronFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0116   // cleanedEcalDrivenGsfElectronsHGC
0117   edm::ParameterSetDescription desc;
0118   desc.add<edm::InputTag>("inputGsfElectrons", edm::InputTag("ecalDrivenGsfElectronsHGC"));
0119   desc.add<bool>("cleanBarrel", false);
0120   desc.add<std::string>("outputCollection", "");
0121   descriptions.add("cleanedEcalDrivenGsfElectronsHGC", desc);
0122 }
0123 
0124 //define this as a plug-in
0125 DEFINE_FWK_MODULE(HGCalElectronFilter);