CMS 3D CMS Logo

PrescaleWeightProvider.h
Go to the documentation of this file.
1 #ifndef CommonTools_TriggerUtils_PrescaleWeightProvider_h
2 #define CommonTools_TriggerUtils_PrescaleWeightProvider_h
3 
4 // -*- C++ -*-
5 //
6 // Package: CommonTools/TriggerUtils
7 // Class: PrescaleWeightProvider
8 //
9 /*
10  \class PrescaleWeightProvider PrescaleWeightProvider.h "CommonTools/TriggerUtils/interface/PrescaleWeightProvider.h"
11  \brief
12 
13  This class takes a vector of HLT paths and returns a weight based on their
14  HLT and L1 prescales. The weight is equal to the lowest combined (L1*HLT) prescale
15  of the selected paths
16 
17  \author Aram Avetisyan
18 */
19 
20 #include <memory>
21 #include <string>
22 #include <vector>
23 #include <type_traits>
24 
26 
34 
36 
38 
39 namespace edm {
40  class ConsumesCollector;
41  class Event;
42  class EventSetup;
43  class ParameterSet;
44  class Run;
45  class TriggerResults;
46 } // namespace edm
47 
50  bool init_;
51  std::unique_ptr<HLTPrescaleProvider> hltPrescaleProvider_;
53 
54  std::vector<std::string> l1SeedPaths_;
55 
56  // configuration parameters
57  unsigned verbosity_; // optional (default: 0)
58  edm::InputTag triggerResultsTag_; // optional (default: "TriggerResults::HLT")
60  edm::InputTag l1GtTriggerMenuLiteTag_; // optional (default: "l1GtTriggerMenuLite")
61  edm::EDGetTokenT<L1GtTriggerMenuLite> l1GtTriggerMenuLiteToken_; // optional (default: "l1GtTriggerMenuLite")
62  std::vector<std::string> hltPaths_;
63 
64 public:
65  // The constructor must be called from the ED module's c'tor
66  template <typename T>
68 
69  template <typename T>
71 
73 
74  // to be called from the ED module's beginRun() method
75  void initRun(const edm::Run& run, const edm::EventSetup& setup);
76 
77  // to be called from the ED module's event loop method
78  template <typename T = int>
80 
81 private:
83 
84  void parseL1Seeds(const std::string& l1Seeds);
85 };
86 
87 template <typename T>
90 
91 template <typename T>
94  hltPrescaleProvider_ = std::make_unique<HLTPrescaleProvider>(config, iC, module);
95 }
96 
97 template <typename T>
99  static_assert(std::is_same_v<T, double> or std::is_same_v<T, FractionalPrescale>,
100  "\n\tPlease use prescaleWeight<double> or prescaleWeight<FractionalPrescale>"
101  "\n\t(other types for HLT prescales are not supported anymore by PrescaleWeightProvider");
102  if (!init_)
103  return 1;
104 
105  // L1
106  L1GtUtils const& l1GtUtils = hltPrescaleProvider_->l1GtUtils();
107 
108  // HLT
109  HLTConfigProvider const& hltConfig = hltPrescaleProvider_->hltConfigProvider();
110 
112  event.getByToken(triggerResultsToken_, triggerResults);
113  if (!triggerResults.isValid()) {
114  if (verbosity_ > 0)
115  edm::LogError("PrescaleWeightProvider::prescaleWeight")
116  << "TriggerResults product not found for InputTag \"" << triggerResultsTag_.encode() << "\"";
117  return 1;
118  }
119 
120  const int SENTINEL(-1);
121  int weight(SENTINEL);
122 
123  for (unsigned ui = 0; ui < hltPaths_.size(); ui++) {
124  const std::string hltPath(hltPaths_.at(ui));
125  unsigned hltIndex(hltConfig.triggerIndex(hltPath));
126  if (hltIndex == hltConfig.size()) {
127  if (verbosity_ > 0)
128  edm::LogError("PrescaleWeightProvider::prescaleWeight") << "HLT path \"" << hltPath << "\" does not exist";
129  continue;
130  }
131  if (!triggerResults->accept(hltIndex))
132  continue;
133 
134  const std::vector<std::pair<bool, std::string> >& level1Seeds = hltConfig.hltL1GTSeeds(hltPath);
135  if (level1Seeds.size() != 1) {
136  if (verbosity_ > 0)
137  edm::LogError("PrescaleWeightProvider::prescaleWeight")
138  << "HLT path \"" << hltPath << "\" provides too many L1 seeds";
139  return 1;
140  }
141  parseL1Seeds(level1Seeds.at(0).second);
142  if (l1SeedPaths_.empty()) {
143  if (verbosity_ > 0)
144  edm::LogWarning("PrescaleWeightProvider::prescaleWeight")
145  << "Failed to parse L1 seeds for HLT path \"" << hltPath << "\"";
146  continue;
147  }
148 
149  int l1Prescale(SENTINEL);
150  for (unsigned uj = 0; uj < l1SeedPaths_.size(); uj++) {
151  int l1TempPrescale(SENTINEL);
152  int errorCode(0);
153  if (level1Seeds.at(0).first) { // technical triggers
154  unsigned techBit(atoi(l1SeedPaths_.at(uj).c_str()));
155  const std::string techName(*(triggerMenuLite_->gtTechTrigName(techBit, errorCode)));
156  if (errorCode != 0)
157  continue;
158  if (!l1GtUtils.decision(event, techName, errorCode))
159  continue;
160  if (errorCode != 0)
161  continue;
162  l1TempPrescale = l1GtUtils.prescaleFactor(event, techName, errorCode);
163  if (errorCode != 0)
164  continue;
165  } else { // algorithmic triggers
166  if (!l1GtUtils.decision(event, l1SeedPaths_.at(uj), errorCode))
167  continue;
168  if (errorCode != 0)
169  continue;
170  l1TempPrescale = l1GtUtils.prescaleFactor(event, l1SeedPaths_.at(uj), errorCode);
171  if (errorCode != 0)
172  continue;
173  }
174  if (l1TempPrescale > 0) {
175  if (l1Prescale == SENTINEL || l1Prescale > l1TempPrescale)
176  l1Prescale = l1TempPrescale;
177  }
178  }
179  if (l1Prescale == SENTINEL) {
180  if (verbosity_ > 0)
181  edm::LogError("PrescaleWeightProvider::prescaleWeight")
182  << "Unable to find the L1 prescale for HLT path \"" << hltPath << "\"";
183  continue;
184  }
185 
186  auto const prescale = l1Prescale * hltPrescaleProvider_->prescaleValue<T>(event, setup, hltPath);
187 
188  if (prescale > 0) {
189  if (weight == SENTINEL || weight > prescale) {
190  weight = prescale;
191  }
192  }
193  }
194 
195  if (weight == SENTINEL) {
196  if (verbosity_ > 0)
197  edm::LogWarning("PrescaleWeightProvider::prescaleWeight")
198  << "No valid weight for any requested HLT path, returning default weight of 1";
199  return 1;
200  }
201 
202  return weight;
203 }
204 
205 #endif // CommonTools_TriggerUtils_PrescaleWeightProvider_h
PrescaleWeightProvider(const edm::ParameterSet &config, edm::ConsumesCollector &&iC, T &module)
std::string encode() const
Definition: InputTag.cc:159
const int prescaleFactor(const edm::Event &iEvent, const std::string &nameAlgoTechTrig, int &errorCode) const
return prescale factor for a given algorithm or technical trigger
Definition: L1GtUtils.cc:1078
edm::Handle< L1GtTriggerMenuLite > triggerMenuLite_
Definition: weight.py:1
Definition: config.py:1
Log< level::Error, false > LogError
std::vector< std::string > l1SeedPaths_
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
std::vector< std::string > hltPaths_
std::unique_ptr< HLTPrescaleProvider > hltPrescaleProvider_
void initRun(const edm::Run &run, const edm::EventSetup &setup)
static std::string const triggerResults
Definition: EdmProvDump.cc:47
edm::EDGetTokenT< L1GtTriggerMenuLite > l1GtTriggerMenuLiteToken_
T prescaleWeight(const edm::Event &event, const edm::EventSetup &setup)
const std::string * gtTechTrigName(const unsigned int bitNumber, int &errorCode) const
const bool decision(const edm::Event &iEvent, const std::string &nameAlgoTechTrig, int &errorCode) const
Definition: L1GtUtils.cc:1066
HLT enums.
dictionary config
Read in AllInOne config in JSON format.
Definition: DiMuonV_cfg.py:30
edm::EDGetTokenT< edm::TriggerResults > triggerResultsToken_
Log< level::Warning, false > LogWarning
long double T
void parseL1Seeds(const std::string &l1Seeds)
Definition: event.py:1
Definition: Run.h:45