CMS 3D CMS Logo

HLTSinglet.cc
Go to the documentation of this file.
1 
12 
19 
21 
23 
24 // extract the candidate type
25 template<typename T>
26 int getObjectType(const T &) {
27  return 0;
28 }
29 
30 // specialize for type l1extra::L1EmParticle
31 template<typename T>
32 int getObjectType(const l1extra::L1EmParticle & candidate) {
33  switch (candidate.type()) {
38  default:
39  return 0;
40  }
41 }
42 
43 // specialize for type l1extra::L1EtMissParticle
44 template<typename T>
45 int getObjectType(const l1extra::L1EtMissParticle & candidate) {
46  switch (candidate.type()) {
48  return trigger::TriggerL1ETM;
50  return trigger::TriggerL1HTM;
51  default:
52  return 0;
53  }
54 }
55 
56 // specialize for type l1extra::L1JetParticle
57 template<typename T>
58 int getObjectType(const l1extra::L1JetParticle & candidate) {
59  switch (candidate.type()) {
66  default:
67  return 0;
68  }
69 }
70 
71 
72 //
73 // constructors and destructor
74 //
75 template<typename T>
77  inputTag_ (iConfig.template getParameter<edm::InputTag>("inputTag")),
78  inputToken_ (consumes<std::vector<T> >(inputTag_)),
79  triggerType_ (iConfig.template getParameter<int>("triggerType")),
80  min_N_ (iConfig.template getParameter<int> ("MinN" )),
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  min_Eta_ (iConfig.template getParameter<double> ("MinEta" )),
85  max_Eta_ (iConfig.template getParameter<double> ("MaxEta" ))
86 {
87  LogDebug("") << "Input/ptcut/etacut/ncut : "
88  << inputTag_.encode() << " "
89  << min_E_ << " " << min_Pt_ << " " << min_Mass_ << " "
90  << min_Eta_ << " " << max_Eta_ << " " << min_N_ ;
91 }
92 
93 template<typename T>
94 HLTSinglet<T>::~HLTSinglet() = default;
95 
96 template<typename T>
97 void
101  desc.add<edm::InputTag>("inputTag",edm::InputTag("hltCollection"));
102  desc.add<int>("triggerType",0);
103  desc.add<double>("MinE",-1.0);
104  desc.add<double>("MinPt",-1.0);
105  desc.add<double>("MinMass",-1.0);
106  desc.add<double>("MinEta",-1.0);
107  desc.add<double>("MaxEta",-1.0);
108  desc.add<int>("MinN",1);
109  descriptions.add(defaultModuleLabel<HLTSinglet<T>>(), desc);
110 }
111 
112 //
113 // member functions
114 //
115 
116 // ------------ method called to produce the data ------------
117 template<typename T>
118 bool
120 {
121  using namespace std;
122  using namespace edm;
123  using namespace reco;
124  using namespace trigger;
125 
126  typedef vector<T> TCollection;
127  typedef Ref<TCollection> TRef;
128 
129  // All HLT filters must create and fill an HLT filter object,
130  // recording any reconstructed physics objects satisfying (or not)
131  // this HLT filter, and place it in the Event.
132 
133  // The filter object
134  if (saveTags()) filterproduct.addCollectionTag(inputTag_);
135 
136  // Ref to Candidate object to be recorded in filter object
137  TRef ref;
138 
139 
140  // get hold of collection of objects
142  iEvent.getByToken(inputToken_,objects);
143 
144  // look at all objects, check cuts and add to filter object
145  int n(0);
146  typename TCollection::const_iterator i ( objects->begin() );
147  for (; i!=objects->end(); i++) {
148  if ( (i->energy() >= min_E_) &&
149  (i->pt() >= min_Pt_) &&
150  (i->mass() >= min_Mass_) &&
151  ( (min_Eta_ < 0.0) || (std::abs(i->eta()) >= min_Eta_) ) &&
152  ( (max_Eta_ < 0.0) || (std::abs(i->eta()) <= max_Eta_) ) ) {
153  n++;
154  ref=TRef(objects,distance(objects->begin(),i));
155  int tid = getObjectType<T>(*i);
156  if (tid == 0)
157  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)
std::string defaultModuleLabel()
int getObjectType(const T &)
Definition: HLTSinglet.cc:26
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
EmType type() const
Definition: L1EmParticle.h:62
const double min_E_
Definition: HLTSinglet.h:39
EtMissType type() const
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
std::string encode() const
Definition: InputTag.cc:165
JetType type() const
Definition: L1JetParticle.h:63
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
virtual bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
Definition: HLTSinglet.cc:119
int iEvent
Definition: GenABIO.cc:230
const double max_Eta_
Definition: HLTSinglet.h:43
const double min_Eta_
Definition: HLTSinglet.h:42
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
def template(fileName, svg, replaceme="REPLACEME")
Definition: svgfig.py:520
const double min_Mass_
Definition: HLTSinglet.h:41
const double min_Pt_
Definition: HLTSinglet.h:40
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
const int min_N_
Definition: HLTSinglet.h:38
const int triggerType_
Definition: HLTSinglet.h:37
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
bool saveTags() const
Definition: HLTFilter.h:45
HLT enums.
const edm::EDGetTokenT< std::vector< T > > inputToken_
Definition: HLTSinglet.h:36
const edm::InputTag inputTag_
Definition: HLTSinglet.h:35
long double T
HLTSinglet(const edm::ParameterSet &)
Definition: HLTSinglet.cc:76
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: HLTSinglet.cc:98