CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTSinglet.cc
Go to the documentation of this file.
1 
14 
21 
23 
24 #include <typeinfo>
25 
26 // extract the candidate type
27 template<typename T>
28 int getObjectType(const T &) {
29  return 0;
30 }
31 
32 // specialize for type l1extra::L1EmParticle
33 template<typename T>
34 int getObjectType(const l1extra::L1EmParticle & candidate) {
35  switch (candidate.type()) {
40  default:
41  return 0;
42  }
43 }
44 
45 // specialize for type l1extra::L1EtMissParticle
46 template<typename T>
47 int getObjectType(const l1extra::L1EtMissParticle & candidate) {
48  switch (candidate.type()) {
50  return trigger::TriggerL1ETM;
52  return trigger::TriggerL1HTM;
53  default:
54  return 0;
55  }
56 }
57 
58 // specialize for type l1extra::L1JetParticle
59 template<typename T>
60 int getObjectType(const l1extra::L1JetParticle & candidate) {
61  switch (candidate.type()) {
68  default:
69  return 0;
70  }
71 }
72 
73 
74 //
75 // constructors and destructor
76 //
77 template<typename T>
79  inputTag_ (iConfig.template getParameter<edm::InputTag>("inputTag")),
80  triggerType_ (iConfig.template getParameter<int>("triggerType")),
81  min_E_ (iConfig.template getParameter<double> ("MinE" )),
82  min_Pt_ (iConfig.template getParameter<double> ("MinPt" )),
83  min_Mass_ (iConfig.template getParameter<double> ("MinMass" )),
84  max_Eta_ (iConfig.template getParameter<double> ("MaxEta" )),
85  min_N_ (iConfig.template getParameter<int> ("MinN" )),
86  tid_ (triggerType_)
87 {
88  LogDebug("") << "Input/ptcut/etacut/ncut : "
89  << inputTag_.encode() << " "
90  << min_E_ << " " << min_Pt_ << " " << min_Mass_ << " "
91  << max_Eta_ << " " << min_N_ ;
92 }
93 
94 template<typename T>
96 {
97 }
98 
99 template<typename T>
100 void
103  makeHLTFilterDescription(desc);
104  desc.add<edm::InputTag>("inputTag",edm::InputTag("hltCollection"));
105  desc.add<int>("triggerType",0);
106  desc.add<double>("MinE",-1.0);
107  desc.add<double>("MinPt",-1.0);
108  desc.add<double>("MinMass",-1.0);
109  desc.add<double>("MaxEta",-1.0);
110  desc.add<int>("MinN",1);
111  descriptions.add(std::string("hlt")+std::string(typeid(HLTSinglet<T>).name()),desc);
112 }
113 
114 //
115 // member functions
116 //
117 
118 // ------------ method called to produce the data ------------
119 template<typename T>
120 bool
122 {
123  using namespace std;
124  using namespace edm;
125  using namespace reco;
126  using namespace trigger;
127 
128  typedef vector<T> TCollection;
129  typedef Ref<TCollection> TRef;
130 
131  // All HLT filters must create and fill an HLT filter object,
132  // recording any reconstructed physics objects satisfying (or not)
133  // this HLT filter, and place it in the Event.
134 
135  // The filter object
136  if (saveTags()) filterproduct.addCollectionTag(inputTag_);
137 
138  // Ref to Candidate object to be recorded in filter object
139  TRef ref;
140 
141 
142  // get hold of collection of objects
143  Handle<TCollection> objects;
144  iEvent.getByLabel (inputTag_,objects);
145 
146  // look at all objects, check cuts and add to filter object
147  int n(0);
148  typename TCollection::const_iterator i ( objects->begin() );
149  for (; i!=objects->end(); i++) {
150  if ( (i->energy() >= min_E_) &&
151  (i->pt() >= min_Pt_) &&
152  (i->mass() >= min_Mass_) &&
153  ( (max_Eta_ < 0.0) || (std::abs(i->eta()) <= max_Eta_) ) ) {
154  n++;
155  ref=TRef(objects,distance(objects->begin(),i));
156  tid_=getObjectType<T>(*i);
157  if (tid_==0) tid_=triggerType_;
158  filterproduct.addObject(tid_,ref);
159  }
160  }
161 
162  // filter decision
163  bool accept(n>=min_N_);
164 
165  return accept;
166 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
double min_E_
Definition: HLTSinglet.h:39
double min_Pt_
Definition: HLTSinglet.h:40
int getObjectType(const T &)
Definition: HLTSinglet.cc:28
#define abs(x)
Definition: mlp_lapack.h:159
EmType type() const
Definition: L1EmParticle.h:63
EtMissType type() const
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:22
std::string encode() const
Definition: InputTag.cc:72
JetType type() const
Definition: L1JetParticle.h:64
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
double min_Mass_
Definition: HLTSinglet.h:41
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
edm::InputTag inputTag_
Definition: HLTSinglet.h:37
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
virtual bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct)
Definition: HLTSinglet.cc:121
int min_N_
Definition: HLTSinglet.h:43
long double T
double max_Eta_
Definition: HLTSinglet.h:42
HLTSinglet(const edm::ParameterSet &)
Definition: HLTSinglet.cc:78
def template
Definition: svgfig.py:520
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: HLTSinglet.cc:101