CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HLTSinglet.cc
Go to the documentation of this file.
1 
12 
19 
20 #include "HLTSinglet.h"
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>
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 // constructors and destructor
73 //
74 template <typename T>
76  : HLTFilter(iConfig),
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  max_Mass_(iConfig.template getParameter<double>("MaxMass")),
85  min_Eta_(iConfig.template getParameter<double>("MinEta")),
86  max_Eta_(iConfig.template getParameter<double>("MaxEta")) {
87  LogDebug("") << "Input/ptcut/etacut/ncut : " << inputTag_.encode() << " " << min_E_ << " " << min_Pt_ << " "
88  << min_Mass_ << " " << max_Mass_ << " " << min_Eta_ << " " << max_Eta_ << " " << min_N_;
89 }
90 
91 template <typename T>
92 HLTSinglet<T>::~HLTSinglet() = default;
93 
94 template <typename T>
97  makeHLTFilterDescription(desc);
98  desc.add<edm::InputTag>("inputTag", edm::InputTag("hltCollection"));
99  desc.add<int>("triggerType", 0);
100  desc.add<double>("MinE", -1.0);
101  desc.add<double>("MinPt", -1.0);
102  desc.add<double>("MinMass", -1.0);
103  desc.add<double>("MaxMass", -1.0);
104  desc.add<double>("MinEta", -1.0);
105  desc.add<double>("MaxEta", -1.0);
106  desc.add<int>("MinN", 1);
107  descriptions.add(defaultModuleLabel<HLTSinglet<T>>(), desc);
108 }
109 
110 //
111 // member functions
112 //
113 
114 // ------------ method called to produce the data ------------
115 template <typename T>
117  const edm::EventSetup& iSetup,
118  trigger::TriggerFilterObjectWithRefs& filterproduct) const {
119  using namespace std;
120  using namespace edm;
121  using namespace reco;
122  using namespace trigger;
123 
124  typedef vector<T> TCollection;
125  typedef Ref<TCollection> TRef;
126 
127  // All HLT filters must create and fill an HLT filter object,
128  // recording any reconstructed physics objects satisfying (or not)
129  // this HLT filter, and place it in the Event.
130 
131  // The filter object
132  if (saveTags())
133  filterproduct.addCollectionTag(inputTag_);
134 
135  // Ref to Candidate object to be recorded in filter object
136  TRef ref;
137 
138  // get hold of collection of objects
140  iEvent.getByToken(inputToken_, objects);
141 
142  // look at all objects, check cuts and add to filter object
143  int n(0);
144  typename TCollection::const_iterator i(objects->begin());
145  for (; i != objects->end(); i++) {
146  if ((i->energy() >= min_E_) && (i->pt() >= min_Pt_) && (i->mass() >= min_Mass_) &&
147  ((max_Mass_ < 0.0) || (i->mass() <= max_Mass_)) && ((min_Eta_ < 0.0) || (std::abs(i->eta()) >= min_Eta_)) &&
148  ((max_Eta_ < 0.0) || (std::abs(i->eta()) <= max_Eta_))) {
149  n++;
150  ref = TRef(objects, distance(objects->begin(), i));
151  int tid = getObjectType<T>(*i);
152  if (tid == 0)
153  tid = triggerType_;
154  filterproduct.addObject(tid, ref);
155  }
156  }
157 
158  // filter decision
159  bool accept(n >= min_N_);
160 
161  return accept;
162 }
int getObjectType(const T &)
Definition: HLTSinglet.cc:26
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
const double max_Mass_
Definition: HLTSinglet.h:43
EmType type() const
Definition: L1EmParticle.h:46
std::string defaultModuleLabel()
const double min_E_
Definition: HLTSinglet.h:40
EtMissType type() const
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:31
std::string encode() const
Definition: InputTag.cc:159
JetType type() const
Definition: L1JetParticle.h:46
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:224
const double max_Eta_
Definition: HLTSinglet.h:45
const double min_Eta_
Definition: HLTSinglet.h:44
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const double min_Mass_
Definition: HLTSinglet.h:42
const double min_Pt_
Definition: HLTSinglet.h:41
ParameterDescriptionBase * add(U const &iLabel, T const &value)
const int min_N_
Definition: HLTSinglet.h:39
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
Definition: HLTSinglet.cc:116
void add(std::string const &label, ParameterSetDescription const &psetDescription)
~HLTSinglet() override
const edm::InputTag inputTag_
Definition: HLTSinglet.h:36
long double T
HLTSinglet(const edm::ParameterSet &)
Definition: HLTSinglet.cc:75
def template
Definition: svgfig.py:521
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: HLTSinglet.cc:95
#define LogDebug(id)