CMS 3D CMS Logo

HLTElectronEtFilter.cc
Go to the documentation of this file.
1 
8 #include "HLTElectronEtFilter.h"
9 
11 
15 
21 
24 
25 //
26 // constructors and destructor
27 //
29  candTag_ = iConfig.getParameter< edm::InputTag > ("candTag");
30  EtEB_ = iConfig.getParameter<double> ("EtCutEB");
31  EtEE_ = iConfig.getParameter<double> ("EtCutEE");
32 
33  ncandcut_ = iConfig.getParameter<int> ("ncandcut");
34 
35  l1EGTag_ = iConfig.getParameter< edm::InputTag > ("l1EGCand");
36 
37  candToken_ = consumes<trigger::TriggerFilterObjectWithRefs> (candTag_);
38 }
39 
41 
42 void
46  desc.add<edm::InputTag>("candTag", edm::InputTag("hltElectronPixelMatchFilter"));
47  desc.add<double>("EtCutEB", 0.0);
48  desc.add<double>("EtCutEE", 0.0);
49  desc.add<int>("ncandcut", 1);
50  desc.add<edm::InputTag>("l1EGCand", edm::InputTag("hltL1IsoRecoEcalCandidate"));
51  descriptions.add("hltElectronEtFilter", desc);
52 }
53 
54 // ------------ method called to produce the data ------------
56 {
57  using namespace trigger;
58  if (saveTags()) {
59  filterproduct.addCollectionTag(l1EGTag_);
60  }
61 
62  // Ref to Candidate object to be recorded in filter object
64 
66  iEvent.getByToken (candToken_, PrevFilterOutput);
67 
68  std::vector<edm::Ref<reco::ElectronCollection> > elecands;
69  PrevFilterOutput->getObjects(TriggerElectron, elecands);
70 
71 
72 
73  // look at all photons, check cuts and add to filter object
74  int n = 0;
75 
76  for (auto & elecand : elecands) {
77 
78  ref = elecand;
79  float Pt = ref->pt();
80  float Eta = fabs(ref->eta());
81 
82  if ( (Eta < 1.479 && Pt > EtEB_) || (Eta >= 1.479 && Pt > EtEE_) ) {
83  n++;
84  filterproduct.addObject(TriggerElectron, ref);
85  }
86 
87  }
88  // filter decision
89  bool accept(n>=ncandcut_);
90 
91  return accept;
92 }
93 
94 // declare this class as a framework plugin
T getParameter(std::string const &) const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > candToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
~HLTElectronEtFilter() override
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
HLTElectronEtFilter(const edm::ParameterSet &)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool saveTags() const
Definition: HLTFilter.h:45