CMS 3D CMS Logo

HLTPrescaleProvider.cc
Go to the documentation of this file.
1 
14 
15 #include <cassert>
16 #include <sstream>
17 
18 static const bool useL1EventSetup(true);
19 static const bool useL1GtTriggerMenuLite(false);
20 static const unsigned char countMax(2);
21 
23  const edm::EventSetup& iSetup,
24  const std::string& processName,
25  bool& changed) {
26  inited_ = true;
27 
28  count_[0] = 0;
29  count_[1] = 0;
30  count_[2] = 0;
31  count_[3] = 0;
32  count_[4] = 0;
33 
34  const bool result(hltConfigProvider_.init(iRun, iSetup, processName, changed));
35 
36  const unsigned int l1tType(hltConfigProvider_.l1tType());
37  if (l1tType == 1) {
40  l1GtUtils_->getL1GtRunCache(iRun, iSetup, useL1EventSetup, useL1GtTriggerMenuLite);
41  } else if (l1tType == 2) {
43  l1tGlobalUtil_->retrieveL1Setup(iSetup);
44  } else {
45  edm::LogError("HLTPrescaleProvider") << " Unknown L1T Type " << l1tType << " - prescales will not be avaiable!";
46  }
47 
48  return result;
49 }
50 
53  return *l1GtUtils_;
54 }
55 
58  return *l1tGlobalUtil_;
59 }
60 
62  if (!inited_) {
63  throw cms::Exception("LogicError") << "HLTPrescaleProvider::prescaleSet,\n"
64  "HLTPrescaleProvider::init was not called at beginRun\n";
65  }
66  const unsigned int l1tType(hltConfigProvider_.l1tType());
67  if (l1tType == 1) {
69 
70  // return hltPrescaleTable_.set();
71  l1GtUtils_->getL1GtRunCache(iEvent, iSetup, useL1EventSetup, useL1GtTriggerMenuLite);
72  int errorTech(0);
73  const int psfsiTech(l1GtUtils_->prescaleFactorSetIndex(iEvent, L1GtUtils::TechnicalTrigger, errorTech));
74  int errorPhys(0);
75  const int psfsiPhys(l1GtUtils_->prescaleFactorSetIndex(iEvent, L1GtUtils::AlgorithmTrigger, errorPhys));
76  assert(psfsiTech == psfsiPhys);
77  if ((errorTech == 0) && (errorPhys == 0) && (psfsiTech >= 0) && (psfsiPhys >= 0) && (psfsiTech == psfsiPhys)) {
78  return psfsiPhys;
79  } else {
81  if (count_[0] < countMax) {
82  count_[0] += 1;
83  edm::LogError("HLTPrescaleProvider")
84  << " Using processName '" << hltConfigProvider_.processName() << "':"
85  << " Error in determining HLT prescale set index from L1 data using L1GtUtils:"
86  << " Tech/Phys error = " << errorTech << "/" << errorPhys << " Tech/Phys psfsi = " << psfsiTech << "/"
87  << psfsiPhys;
88  }
89  return -1;
90  }
91  } else if (l1tType == 2) {
93  l1tGlobalUtil_->retrieveL1Event(iEvent, iSetup);
94  return static_cast<int>(l1tGlobalUtil_->prescaleColumn());
95  } else {
96  if (count_[0] < countMax) {
97  count_[0] += 1;
98  edm::LogError("HLTPrescaleProvider")
99  << " Using processName '" << hltConfigProvider_.processName() << "':"
100  << " Unknown L1T Type " << l1tType << " - can not determine prescale set index!";
101  }
102  return -1;
103  }
104 }
105 
106 template <>
108  int numer = static_cast<int>(val * kL1PrescaleDenominator_ + 0.5);
109  static constexpr double kL1RoundingEpsilon = 0.001;
110  if (std::abs(numer - val * kL1PrescaleDenominator_) > kL1RoundingEpsilon) {
111  edm::LogWarning("ValueError") << " Error, L1 prescale val " << val
112  << " does not appear to precisely expressable as int / " << kL1PrescaleDenominator_
113  << ", using a FractionalPrescale is a loss of precision";
114  }
115 
116  return {numer, kL1PrescaleDenominator_};
117 }
118 
120  const edm::EventSetup& iSetup,
121  const std::string& trigger) {
122  // get L1T prescale - works only for those hlt trigger paths with
123  // exactly one L1GT seed module which has exactly one L1T name as seed
124 
125  double result = -1;
126 
127  const unsigned int l1tType(hltConfigProvider_.l1tType());
128  if (l1tType == 1) {
129  checkL1GtUtils();
130  const unsigned int nL1GTSeedModules(hltConfigProvider_.hltL1GTSeeds(trigger).size());
131  if (nL1GTSeedModules == 0) {
132  // no L1 seed module on path hence no L1 seed hence formally no L1 prescale
133  result = 1;
134  } else if (nL1GTSeedModules == 1) {
135  l1GtUtils_->getL1GtRunCache(iEvent, iSetup, useL1EventSetup, useL1GtTriggerMenuLite);
136  const std::string l1tname(hltConfigProvider_.hltL1GTSeeds(trigger).at(0).second);
137  int l1error(0);
138  result = l1GtUtils_->prescaleFactor(iEvent, l1tname, l1error);
139  if (l1error != 0) {
140  if (count_[1] < countMax) {
141  count_[1] += 1;
142  edm::LogError("HLTPrescaleProvider")
143  << " Error in determining L1T prescale for HLT path: '" << trigger << "' with L1T seed: '" << l1tname
144  << "' using L1GtUtils: error code = " << l1error << "." << std::endl
145  << " Note: for this method ('prescaleValues'), only a single L1T name (and not a bit number) is allowed "
146  "as seed!"
147  << std::endl
148  << " For seeds being complex logical expressions, try the new method 'prescaleValuesInDetail'."
149  << std::endl;
150  }
151  result = -1;
152  }
153  } else {
155  if (count_[2] < countMax) {
156  count_[2] += 1;
157  std::string dump("'" + hltConfigProvider_.hltL1GTSeeds(trigger).at(0).second + "'");
158  for (unsigned int i = 1; i != nL1GTSeedModules; ++i) {
159  dump += " * '" + hltConfigProvider_.hltL1GTSeeds(trigger).at(i).second + "'";
160  }
161  edm::LogError("HLTPrescaleProvider")
162  << " Error in determining L1T prescale for HLT path: '" << trigger << "' has multiple L1GTSeed modules, "
163  << nL1GTSeedModules << ", with L1 seeds: " << dump
164  << ". (Note: at most one L1GTSeed module is allowed for a proper determination of the L1T prescale!)";
165  }
166  result = -1;
167  }
168  } else if (l1tType == 2) {
170  const unsigned int nL1TSeedModules(hltConfigProvider_.hltL1TSeeds(trigger).size());
171  if (nL1TSeedModules == 0) {
172  // no L1 seed module on path hence no L1 seed hence formally no L1 prescale
173  result = 1;
174  } else if (nL1TSeedModules == 1) {
175  // l1tGlobalUtil_->retrieveL1Event(iEvent,iSetup);
176  const std::string l1tname(hltConfigProvider_.hltL1TSeeds(trigger).at(0));
177  bool l1error(!l1tGlobalUtil_->getPrescaleByName(l1tname, result));
178  if (l1error) {
179  if (count_[1] < countMax) {
180  count_[1] += 1;
181  edm::LogError("HLTPrescaleProvider")
182  << " Error in determining L1T prescale for HLT path: '" << trigger << "' with L1T seed: '" << l1tname
183  << "' using L1TGlobalUtil: error cond = " << l1error << "." << std::endl
184  << " Note: for this method ('prescaleValues'), only a single L1T name (and not a bit number) is allowed "
185  "as seed!"
186  << std::endl
187  << " For seeds being complex logical expressions, try the new method 'prescaleValuesInDetail'."
188  << std::endl;
189  }
190  result = -1;
191  }
192  } else {
194  if (count_[2] < countMax) {
195  count_[2] += 1;
197  for (unsigned int i = 1; i != nL1TSeedModules; ++i) {
198  dump += " * '" + hltConfigProvider_.hltL1TSeeds(trigger).at(i) + "'";
199  }
200  edm::LogError("HLTPrescaleProvider")
201  << " Error in determining L1T prescale for HLT path: '" << trigger << "' has multiple L1TSeed modules, "
202  << nL1TSeedModules << ", with L1T seeds: " << dump
203  << ". (Note: at most one L1TSeed module is allowed for a proper determination of the L1T prescale!)";
204  }
205  result = -1;
206  }
207  } else {
208  if (count_[1] < countMax) {
209  count_[1] += 1;
210  edm::LogError("HLTPrescaleProvider") << " Unknown L1T Type " << l1tType << " - can not determine L1T prescale! ";
211  }
212  result = -1;
213  }
214 
215  return result;
216 }
217 
218 std::vector<std::pair<std::string, double>> HLTPrescaleProvider::getL1PrescaleValueInDetail(
219  const edm::Event& iEvent, const edm::EventSetup& iSetup, const std::string& trigger) {
220  std::vector<std::pair<std::string, double>> result;
221 
222  const unsigned int l1tType(hltConfigProvider_.l1tType());
223  if (l1tType == 1) {
224  checkL1GtUtils();
225 
226  const unsigned int nL1GTSeedModules(hltConfigProvider_.hltL1GTSeeds(trigger).size());
227  if (nL1GTSeedModules == 0) {
228  // no L1 seed module on path hence no L1 seed hence formally no L1 prescale
229  result.clear();
230  } else if (nL1GTSeedModules == 1) {
231  l1GtUtils_->getL1GtRunCache(iEvent, iSetup, useL1EventSetup, useL1GtTriggerMenuLite);
232  const std::string l1tname(hltConfigProvider_.hltL1GTSeeds(trigger).at(0).second);
234  l1Logical.logicalExpressionRunUpdate(iEvent.getRun(), iSetup, l1tname);
235  const std::vector<std::pair<std::string, int>>& errorCodes(l1Logical.errorCodes(iEvent));
236  auto resultInt = l1Logical.prescaleFactors();
237  result.clear();
238  for (const auto& entry : resultInt) {
239  result.push_back(entry);
240  }
241 
242  int l1error(l1Logical.isValid() ? 0 : 1);
243  for (auto const& errorCode : errorCodes) {
244  l1error += std::abs(errorCode.second);
245  }
246  if (l1error != 0) {
247  if (count_[3] < countMax) {
248  count_[3] += 1;
249  std::ostringstream message;
250  message << " Error in determining L1T prescales for HLT path: '" << trigger << "' with complex L1T seed: '"
251  << l1tname << "' using L1GtUtils: " << std::endl
252  << " isValid=" << l1Logical.isValid() << " l1tname/error/prescale " << errorCodes.size() << std::endl;
253  for (unsigned int i = 0; i < errorCodes.size(); ++i) {
254  message << " " << i << ":" << errorCodes[i].first << "/" << errorCodes[i].second << "/" << result[i].second;
255  }
256  message << ".";
257  edm::LogError("HLTPrescaleProvider") << message.str();
258  }
259  result.clear();
260  }
261  } else {
263  if (count_[4] < countMax) {
264  count_[4] += 1;
265  std::string dump("'" + hltConfigProvider_.hltL1GTSeeds(trigger).at(0).second + "'");
266  for (unsigned int i = 1; i != nL1GTSeedModules; ++i) {
267  dump += " * '" + hltConfigProvider_.hltL1GTSeeds(trigger).at(i).second + "'";
268  }
269  edm::LogError("HLTPrescaleProvider")
270  << " Error in determining L1T prescale for HLT path: '" << trigger << "' has multiple L1GTSeed modules, "
271  << nL1GTSeedModules << ", with L1 seeds: " << dump
272  << ". (Note: at most one L1GTSeed module is allowed for a proper determination of the L1T prescale!)";
273  }
274  result.clear();
275  }
276  } else if (l1tType == 2) {
278  const unsigned int nL1TSeedModules(hltConfigProvider_.hltL1TSeeds(trigger).size());
279  if (nL1TSeedModules == 0) {
280  // no L1 seed module on path hence no L1 seed hence formally no L1 prescale
281  result.clear();
282  } else if (nL1TSeedModules == 1) {
283  // l1tGlobalUtil_->retrieveL1Event(iEvent,iSetup);
285  GlobalLogicParser l1tGlobalLogicParser = GlobalLogicParser(l1tname);
286  const std::vector<GlobalLogicParser::OperandToken> l1tSeeds = l1tGlobalLogicParser.expressionSeedsOperandList();
287  int l1error(0);
288  double l1tPrescale(-1);
289  for (auto const& i : l1tSeeds) {
290  const string& l1tSeed = i.tokenName;
291  if (!l1tGlobalUtil_->getPrescaleByName(l1tSeed, l1tPrescale)) {
292  l1error += 1;
293  }
294  result.push_back(std::pair<std::string, double>(l1tSeed, l1tPrescale));
295  }
296  if (l1error != 0) {
297  if (count_[3] < countMax) {
298  count_[3] += 1;
299  string l1name = l1tname;
300  std::ostringstream message;
301  message << " Error in determining L1T prescales for HLT path: '" << trigger << "' with complex L1T seed: '"
302  << l1tname << "' using L1TGlobalUtil: " << std::endl
303  << " isValid=" << l1tGlobalLogicParser.checkLogicalExpression(l1name) << " l1tname/error/prescale "
304  << l1tSeeds.size() << std::endl;
305  for (unsigned int i = 0; i < l1tSeeds.size(); ++i) {
306  const string& l1tSeed = l1tSeeds[i].tokenName;
307  message << " " << i << ":" << l1tSeed << "/" << l1tGlobalUtil_->getPrescaleByName(l1tSeed, l1tPrescale)
308  << "/" << result[i].second;
309  }
310  message << ".";
311  edm::LogError("HLTPrescaleProvider") << message.str();
312  }
313  result.clear();
314  }
315  } else {
317  if (count_[4] < countMax) {
318  count_[4] += 1;
320  for (unsigned int i = 1; i != nL1TSeedModules; ++i) {
321  dump += " * '" + hltConfigProvider_.hltL1TSeeds(trigger).at(i) + "'";
322  }
323  edm::LogError("HLTPrescaleProvider")
324  << " Error in determining L1T prescale for HLT path: '" << trigger << "' has multiple L1TSeed modules, "
325  << nL1TSeedModules << ", with L1T seeds: " << dump
326  << ". (Note: at most one L1TSeed module is allowed for a proper determination of the L1T prescale!)";
327  }
328  result.clear();
329  }
330  } else {
331  if (count_[3] < countMax) {
332  count_[3] += 1;
333  edm::LogError("HLTPrescaleProvider") << " Unknown L1T Type " << l1tType << " - can not determine L1T prescale! ";
334  }
335  result.clear();
336  }
337 
338  return result;
339 }
340 
342  return hltConfigProvider_.moduleType(hltConfigProvider_.moduleLabel(i, triggerResults.index(i))) == "HLTPrescaler";
343 }
344 
346  if (!l1GtUtils_) {
347  throw cms::Exception("Configuration") << "HLTPrescaleProvider::checkL1GtUtils(),\n"
348  "Attempt to use L1GtUtils object when none was constructed.\n"
349  "Possibly the proper era is not configured or\n"
350  "the module configuration does not use the era properly\n"
351  "or input is from mixed eras";
352  }
353 }
354 
356  if (!l1tGlobalUtil_) {
357  throw cms::Exception("Configuration") << "HLTPrescaleProvider:::checkL1TGlobalUtil(),\n"
358  "Attempt to use L1TGlobalUtil object when none was constructed.\n"
359  "Possibly the proper era is not configured or\n"
360  "the module configuration does not use the era properly\n"
361  "or input is from mixed eras";
362  }
363 }
364 
366  unsigned int stageL1Trigger,
369  bool readPrescalesFromFile) {
370  desc.add<unsigned int>("stageL1Trigger", stageL1Trigger);
373 }
HLTConfigProvider hltConfigProvider_
static void fillPSetDescription(edm::ParameterSetDescription &desc, unsigned int stageL1Trigger, edm::InputTag const &l1tAlgBlkInputTag, edm::InputTag const &l1tExtBlkInputTag, bool readPrescalesFromFile)
std::unique_ptr< L1GtUtils > l1GtUtils_
const std::vector< std::vector< std::string > > & hltL1TSeeds() const
const std::string moduleType(const std::string &module) const
C++ class name of module.
Log< level::Error, false > LogError
assert(be >=bs)
l1t::L1TGlobalUtil const & l1tGlobalUtil() const
double getL1PrescaleValue(const edm::Event &iEvent, const edm::EventSetup &iSetup, const std::string &trigger)
static const unsigned char countMax(2)
unsigned int l1tType() const
L1T type (0=unknown, 1=legacy/stage-1 or 2=stage-2)
boost::rational< int > FractionalPrescale
int iEvent
Definition: GenABIO.cc:224
const std::string & moduleLabel(unsigned int trigger, unsigned int module) const
const std::vector< std::vector< std::pair< bool, std::string > > > & hltL1GTSeeds() const
static const bool useL1GtTriggerMenuLite(false)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static std::string const triggerResults
Definition: EdmProvDump.cc:47
static void fillDescription(edm::ParameterSetDescription &desc, edm::InputTag const &iAlg, edm::InputTag const &iExt, bool readPrescalesFromFile)
int prescaleSet(const edm::Event &iEvent, const edm::EventSetup &iSetup)
bool rejectedByHLTPrescaler(const edm::TriggerResults &triggerResults, unsigned int i) const
static void fillDescription(edm::ParameterSetDescription &desc)
Definition: L1GtUtils.h:137
const std::string & processName() const
process name
bool init(const edm::Run &iRun, const edm::EventSetup &iSetup, const std::string &processName, bool &changed)
d&#39;tor
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)
unsigned char count_[5]
L1GtUtils const & l1GtUtils() const
Log< level::Warning, false > LogWarning
T convertL1PS(double val) const
static const bool useL1EventSetup(true)
bool init(const edm::Run &iRun, const edm::EventSetup &iSetup, const std::string &processName, bool &changed)
static constexpr int kL1PrescaleDenominator_
Definition: Run.h:45