CMS 3D CMS Logo

PackedTriggerPrescales.h
Go to the documentation of this file.
1 #ifndef DataFormats_PatCandidates_PackedTriggerPrescales_h
2 #define DataFormats_PatCandidates_PackedTriggerPrescales_h
3 
4 #include <cstring>
5 #include <type_traits>
6 
10 
11 namespace pat {
12 
14  public:
17  ~PackedTriggerPrescales() = default;
18 
19  // get prescale by index
20  // method templated to force correct choice of output type
21  // (as part of deprecating integer types for trigger prescales)
22  template <typename T = int>
23  T getPrescaleForIndex(int index) const;
24 
25  // get prescale by name or name prefix (if setTriggerNames was called)
26  // method templated to force correct choice of output type
27  // (as part of deprecating integer types for trigger prescales)
28  template <typename T = int>
29  T getPrescaleForName(const std::string &name, bool prefixOnly = false) const;
30 
31  // return the TriggerResults associated with this
32  const edm::TriggerResults &triggerResults() const { return *edm::getProduct<edm::TriggerResults>(triggerResults_); }
33 
34  // use this method first if you want to be able to access the prescales by name
35  // you can get the TriggerNames from the TriggerResults and the Event (edm or fwlite)
37 
38  // set that the trigger of given index has a given prescale
39  void addPrescaledTrigger(int index, double prescale);
40 
41  protected:
42  std::vector<double> prescaleValues_;
45  };
46 
47  template <typename T>
49  static_assert(std::is_same_v<T, double>,
50  "\n\n\tPlease use getPrescaleForIndex<double> "
51  "(other types for trigger prescales are not supported anymore by PackedTriggerPrescales)");
52  if (unsigned(index) >= triggerResults().size())
53  throw cms::Exception("InvalidReference", "Index out of bounds");
54  return prescaleValues_[index];
55  }
56 
57  template <typename T>
59  static_assert(std::is_same_v<T, double>,
60  "\n\n\tPlease use getPrescaleForName<double> "
61  "(other types for trigger prescales are not supported anymore by PackedTriggerPrescales)");
62  if (triggerNames_ == nullptr)
63  throw cms::Exception("LogicError", "getPrescaleForName called without having called setTriggerNames first");
64  if (prefixOnly) {
65  auto const &names = triggerNames_->triggerNames();
66  if (name.empty())
67  throw cms::Exception("EmptyName",
68  "getPrescaleForName called with invalid arguments (name is empty, prefixOnly=true");
69  size_t siz = name.length() - 1;
70  while (siz > 0 && (name[siz] == '*' || name[siz] == '\0'))
71  siz--;
72  for (unsigned int i = 0, n = names.size(); i < n; ++i) {
73  if (strncmp(name.c_str(), names[i].c_str(), siz) == 0) {
74  return getPrescaleForIndex<T>(i);
75  }
76  }
77  throw cms::Exception("InvalidReference", "Index out of bounds");
78  } else {
80  return getPrescaleForIndex<T>(index);
81  }
82  }
83 
84 } // namespace pat
85 
86 #endif // DataFormats_PatCandidates_PackedTriggerPrescales_h
size
Write out results.
T getPrescaleForName(const std::string &name, bool prefixOnly=false) const
Strings const & triggerNames() const
Definition: TriggerNames.cc:48
const std::string names[nVars_]
const edm::TriggerResults & triggerResults() const
Definition: HeavyIon.h:7
unsigned int triggerIndex(std::string_view name) const
Definition: TriggerNames.cc:52
std::vector< double > prescaleValues_
T getPrescaleForIndex(int index) const
const edm::TriggerNames * triggerNames_
void setTriggerNames(const edm::TriggerNames &names)
void addPrescaledTrigger(int index, double prescale)
long double T