CMS 3D CMS Logo

HLTPrescaleProvider.h
Go to the documentation of this file.
1 #ifndef HLTrigger_HLTcore_HLTPrescaleProvider_h
2 #define HLTrigger_HLTcore_HLTPrescaleProvider_h
3 
23 
24 #include <memory>
25 #include <string>
26 #include <utility>
27 #include <vector>
28 #include <type_traits>
29 
30 namespace edm {
31  class ConsumesCollector;
32  class Event;
33  class EventSetup;
34  class ParameterSet;
35  class Run;
37 } // namespace edm
38 
40 public:
41  template <typename T>
43 
44  template <typename T>
46 
51  bool init(const edm::Run& iRun, const edm::EventSetup& iSetup, const std::string& processName, bool& changed);
52 
54  L1GtUtils const& l1GtUtils() const;
55  l1t::L1TGlobalUtil const& l1tGlobalUtil() const;
56 
59  int prescaleSet(const edm::Event& iEvent, const edm::EventSetup& iSetup);
60  // negative == error
61 
63  template <typename T = unsigned int>
65  const int set(prescaleSet(iEvent, iSetup));
66  //there is a template specialisation for unsigned in which returns +1 which
67  //emulates old behaviour
68  return set < 0 ? -1 : hltConfigProvider_.prescaleValue<T>(static_cast<unsigned int>(set), trigger);
69  }
70 
72  template <typename TL1 = int, typename THLT = TL1>
73  std::pair<TL1, THLT> prescaleValues(const edm::Event& iEvent,
74  const edm::EventSetup& iSetup,
75  const std::string& trigger) {
76  return {convertL1PS<TL1>(getL1PrescaleValue(iEvent, iSetup, trigger)),
77  prescaleValue<THLT>(iEvent, iSetup, trigger)};
78  }
79  // any one negative => error in retrieving this (L1T or HLT) prescale
80 
81  // In case of a complex Boolean expression as L1 seed
82  template <typename TL1 = int, typename THLT = TL1>
83  std::pair<std::vector<std::pair<std::string, TL1>>, THLT> prescaleValuesInDetail(const edm::Event& iEvent,
84  const edm::EventSetup& iSetup,
85  const std::string& trigger) {
86  std::pair<std::vector<std::pair<std::string, TL1>>, THLT> retval;
87  for (auto& entry : getL1PrescaleValueInDetail(iEvent, iSetup, trigger)) {
88  retval.first.emplace_back(std::move(entry.first), convertL1PS<TL1>(entry.second));
89  }
90  retval.second = prescaleValue<THLT>(iEvent, iSetup, trigger);
91  return retval;
92  }
93  // Event rejected by HLTPrescaler on ith HLT path?
94  bool rejectedByHLTPrescaler(const edm::TriggerResults& triggerResults, unsigned int i) const;
96 
98  unsigned int stageL1Trigger,
101  bool readPrescalesFromFile);
102 
103 private:
104  static constexpr const char* l1tGlobalDecisionKeyword_ = "L1GlobalDecision";
105 
106  void checkL1GtUtils() const;
107  void checkL1TGlobalUtil() const;
108 
109  template <typename T>
110  T convertL1PS(double val) const {
111  static_assert(std::is_same_v<T, double> or std::is_same_v<T, FractionalPrescale>,
112  "\n\n\tPlease use convertL1PS<double> or convertL1PS<FractionalPrescale>"
113  " (other types for L1T prescales are not supported anymore by HLTPrescaleProvider)"
114  "\n\tconvertL1PS is used inside prescaleValues and prescaleValuesInDetail,"
115  " so it might be necessary to specify template arguments for those calls,"
116  "\n\te.g. prescaleValues<double, FractionalPrescale>"
117  " (the 1st argument applies to L1T prescales, the 2nd to HLT prescales)\n");
118  return T(val);
119  }
120 
121  double getL1PrescaleValue(const edm::Event& iEvent, const edm::EventSetup& iSetup, const std::string& trigger);
122 
123  std::vector<std::pair<std::string, double>> getL1PrescaleValueInDetail(const edm::Event& iEvent,
124  const edm::EventSetup& iSetup,
125  const std::string& trigger);
126 
128 
130 
131  std::unique_ptr<L1GtUtils> l1GtUtils_;
132  std::unique_ptr<l1t::L1TGlobalUtil> l1tGlobalUtil_;
133 
134  unsigned char count_[5] = {0, 0, 0, 0, 0};
135  bool inited_ = false;
136 };
137 
138 template <typename T>
140  : HLTPrescaleProvider(pset, iC, module) {}
141 
142 template <typename T>
144  unsigned int stageL1Trigger = pset.getParameter<unsigned int>("stageL1Trigger");
145  if (stageL1Trigger <= 1) {
146  l1GtUtils_ = std::make_unique<L1GtUtils>(pset, iC, false, module, L1GtUtils::UseEventSetupIn::Run);
147  } else {
148  l1tGlobalUtil_ = std::make_unique<l1t::L1TGlobalUtil>(pset, iC, module, l1t::UseEventSetupIn::Run);
149  }
150 }
151 
152 template <>
154 
155 #endif // HLTrigger_HLTcore_HLTPrescaleProvider_h
HLTConfigProvider hltConfigProvider_
static void fillPSetDescription(edm::ParameterSetDescription &desc, unsigned int stageL1Trigger, edm::InputTag const &l1tAlgBlkInputTag, edm::InputTag const &l1tExtBlkInputTag, bool readPrescalesFromFile)
static int l1PrescaleDenominator()
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)
HLTConfigProvider const & hltConfigProvider() const
l1t::L1TGlobalUtil const & l1tGlobalUtil() const
double getL1PrescaleValue(const edm::Event &iEvent, const edm::EventSetup &iSetup, const std::string &trigger)
T prescaleValue(unsigned int set, const std::string &trigger) const
HLT prescale value in specific prescale set for a specific trigger path.
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.
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
static constexpr const char * l1tGlobalDecisionKeyword_
static std::string const triggerResults
Definition: EdmProvDump.cc:47
int prescaleSet(const edm::Event &iEvent, const edm::EventSetup &iSetup)
HLTPrescaleProvider(edm::ParameterSet const &pset, edm::ConsumesCollector &&iC, T &module)
bool rejectedByHLTPrescaler(const edm::TriggerResults &triggerResults, unsigned int i) const
std::unique_ptr< l1t::L1TGlobalUtil > l1tGlobalUtil_
std::vector< std::pair< std::string, double > > getL1PrescaleValueInDetail(const edm::Event &iEvent, const edm::EventSetup &iSetup, const std::string &trigger)
HLT enums.
unsigned char count_[5]
L1GtUtils const & l1GtUtils() const
T convertL1PS(double val) const
long double T
bool init(const edm::Run &iRun, const edm::EventSetup &iSetup, const std::string &processName, bool &changed)
def move(src, dest)
Definition: eostools.py:511
static constexpr int kL1PrescaleDenominator_
Definition: Run.h:45