CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTElectronGenericFilter.cc
Go to the documentation of this file.
1 
10 
12 
14 
16 
22 
25 
26 //
27 // constructors and destructor
28 //
30  candTag_ = iConfig.getParameter< edm::InputTag > ("candTag");
31  isoTag_ = iConfig.getParameter< edm::InputTag > ("isoTag");
32  nonIsoTag_ = iConfig.getParameter< edm::InputTag > ("nonIsoTag");
33  lessThan_ = iConfig.getParameter<bool> ("lessThan");
34  thrRegularEB_ = iConfig.getParameter<double> ("thrRegularEB");
35  thrRegularEE_ = iConfig.getParameter<double> ("thrRegularEE");
36  thrOverPtEB_ = iConfig.getParameter<double> ("thrOverPtEB");
37  thrOverPtEE_ = iConfig.getParameter<double> ("thrOverPtEE");
38  thrTimesPtEB_ = iConfig.getParameter<double> ("thrTimesPtEB");
39  thrTimesPtEE_ = iConfig.getParameter<double> ("thrTimesPtEE");
40  ncandcut_ = iConfig.getParameter<int> ("ncandcut");
41  doIsolated_ = iConfig.getParameter<bool> ("doIsolated");
42  L1IsoCollTag_= iConfig.getParameter< edm::InputTag > ("L1IsoCand");
43  L1NonIsoCollTag_= iConfig.getParameter< edm::InputTag > ("L1NonIsoCand");
44 }
45 
47 
48 
49 // ------------ method called to produce the data ------------
50 bool
52 {
53  using namespace trigger;
54 
55  if (saveTags()) {
56  filterproduct.addCollectionTag(L1IsoCollTag_);
57  if (not doIsolated_) filterproduct.addCollectionTag(L1NonIsoCollTag_);
58  }
59 
60  // Ref to Candidate object to be recorded in filter object
62 
64 
65  iEvent.getByLabel (candTag_,PrevFilterOutput);
66 
67  std::vector<edm::Ref<reco::ElectronCollection> > elecands;
68  PrevFilterOutput->getObjects(TriggerElectron, elecands);
69 
70 
71  //get hold of isolated association map
73  iEvent.getByLabel (isoTag_,depMap);
74 
75  //get hold of non-isolated association map
77  if(!doIsolated_) iEvent.getByLabel (nonIsoTag_,depNonIsoMap);
78 
79  // look at all photons, check cuts and add to filter object
80  int n = 0;
81 
82  for (unsigned int i=0; i<elecands.size(); i++) {
83 
84  ref = elecands[i];
85  reco::ElectronIsolationMap::const_iterator mapi = (*depMap).find( ref );
86  if (mapi==(*depMap).end() && !doIsolated_) mapi = (*depNonIsoMap).find( ref );
87 
88  float vali = mapi->val;
89  float Pt = ref->pt();
90  float Eta = fabs(ref->eta());
91 
92  if ( lessThan_ ) {
93  if ( (Eta < 1.479 && vali <= thrRegularEB_) || (Eta >= 1.479 && vali <= thrRegularEE_) ) {
94  n++;
95  filterproduct.addObject(TriggerElectron, ref);
96  continue;
97  }
98  if (Pt > 0. && (thrOverPtEB_ > 0. || thrOverPtEE_ > 0. || thrTimesPtEB_ > 0. || thrTimesPtEE_ > 0.) ) {
99  if ((Eta < 1.479 && vali/Pt <= thrOverPtEB_) || (Eta >= 1.479 && vali/Pt <= thrOverPtEE_) ) {
100  n++;
101  filterproduct.addObject(TriggerElectron, ref);
102  continue;
103  }
104  if ((Eta < 1.479 && vali*Pt <= thrTimesPtEB_) || (Eta >= 1.479 && vali*Pt <= thrTimesPtEE_) ) {
105  n++;
106  filterproduct.addObject(TriggerElectron, ref);
107  }
108  }
109  } else {
110  if ( (Eta < 1.479 && vali >= thrRegularEB_) || (Eta >= 1.479 && vali >= thrRegularEE_) ) {
111  n++;
112  filterproduct.addObject(TriggerElectron, ref);
113  continue;
114  }
115  if (Pt > 0. && (thrOverPtEB_ > 0. || thrOverPtEE_ > 0. || thrTimesPtEB_ > 0. || thrTimesPtEE_ > 0.) ) {
116  if ((Eta < 1.479 && vali/Pt >= thrOverPtEB_) || (Eta >= 1.479 && vali/Pt >= thrOverPtEE_) ) {
117  n++;
118  filterproduct.addObject(TriggerElectron, ref);
119  continue;
120  }
121  if ((Eta < 1.479 && vali*Pt >= thrTimesPtEB_) || (Eta >= 1.479 && vali*Pt >= thrTimesPtEE_) ) {
122  n++;
123  filterproduct.addObject(TriggerElectron, ref);
124  }
125  }
126  }
127  }
128 
129  // filter decision
130  bool accept(n>=ncandcut_);
131 
132  return accept;
133 }
134 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:22
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref&lt;C&gt;)
int iEvent
Definition: GenABIO.cc:243
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
bool saveTags() const
Definition: HLTFilter.h:45
HLTElectronGenericFilter(const edm::ParameterSet &)
virtual bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct)