CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HLTPrescaleProvider.h
Go to the documentation of this file.
1 #ifndef HLTcore_HLTPrescaleProvider_h
2 #define HLTcore_HLTPrescaleProvider_h
3 
24 
25 #include <memory>
26 #include <string>
27 #include <utility>
28 #include <vector>
29 
30 namespace edm {
31  class ConsumesCollector;
32  class Event;
33  class EventSetup;
34  class ParameterSet;
35  class Run;
36 } // namespace edm
37 
39 public:
40  template <typename T>
42 
43  template <typename T>
45 
50  bool init(const edm::Run& iRun, const edm::EventSetup& iSetup, const std::string& processName, bool& changed);
51 
53  L1GtUtils const& l1GtUtils() const;
54  l1t::L1TGlobalUtil const& l1tGlobalUtil() const;
55 
58  int prescaleSet(const edm::Event& iEvent, const edm::EventSetup& iSetup);
59  // negative == error
60 
62  template <typename T = unsigned int>
63  T prescaleValue(const edm::Event& iEvent, const edm::EventSetup& iSetup, const std::string& trigger) {
64  const int set(prescaleSet(iEvent, iSetup));
65  //there is a template specialisation for unsigned in which returns +1 which
66  //emulates old behaviour
67  return set < 0 ? -1 : hltConfigProvider_.prescaleValue<T>(static_cast<unsigned int>(set), trigger);
68  }
69 
71  template <typename TL1 = int, typename THLT = TL1>
72  std::pair<TL1, THLT> prescaleValues(const edm::Event& iEvent,
73  const edm::EventSetup& iSetup,
74  const std::string& trigger) {
75  return {convertL1PS<TL1>(getL1PrescaleValue(iEvent, iSetup, trigger)),
76  prescaleValue<THLT>(iEvent, iSetup, trigger)};
77  }
78  // any one negative => error in retrieving this (L1T or HLT) prescale
79 
80  // In case of a complex Boolean expression as L1 seed
81  template <typename TL1 = int, typename THLT = TL1>
82  std::pair<std::vector<std::pair<std::string, TL1> >, THLT> prescaleValuesInDetail(const edm::Event& iEvent,
83  const edm::EventSetup& iSetup,
84  const std::string& trigger) {
85  std::pair<std::vector<std::pair<std::string, TL1> >, THLT> retval;
86  for (auto& entry : getL1PrescaleValueInDetail(iEvent, iSetup, trigger)) {
87  retval.first.emplace_back(std::move(entry.first), convertL1PS<TL1>(entry.second));
88  }
89  retval.second = prescaleValue<THLT>(iEvent, iSetup, trigger);
90  return retval;
91  }
92  // Event rejected by HLTPrescaler on ith HLT path?
93  bool rejectedByHLTPrescaler(const edm::TriggerResults& triggerResults, unsigned int i) const;
95 
96 private:
97  void checkL1GtUtils() const;
98  void checkL1TGlobalUtil() const;
99  template <typename T>
100  T convertL1PS(double val) const {
101  return T(val);
102  }
103 
104  double getL1PrescaleValue(const edm::Event& iEvent, const edm::EventSetup& iSetup, const std::string& trigger);
105  std::vector<std::pair<std::string, double> > getL1PrescaleValueInDetail(const edm::Event& iEvent,
106  const edm::EventSetup& iSetup,
107  const std::string& trigger);
108  static constexpr int kL1PrescaleDenominator_ = 100;
110  std::unique_ptr<L1GtUtils> l1GtUtils_;
111  std::unique_ptr<l1t::L1TGlobalUtil> l1tGlobalUtil_;
112  unsigned char count_[5] = {0, 0, 0, 0, 0};
113  bool inited_ = false;
114 };
115 
116 template <typename T>
118  : HLTPrescaleProvider(pset, iC, module) {}
119 
120 template <typename T>
122  unsigned int stageL1Trigger = pset.getParameter<unsigned int>("stageL1Trigger");
123  if (stageL1Trigger <= 1) {
124  l1GtUtils_ = std::make_unique<L1GtUtils>(pset, iC, false, module, L1GtUtils::UseEventSetupIn::Run);
125  } else {
126  l1tGlobalUtil_ = std::make_unique<l1t::L1TGlobalUtil>(pset, iC, module, l1t::UseEventSetupIn::Run);
127  }
128 }
129 
130 template <>
131 FractionalPrescale HLTPrescaleProvider::convertL1PS<FractionalPrescale>(double val) const;
132 
133 template <>
134 unsigned int HLTPrescaleProvider::prescaleValue<unsigned int>(const edm::Event& iEvent,
135  const edm::EventSetup& iSetup,
136  const std::string& trigger);
137 
138 #endif
HLTConfigProvider hltConfigProvider_
static int l1PrescaleDenominator()
T convertL1PS(double val) const
std::unique_ptr< L1GtUtils > l1GtUtils_
T prescaleValue(const edm::Event &iEvent, const edm::EventSetup &iSetup, const std::string &trigger)
combining the two methods above
std::pair< std::vector< std::pair< std::string, TL1 > >, THLT > prescaleValuesInDetail(const edm::Event &iEvent, const edm::EventSetup &iSetup, const std::string &trigger)
double getL1PrescaleValue(const edm::Event &iEvent, const edm::EventSetup &iSetup, const std::string &trigger)
boost::rational< int > FractionalPrescale
int iEvent
Definition: GenABIO.cc:224
std::pair< TL1, THLT > prescaleValues(const edm::Event &iEvent, const edm::EventSetup &iSetup, const std::string &trigger)
Combined L1T (pair.first) and HLT (pair.second) prescales per HLT path.
def move
Definition: eostools.py:511
static std::string const triggerResults
Definition: EdmProvDump.cc:44
void checkL1TGlobalUtil() const
bool rejectedByHLTPrescaler(const edm::TriggerResults &triggerResults, unsigned int i) const
int prescaleSet(const edm::Event &iEvent, const edm::EventSetup &iSetup)
l1t::L1TGlobalUtil const & l1tGlobalUtil() const
HLTPrescaleProvider(edm::ParameterSet const &pset, edm::ConsumesCollector &&iC, T &module)
std::unique_ptr< l1t::L1TGlobalUtil > l1tGlobalUtil_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::vector< std::pair< std::string, double > > getL1PrescaleValueInDetail(const edm::Event &iEvent, const edm::EventSetup &iSetup, const std::string &trigger)
HLTConfigProvider const & hltConfigProvider() const
unsigned char count_[5]
list entry
Definition: mps_splice.py:68
T prescaleValue(unsigned int set, const std::string &trigger) const
HLT prescale value in specific prescale set for a specific trigger path.
long double T
L1GtUtils const & l1GtUtils() const
bool init(const edm::Run &iRun, const edm::EventSetup &iSetup, const std::string &processName, bool &changed)
static constexpr int kL1PrescaleDenominator_
Definition: Run.h:45
tuple module
Definition: callgraph.py:69