CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TriggerBxMonitor.cc
Go to the documentation of this file.
1 // C++ headers
2 #include <cstring>
3 #include <iterator>
4 #include <string>
5 
6 // {fmt} headers
7 #include <fmt/printf.h>
8 
9 // CMSSW headers
28 
29 namespace {
30 
31  struct RunBasedHistograms {
32  public:
34  RunBasedHistograms()
35  : // L1T and HLT configuration
36  hltConfig(),
37  // L1T and HLT results
38  tcds_bx_all(nullptr),
39  l1t_bx_all(nullptr),
40  hlt_bx_all(nullptr),
41  tcds_bx(),
42  l1t_bx(),
43  hlt_bx(),
44  tcds_bx_2d(),
45  l1t_bx_2d(),
46  hlt_bx_2d() {}
47 
48  public:
49  // HLT configuration
51 
52  // L1T and HLT results
53  dqm::reco::MonitorElement* tcds_bx_all;
54  dqm::reco::MonitorElement* l1t_bx_all;
55  dqm::reco::MonitorElement* hlt_bx_all;
56  std::vector<dqm::reco::MonitorElement*> tcds_bx;
57  std::vector<dqm::reco::MonitorElement*> l1t_bx;
58  std::vector<dqm::reco::MonitorElement*> hlt_bx;
59  std::vector<dqm::reco::MonitorElement*> tcds_bx_2d;
60  std::vector<dqm::reco::MonitorElement*> l1t_bx_2d;
61  std::vector<dqm::reco::MonitorElement*> hlt_bx_2d;
62  };
63 
64 } // namespace
65 
66 class TriggerBxMonitor : public DQMGlobalEDAnalyzer<RunBasedHistograms> {
67 public:
68  explicit TriggerBxMonitor(edm::ParameterSet const&);
69  ~TriggerBxMonitor() override = default;
70 
71  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
72 
73 private:
74  void dqmBeginRun(edm::Run const&, edm::EventSetup const&, RunBasedHistograms&) const override;
75  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&, RunBasedHistograms&) const override;
76  void dqmAnalyze(edm::Event const&, edm::EventSetup const&, RunBasedHistograms const&) const override;
77 
78  // number of bunch crossings
79  static const unsigned int s_bx_range = 3564;
80 
81  // TCDS trigger types
82  // see https://twiki.cern.ch/twiki/bin/viewauth/CMS/TcdsEventRecord
83  static constexpr const char* s_tcds_trigger_types[] = {
84  "Empty", // 0 - No trigger
85  "Physics", // 1 - GT trigger
86  "Calibration", // 2 - Sequence trigger (calibration)
87  "Random", // 3 - Random trigger
88  "Auxiliary", // 4 - Auxiliary (CPM front panel NIM input) trigger
89  nullptr, // 5 - reserved
90  nullptr, // 6 - reserved
91  nullptr, // 7 - reserved
92  "Cyclic", // 8 - Cyclic trigger
93  "Bunch-pattern", // 9 - Bunch-pattern trigger
94  "Software", // 10 - Software trigger
95  "TTS", // 11 - TTS-sourced trigger
96  nullptr, // 12 - reserved
97  nullptr, // 13 - reserved
98  nullptr, // 14 - reserved
99  nullptr // 15 - reserved
100  };
101 
102  // module configuration
107  const bool m_make_1d_plots;
108  const bool m_make_2d_plots;
109  const uint32_t m_ls_range;
110 };
111 
112 // definition
113 constexpr const char* TriggerBxMonitor::s_tcds_trigger_types[];
114 
117  desc.addUntracked<edm::InputTag>("l1tResults", edm::InputTag("gtStage2Digis"));
118  desc.addUntracked<edm::InputTag>("hltResults", edm::InputTag("TriggerResults"));
119  desc.addUntracked<std::string>("dqmPath", "HLT/TriggerBx");
120  desc.addUntracked<bool>("make1DPlots", true);
121  desc.addUntracked<bool>("make2DPlots", false);
122  desc.addUntracked<uint32_t>("lsRange", 4000);
123  descriptions.add("triggerBxMonitor", desc);
124 }
125 
127  : // module configuration
128  m_l1tMenuToken{esConsumes<edm::Transition::BeginRun>()},
129  m_l1t_results(consumes<GlobalAlgBlkBxCollection>(config.getUntrackedParameter<edm::InputTag>("l1tResults"))),
130  m_hlt_results(consumes<edm::TriggerResults>(config.getUntrackedParameter<edm::InputTag>("hltResults"))),
131  m_dqm_path(config.getUntrackedParameter<std::string>("dqmPath")),
132  m_make_1d_plots(config.getUntrackedParameter<bool>("make1DPlots")),
133  m_make_2d_plots(config.getUntrackedParameter<bool>("make2DPlots")),
134  m_ls_range(config.getUntrackedParameter<uint32_t>("lsRange")) {}
135 
137  edm::EventSetup const& setup,
138  RunBasedHistograms& histograms) const {
139  // initialise the TCDS vector
140  if (m_make_1d_plots) {
141  histograms.tcds_bx.clear();
142  histograms.tcds_bx.resize(std::size(s_tcds_trigger_types));
143  }
144  if (m_make_2d_plots) {
145  histograms.tcds_bx_2d.clear();
146  histograms.tcds_bx_2d.resize(std::size(s_tcds_trigger_types));
147  }
148 
149  // cache the L1 trigger menu
150  if (m_make_1d_plots) {
151  histograms.l1t_bx.clear();
152  histograms.l1t_bx.resize(GlobalAlgBlk::maxPhysicsTriggers);
153  }
154  if (m_make_2d_plots) {
155  histograms.l1t_bx_2d.clear();
156  histograms.l1t_bx_2d.resize(GlobalAlgBlk::maxPhysicsTriggers);
157  }
158 
159  // initialise the HLTConfigProvider
160  bool changed = true;
162  labelsForToken(m_hlt_results, labels);
163  if (histograms.hltConfig.init(run, setup, labels.process, changed)) {
164  if (m_make_1d_plots) {
165  histograms.hlt_bx.clear();
166  histograms.hlt_bx.resize(histograms.hltConfig.size());
167  }
168  if (m_make_2d_plots) {
169  histograms.hlt_bx_2d.clear();
170  histograms.hlt_bx_2d.resize(histograms.hltConfig.size());
171  }
172  } else {
173  // HLTConfigProvider not initialised, skip the the HLT monitoring
174  edm::LogError("TriggerBxMonitor")
175  << "failed to initialise HLTConfigProvider, the HLT bx distribution will not be monitored";
176  }
177 }
178 
180  edm::Run const& run,
181  edm::EventSetup const& setup,
182  RunBasedHistograms& histograms) const {
183  // TCDS trigger type plots
184  {
186 
187  // book 2D histogram to monitor all TCDS trigger types in a single plot
189  histograms.tcds_bx_all = booker.book2D("TCDS Trigger Types",
190  "TCDS Trigger Types vs. bunch crossing",
191  s_bx_range + 1,
192  -0.5,
193  s_bx_range + 0.5,
194  size,
195  -0.5,
196  size - 0.5);
197 
198  // book the individual histograms for the known TCDS trigger types
199  booker.setCurrentFolder(m_dqm_path + "/TCDS");
200  for (unsigned int i = 0; i < size; ++i) {
201  if (s_tcds_trigger_types[i]) {
202  if (m_make_1d_plots) {
203  histograms.tcds_bx.at(i) =
204  booker.book1D(s_tcds_trigger_types[i], s_tcds_trigger_types[i], s_bx_range + 1, -0.5, s_bx_range + 0.5);
205  }
206  if (m_make_2d_plots) {
207  std::string const& name_ls = std::string(s_tcds_trigger_types[i]) + " vs LS";
208  histograms.tcds_bx_2d.at(i) = booker.book2D(
209  name_ls, name_ls, s_bx_range + 1, -0.5, s_bx_range + 0.5, m_ls_range, 0.5, m_ls_range + 0.5);
210  }
211  histograms.tcds_bx_all->setBinLabel(i + 1, s_tcds_trigger_types[i], 2); // Y axis
212  }
213  }
214  }
215 
216  // L1T plots
217  {
218  // book 2D histogram to monitor all L1 triggers in a single plot
220  histograms.l1t_bx_all = booker.book2D("Level 1 Triggers",
221  "Level 1 Triggers vs. bunch crossing",
222  s_bx_range + 1,
223  -0.5,
224  s_bx_range + 0.5,
226  -0.5,
228 
229  // book the individual histograms for the L1 triggers that are included in the L1 menu
230  booker.setCurrentFolder(m_dqm_path + "/L1T");
231  auto const& l1tMenu = setup.getData(m_l1tMenuToken);
232  for (auto const& keyval : l1tMenu.getAlgorithmMap()) {
233  unsigned int bit = keyval.second.getIndex();
234  std::string const& name = fmt::sprintf("%s (bit %d)", keyval.first, bit);
235  if (m_make_1d_plots) {
236  histograms.l1t_bx.at(bit) = booker.book1D(name, name, s_bx_range + 1, -0.5, s_bx_range + 0.5);
237  }
238  if (m_make_2d_plots) {
239  std::string const& name_ls = name + " vs LS";
240  histograms.l1t_bx_2d.at(bit) =
241  booker.book2D(name_ls, name_ls, s_bx_range + 1, -0.5, s_bx_range + 0.5, m_ls_range, 0.5, m_ls_range + 0.5);
242  }
243  histograms.l1t_bx_all->setBinLabel(bit + 1, keyval.first, 2); // Y axis
244  }
245  }
246 
247  // HLT plots
248  if (histograms.hltConfig.inited()) {
249  // book 2D histogram to monitor all HLT paths in a single plot
251  histograms.hlt_bx_all = booker.book2D("High Level Triggers",
252  "High Level Triggers vs. bunch crossing",
253  s_bx_range + 1,
254  -0.5,
255  s_bx_range + 0.5,
256  histograms.hltConfig.size(),
257  -0.5,
258  histograms.hltConfig.size() - 0.5);
259 
260  // book the individual HLT triggers histograms
261  booker.setCurrentFolder(m_dqm_path + "/HLT");
262  for (unsigned int i = 0; i < histograms.hltConfig.size(); ++i) {
263  std::string const& name = histograms.hltConfig.triggerName(i);
264  if (m_make_1d_plots) {
265  histograms.hlt_bx[i] = booker.book1D(name, name, s_bx_range + 1, -0.5, s_bx_range + 0.5);
266  }
267  if (m_make_2d_plots) {
268  std::string const& name_ls = name + " vs LS";
269  histograms.hlt_bx_2d[i] =
270  booker.book2D(name_ls, name_ls, s_bx_range + 1, -0.5, s_bx_range + 0.5, m_ls_range, 0.5, m_ls_range + 0.5);
271  }
272  histograms.hlt_bx_all->setBinLabel(i + 1, name, 2); // Y axis
273  }
274  }
275 }
276 
278  edm::EventSetup const& setup,
279  RunBasedHistograms const& histograms) const {
280  unsigned int bx = event.bunchCrossing();
281  unsigned int ls = event.luminosityBlock();
282 
283  // monitor the bx distribution for the TCDS trigger types
284  {
286  unsigned int type = event.experimentType();
287  if (type < size) {
288  if (m_make_1d_plots and histograms.tcds_bx.at(type))
289  histograms.tcds_bx[type]->Fill(bx);
290  if (m_make_2d_plots and histograms.tcds_bx_2d.at(type))
291  histograms.tcds_bx_2d[type]->Fill(bx, ls);
292  }
293  histograms.tcds_bx_all->Fill(bx, type);
294  }
295 
296  // monitor the bx distribution for the L1 triggers
297  {
298  auto const& bxvector = event.get(m_l1t_results);
299  if (not bxvector.isEmpty(0)) {
300  auto const& results = bxvector.at(0, 0);
301  for (unsigned int i = 0; i < GlobalAlgBlk::maxPhysicsTriggers; ++i)
302  if (results.getAlgoDecisionFinal(i)) {
303  if (m_make_1d_plots and histograms.l1t_bx.at(i))
304  histograms.l1t_bx[i]->Fill(bx);
305  if (m_make_2d_plots and histograms.l1t_bx_2d.at(i))
306  histograms.l1t_bx_2d[i]->Fill(bx, ls);
307  histograms.l1t_bx_all->Fill(bx, i);
308  }
309  }
310  }
311 
312  // monitor the bx distribution for the HLT triggers
313  if (histograms.hltConfig.inited()) {
314  auto const& hltResults = event.get(m_hlt_results);
315  for (unsigned int i = 0; i < hltResults.size(); ++i) {
316  if (hltResults.at(i).accept()) {
317  if (m_make_1d_plots and histograms.hlt_bx.at(i))
318  histograms.hlt_bx[i]->Fill(bx);
319  if (m_make_2d_plots and histograms.hlt_bx_2d.at(i))
320  histograms.hlt_bx_2d[i]->Fill(bx, ls);
321  histograms.hlt_bx_all->Fill(bx, i);
322  }
323  }
324  }
325 }
326 
327 //define this as a plug-in
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void dqmAnalyze(edm::Event const &, edm::EventSetup const &, RunBasedHistograms const &) const override
T getUntrackedParameter(std::string const &, T const &) const
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
~TriggerBxMonitor() override=default
const edm::EDGetTokenT< GlobalAlgBlkBxCollection > m_l1t_results
dictionary results
def ls
Definition: eostools.py:349
void dqmBeginRun(edm::Run const &, edm::EventSetup const &, RunBasedHistograms &) const override
Log< level::Error, false > LogError
const edm::EDGetTokenT< edm::TriggerResults > m_hlt_results
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &, RunBasedHistograms &) const override
char const * process
Definition: ProductLabels.h:7
bool getData(T &iHolder) const
Definition: EventSetup.h:128
static constexpr const char * s_tcds_trigger_types[]
const uint32_t m_ls_range
const std::string m_dqm_path
TriggerBxMonitor(edm::ParameterSet const &)
dqm::legacy::MonitorElement MonitorElement
const bool m_make_1d_plots
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:177
static const unsigned int s_bx_range
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void labelsForToken(EDGetToken iToken, Labels &oLabels) const
tuple config
parse the configuration file
static constexpr unsigned int maxPhysicsTriggers
Definition: GlobalAlgBlk.h:52
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
const bool m_make_2d_plots
tuple size
Write out results.
Definition: Run.h:45
const edm::ESGetToken< L1TUtmTriggerMenu, L1TUtmTriggerMenuRcd > m_l1tMenuToken