CMS 3D CMS Logo

L1TGlobalPrescaler.cc
Go to the documentation of this file.
1 #include <vector>
2 #include <array>
3 #include <memory>
4 #include <cassert>
5 #include <string>
6 #include <stdexcept>
7 #include <cstring>
8 
9 #include <boost/format.hpp>
10 
11 template <class T, std::size_t N>
12 std::array<T, N> make_array(std::vector<T> const & values) {
13  assert(N == values.size());
14  std::array<T, N> ret;
15  std::copy(values.begin(), values.end(), ret.begin());
16  return ret;
17 }
18 
19 template <class T>
20 bool empty(T const& container) {
21  return container.empty();
22 }
23 
24 bool empty(const char * container) {
25  return container != nullptr;
26 }
27 
28 
29 using namespace std::literals;
30 
31 namespace {
32 
33  template <class T>
34  struct Entry {
35  T value;
36  const char * tag;
37  const char * description;
38  };
39 
40  template <typename S, typename T, unsigned int N>
41  std::string build_comment_from_entries(S pre, const Entry<T> (&entries)[N]) {
42  std::string comment{ pre };
43  size_t length = 0;
44  for (auto entry: entries)
45  if (entry.tag)
46  length = std::max(std::strlen(entry.tag), length);
47  for (auto entry: entries)
48  if (entry.tag) {
49  comment.reserve(comment.size() + length + std::strlen(entry.description) + 8);
50  comment += "\n \"";
51  comment += entry.tag;
52  comment += "\": ";
53  for (unsigned int i = 0; i < length-std::strlen(entry.tag); ++i) comment += ' ';
54  comment += entry.description;
55  }
56  return comment;
57  }
58 
59  template <typename S1, typename S2, typename T, unsigned int N>
60  std::string build_comment_from_entries(S1 pre, const Entry<T> (&entries)[N], S2 post) {
61  std::string comment = build_comment_from_entries(pre, entries);
62  comment += '\n';
63  comment += post;
64  return comment;
65  }
66 
67 
68  template <class T>
69 #if __cplusplus > 201400
70  // extended constexpr support in C++14 and later
71  constexpr
72 #endif
73  T get_enum_value(Entry<T> const * entries, const char * tag) {
74  for (; entries->tag; ++entries)
75  if (std::strcmp(entries->tag, tag) == 0)
76  return entries->value;
77  throw std::logic_error("invalid tag "s + tag);
78  }
79 
80  template <class T>
81 #if __cplusplus > 201400
82  // extended constexpr support in C++14 and later
83  constexpr
84 #endif
85  T get_enum_value(Entry<T> const * entries, const char * tag, T default_value) {
86  for (; entries->tag; ++entries)
87  if (std::strcmp(entries->tag, tag) == 0)
88  return entries->value;
89  return default_value;
90  }
91 
92 } // anonymous
93 
94 // ############################################################################
95 
96 
109 
110 
111 namespace {
112 
113  template <typename T, std::size_t N, typename S>
114  std::array<T, N>
115  getParameterArray(edm::ParameterSet const & config, S const & name)
116  {
117  std::vector<T> values = config.getParameter<std::vector<T>>(name);
118  if (values.size() != N)
120  << "Parameter \"" << name << "\" should have " << N << " elements.\n"
121  << "The number of elements in the configuration is incorrect.";
122  std::array<T, N> ret;
123  std::copy(values.begin(), values.end(), ret.begin());
124  return ret;
125  }
126 
127 } // anonymous
128 
130 public:
132 
133  bool filter(edm::Event& event, edm::EventSetup const& setup) override;
134 
135  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
136 
137 private:
138  enum class Mode {
139  ApplyPrescaleValues, // apply the given prescale values
140  ApplyPrescaleRatios, // apply prescales equal to ratio between the given values and the ones read from the EventSetup
141  ApplyColumnValues, // apply the prescale values from the EventSetup corresponding to the given column index
142  ApplyColumnRatios, // apply prescales equal to ratio between the values corresponsing to the given column index, and the ones read from the EventSetup
143  ForcePrescaleValues, // apply the given prescale values, ignoring the prescales and masks already applied
144  ForceColumnValues, // apply the prescale values from the EventSetup corresponding to the given column index, ignoring the prescales and masks already applied
145  Invalid = -1
146  };
147 
148  static const
149  constexpr Entry<Mode> s_modes[] {
150  { Mode::ApplyPrescaleValues, "applyPrescaleValues", "apply the given prescale values" },
151  { Mode::ApplyPrescaleRatios, "applyPrescaleRatios", "apply prescales equal to ratio between the given values and the ones read from the EventSetup" },
152  { Mode::ApplyColumnValues, "applyColumnValues", "apply the prescale values from the EventSetup corresponding to the given column index" },
153  { Mode::ApplyColumnRatios, "applyColumnRatios", "apply prescales equal to ratio between the values corresponsing to the given column index, and the ones read from the EventSetup" },
154  { Mode::ForcePrescaleValues, "forcePrescaleValues", "apply the given prescale values, ignoring the prescales and masks already applied" },
155  { Mode::ForceColumnValues, "forceColumnValues", "apply the prescale values from the EventSetup corresponding to the given column index, ignoring the prescales and masks already applied" },
156  { Mode::Invalid, nullptr, nullptr }
157  };
158 
159  const Mode m_mode;
161  const std::array<double, GlobalAlgBlk::maxPhysicsTriggers> m_l1tPrescales;
162  std::array<double, GlobalAlgBlk::maxPhysicsTriggers> m_prescales;
163  std::array<unsigned int, GlobalAlgBlk::maxPhysicsTriggers> m_counters;
166 
167 };
168 
169 const
170 constexpr Entry<L1TGlobalPrescaler::Mode> L1TGlobalPrescaler::s_modes[];
171 
172 
174  m_mode( get_enum_value(s_modes, config.getParameter<std::string>("mode").c_str(), Mode::Invalid) ),
175  m_l1tResultsToken( consumes<GlobalAlgBlkBxCollection>(config.getParameter<edm::InputTag>("l1tResults")) ),
176  m_l1tPrescales( m_mode == Mode::ApplyPrescaleValues or m_mode == Mode::ApplyPrescaleRatios or m_mode == Mode::ForcePrescaleValues ?
177  getParameterArray<double, GlobalAlgBlk::maxPhysicsTriggers>(config, "l1tPrescales") :
178  std::array<double, GlobalAlgBlk::maxPhysicsTriggers>() ),
179  m_l1tPrescaleColumn( m_mode == Mode::ApplyColumnValues or m_mode == Mode::ApplyColumnRatios or m_mode == Mode::ForceColumnValues ?
180  config.getParameter<uint32_t>("l1tPrescaleColumn") : 0 ),
181  m_oldIndex(-1)
182 {
183  switch (m_mode) {
184  // if the mode is "applyPrescaleValues", use the given values
188  break;
189 
190  // otherwise we need to read the prescale values from the EventSetup
195  break;
196 
197  // this should never happen
198  case Mode::Invalid:
199  throw edm::Exception(edm::errors::Configuration) << "invalid mode \"" << config.getParameter<std::string>("mode") << "\"";
200  }
201 
202  m_counters.fill(0);
203  produces<GlobalAlgBlkBxCollection>();
204 }
205 
208  event.getByToken(m_l1tResultsToken, handle);
209 
210  // if the input collection does not have any information for bx 0,
211  // produce an empty collection, and fail
212  if (handle->isEmpty(0)) {
213  std::unique_ptr<GlobalAlgBlkBxCollection> result(new GlobalAlgBlkBxCollection());
214  event.put(std::move(result));
215  return false;
216  }
217 
218  // read the prescale index
219  int index = handle->at(0,0).getPreScColumn();
220  assert(index >= 0);
221 
222  // Mode::ApplyPrescaleRatios
223  // apply prescales equal to ratio between the given values and the ones read from the EventSetup
224  if (m_mode == Mode::ApplyPrescaleRatios and m_oldIndex != index) {
226  setup.get<L1TGlobalPrescalesVetosRcd>().get(h);
227  auto const & prescaleTable = h->prescale_table_;
228  if (index >= (int) prescaleTable.size())
229  throw edm::Exception(edm::errors::LogicError) << boost::format("The prescale index %d is invalid, it should be smaller than the prescale table size %d.") % index % prescaleTable.size();
230  auto const & prescales = prescaleTable[index];
231  unsigned long i = 0;
232  for (; i < std::min(prescales.size(), (unsigned long) GlobalAlgBlk::maxPhysicsTriggers); ++i)
233  if (m_l1tPrescales[i] == 0) {
234  // if the trigger is requested to be disabled, just do it
235  m_prescales[i] = 0.;
236  } else if (prescales[i] == 0) {
237  // othersie, if the trigger was originally disabled, warn the user and keep it that way
238  m_prescales[i] = 0.;
239  edm::LogWarning("L1TGlobalPrescaler") << "Request to enable the trigger " << i << " which was originally disabled\nIt will be kept disabled.";
240  } else if (m_l1tPrescales[i] < prescales[i]) {
241  // if the target prescale is lower than the original prescale, keep the trigger unprescaled
242  m_prescales[i] = 1.;
243  edm::LogWarning("L1TGlobalPrescaler") << "Request to prescale the trigger " << i << " less than it was originally prescaled\nNo further prescale will be applied.";
244  } else {
245  // apply the ratio of the new and old prescales
246  m_prescales[i] = (double) m_l1tPrescales[i] / prescales[i];
247  }
248  for (; i < (unsigned long) GlobalAlgBlk::maxPhysicsTriggers; ++i)
249  // disable the triggers not included in the prescale table
250  m_prescales[i] = 0.;
251  // reset the prescales
252  m_counters.fill(0);
253  m_oldIndex = index;
254  }
255 
256  // Mode::ApplyColumnValues and Mode::ForceColumnValues
257  // apply the prescale values from the EventSetup corresponding to the given column index
260  setup.get<L1TGlobalPrescalesVetosRcd>().get(h);
261  auto const & prescaleTable = h->prescale_table_;
262  if (m_l1tPrescaleColumn >= (int) prescaleTable.size())
263  throw edm::Exception(edm::errors::Configuration) << boost::format("The prescale index %d is invalid, it should be smaller than the prescale table size %d.") % m_l1tPrescaleColumn % prescaleTable.size();
264  auto const & targets = prescaleTable[m_l1tPrescaleColumn];
265  unsigned long i = 0;
266  for (; i < std::min(targets.size(), (unsigned long) GlobalAlgBlk::maxPhysicsTriggers); ++i)
267  // read the prescales from the EventSetup
268  m_prescales[i] = targets[i];
269  for (; i < (unsigned long) GlobalAlgBlk::maxPhysicsTriggers; ++i)
270  // disable the triggers not included in the prescale table
271  m_prescales[i] = 0.;
272  // reset the prescales
273  m_counters.fill(0);
275  }
276 
277  // Mode::ApplyColumnRatios
278  // apply prescales equal to ratio between the values corresponsing to the given column index, and the ones read from the EventSetup
279  if (m_mode == Mode::ApplyColumnRatios and m_oldIndex != index) {
281  setup.get<L1TGlobalPrescalesVetosRcd>().get(h);
282  auto const & prescaleTable = h->prescale_table_;
283  if (index >= (int) prescaleTable.size())
284  throw edm::Exception(edm::errors::LogicError) << boost::format("The prescale index %d is invalid, it should be smaller than the prescale table size %d.") % index % prescaleTable.size();
285  if (m_l1tPrescaleColumn >= (int) prescaleTable.size())
286  throw edm::Exception(edm::errors::Configuration) << boost::format("The prescale index %d is invalid, it should be smaller than the prescale table size %d.") % m_l1tPrescaleColumn % prescaleTable.size();
287  auto const & prescales = prescaleTable[index];
288  auto const & targets = prescaleTable[m_l1tPrescaleColumn];
289  unsigned long i = 0;
290  for (; i < std::min({prescales.size(), targets.size(), (unsigned long) GlobalAlgBlk::maxPhysicsTriggers}); ++i)
291  if (prescales[i] == 0)
292  // if the trigger was disabled, keep it disabled
293  m_prescales[i] = 0.;
294  else
295  // if the target prescale is lower than the original prescale, keep the trigger unprescaled
296  m_prescales[i] = targets[i] < prescales[i] ? 1. : (double) targets[i] / prescales[i];
297  for (; i < (unsigned long) GlobalAlgBlk::maxPhysicsTriggers; ++i)
298  // disable the triggers not included in the prescale table
299  m_prescales[i] = 0.;
300  // reset the prescales
301  m_counters.fill(0);
302  m_oldIndex = index;
303  }
304 
305  // make a copy of the GlobalAlgBlk for bx 0
306  GlobalAlgBlk algoBlock = handle->at(0,0);
307 
308  bool finalOr = false;
309  std::vector<bool> const& decision = (m_mode == Mode::ForceColumnValues or m_mode == Mode::ForcePrescaleValues) ?
310  algoBlock.getAlgoDecisionInitial() :
311  algoBlock.getAlgoDecisionFinal();
312 
313  for (unsigned int i = 0; i < GlobalAlgBlk::maxPhysicsTriggers; ++i) {
314  if (m_prescales[i] == 0) {
315  // mask this trigger: reset the bit
316  algoBlock.setAlgoDecisionFinal(i, false);
317  } else if (decision[i]) {
318  // prescale this trigger
319  ++m_counters[i];
320  if (std::fmod(m_counters[i], m_prescales[i]) < 1) {
321  // the prescale is successful, set the bit
322  algoBlock.setAlgoDecisionFinal(i, true);
323  finalOr = true;
324  } else {
325  // the prescale failed, reset the bit
326  algoBlock.setAlgoDecisionFinal(i, false);
327  }
328  }
329  }
330 
331  // set the final OR
332  algoBlock.setFinalORPreVeto(finalOr);
333  if (algoBlock.getFinalORVeto())
334  finalOr = false;
335  algoBlock.setFinalOR(finalOr);
336 
337  // set the new prescale column
340 
341  // create a new GlobalAlgBlkBxCollection, and set the new prescaled decisions for bx 0
342  std::unique_ptr<GlobalAlgBlkBxCollection> result(new GlobalAlgBlkBxCollection());
343  result->push_back(0, algoBlock);
344  event.put(std::move(result));
345 
346  return finalOr;
347 }
348 
349 
351  // collection with the original uGT results
352  edm::ParameterDescription<edm::InputTag> l1tResults("l1tResults", edm::InputTag("gtStage2Digis"), true);
353  l1tResults.setComment("Collection with the original uGT results");
354 
355  // define how to apply the prescale values
356  edm::ParameterDescription<std::string> mode("mode", "applyPrescaleValues", true);
357  mode.setComment(build_comment_from_entries("Define how to apply the prescale values:", s_modes));
358 
359  // target prescale values (for modes "applyPrescaleValues" or "applyPrescaleRatios")
360  edm::ParameterDescription<std::vector<double>> l1tPrescales("l1tPrescales", std::vector<double>(GlobalAlgBlk::maxPhysicsTriggers, 1.), true);
361  l1tPrescales.setComment("Target prescale values (for modes \"applyPrescaleValues\", \"applyPrescaleRatios\" or \"forcePrescaleValues\")");
362 
363  // target prescale column (for modes "applyColumnValues" or "applyColumnRatios")
364  edm::ParameterDescription<uint32_t> l1tPrescaleColumn("l1tPrescaleColumn", 0, true);
365  l1tPrescaleColumn.setComment("Target prescale column (for modes \"applyColumnValues\", \"applyColumnRatios\" or \"forceColumnValues\")");
366 
367  // validaton of all possible configurations and applyPrescaleValues example
368  {
370  desc.addNode(l1tResults);
371  desc.ifValue(mode,
372  // if mode is "applyPrescaleValues", "applyPrescaleRatios" or "forcePrescaleValues", read the target prescales
373  "applyPrescaleValues" >> l1tPrescales or
374  "applyPrescaleRatios" >> l1tPrescales or
375  "forcePrescaleValues" >> l1tPrescales or
376  // if mode is "applyColumnValues", "applyColumnRatios" or "forceColumnValues", read the target column
377  "applyColumnValues" >> l1tPrescaleColumn or
378  "applyColumnRatios" >> l1tPrescaleColumn or
379  "forceColumnValues" >> l1tPrescaleColumn );
380  descriptions.add("l1tGlobalPrescaler", desc);
381  }
382 
383  // applyColumnRatios example
384  {
386  desc.addNode(l1tResults);
387  desc.add<std::string>("mode", "applyColumnRatios")->setComment("apply prescales equal to ratio between the values corresponsing to the given column index, and the ones read from the EventSetup");
388  desc.addNode(l1tPrescaleColumn);
389  descriptions.add("l1tGlobalPrescalerTargetColumn", desc);
390  }
391 }
392 
393 
394 // register as a framework plugin
T getParameter(std::string const &) const
void setComment(std::string const &value)
ParameterDescriptionNode * ifValue(ParameterDescription< T > const &switchParameter, std::unique_ptr< ParameterDescriptionCases< T >> cases)
std::vector< bool > const & getAlgoDecisionInitial() const
Get decision bits.
Definition: GlobalAlgBlk.h:84
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
def copy(args, dbName)
void setFinalORPreVeto(bool fOR)
Definition: GlobalAlgBlk.h:61
std::vector< bool > const & getAlgoDecisionFinal() const
Definition: GlobalAlgBlk.h:90
const std::array< double, GlobalAlgBlk::maxPhysicsTriggers > m_l1tPrescales
std::array< double, GlobalAlgBlk::maxPhysicsTriggers > m_prescales
ParameterDescriptionNode * addNode(ParameterDescriptionNode const &node)
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:2
std::array< unsigned int, GlobalAlgBlk::maxPhysicsTriggers > m_counters
bool isEmpty(int bx) const
BXVector< GlobalAlgBlk > GlobalAlgBlkBxCollection
Definition: GlobalAlgBlk.h:32
Definition: config.py:1
const bool getFinalORVeto() const
Definition: GlobalAlgBlk.h:71
std::array< T, N > make_array(std::vector< T > const &values)
Entry(const ALIstring &type)
Definition: Entry.cc:22
bool filter(edm::Event &event, edm::EventSetup const &setup) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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
bool empty(T const &container)
format
Some error handling for the usage.
T min(T a, T b)
Definition: MathUtil.h:58
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void setFinalOR(bool fOR)
Definition: GlobalAlgBlk.h:62
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const edm::EDGetTokenT< GlobalAlgBlkBxCollection > m_l1tResultsToken
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
static unsigned int maxPhysicsTriggers
Definition: GlobalAlgBlk.h:54
#define N
Definition: blowfish.cc:9
std::vector< std::vector< int > > prescale_table_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
HLT enums.
T get() const
Definition: EventSetup.h:71
void setPreScColumn(int psC)
Definition: GlobalAlgBlk.h:63
long double T
static const Entry< Mode > s_modes[]
def move(src, dest)
Definition: eostools.py:511
#define constexpr
void setAlgoDecisionFinal(unsigned int bit, bool val)
L1TGlobalPrescaler(edm::ParameterSet const &config)
Definition: event.py:1
#define comment(par)
Definition: vmac.h:163
const T & at(int bx, unsigned i) const