CMS 3D CMS Logo

HGCalElectronFilter.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: EgammaTools/HGCalElectronFilter
4 // Class: HGCalElectronFilter
5 //
13 //
14 // Original Author: Florian beaudette
15 // Created: Mon, 06 Nov 2017 21:49:54 GMT
16 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
26 
29 
32 
36 
37 //
38 // class declaration
39 //
40 
42  public:
43  explicit HGCalElectronFilter(const edm::ParameterSet&);
44  ~HGCalElectronFilter() override;
45 
46  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
47 
48  private:
49  void beginStream(edm::StreamID) override;
50  void produce(edm::Event&, const edm::EventSetup&) override;
51  void endStream() override;
52 
53  // ----------member data ---------------------------
57 };
58 
59 
61  electronsToken_(consumes<edm::View<reco::GsfElectron>>(iConfig.getParameter<edm::InputTag>("inputGsfElectrons"))),
62  outputCollection_(iConfig.getParameter<std::string>("outputCollection")),
63  cleanBarrel_(iConfig.getParameter<bool>("cleanBarrel"))
64 {
65  produces<reco::GsfElectronCollection>(outputCollection_);
66 }
67 
68 
70 }
71 
72 
73 //
74 // member functions
75 //
76 
77 // ------------ method called to produce the data ------------
78 void
80 
81  auto gsfElectrons_p = std::make_unique<reco::GsfElectronCollection>();
82 
84  iEvent.getByToken(electronsToken_, electronsH);
85 
86  for(const auto& electron1 : *electronsH) {
87  bool isBest = true;
88  if (!cleanBarrel_ && electron1.isEB()) { // keep all barrel electrons
89  gsfElectrons_p->push_back(electron1);
90  continue;
91  } else {
92  for (const auto& electron2 : *electronsH) {
93  if (&electron1 == &electron2) continue;
94  if (electron1.superCluster() != electron2.superCluster()) continue;
95  if (electron1.electronCluster() != electron2.electronCluster()) {
96  if (electron1.electronCluster()->energy() < electron2.electronCluster()->energy()) {
97  isBest=false;
98  break;
99  }
100  } else {
101  if (fabs(electron1.eEleClusterOverPout()-1.) > fabs(electron2.eEleClusterOverPout()-1.)) {
102  isBest=false;
103  break;
104  }
105  }
106  }
107  if (isBest) gsfElectrons_p->push_back(electron1);
108  }
109  }
110 
111  iEvent.put(std::move(gsfElectrons_p), outputCollection_);
112 }
113 
114 // ------------ method called once each stream before processing any runs, lumis or events ------------
115 void
117 {
118 }
119 
120 // ------------ method called once each stream after processing all runs, lumis and events ------------
121 void
123 }
124 
125 
126 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
127 void
129  // cleanedEcalDrivenGsfElectronsFromMultiCl
131  desc.add<edm::InputTag>("inputGsfElectrons", edm::InputTag("ecalDrivenGsfElectronsFromMultiCl"));
132  desc.add<bool>("cleanBarrel", false);
133  desc.add<std::string>("outputCollection", "");
134  descriptions.add("cleanedEcalDrivenGsfElectronsFromMultiCl", desc);
135 }
136 
137 //define this as a plug-in
void beginStream(edm::StreamID) override
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
HGCalElectronFilter(const edm::ParameterSet &)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< edm::View< reco::GsfElectron > > electronsToken_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
HLT enums.
void produce(edm::Event &, const edm::EventSetup &) override
def move(src, dest)
Definition: eostools.py:511