CMS 3D CMS Logo

TriggerBxMonitor.cc
Go to the documentation of this file.
1 // C++ headers
2 #include <cstring>
3 #include <iterator>
4 #include <string>
5 
6 // boost headers
7 #include <boost/format.hpp>
8 
9 // CMSSW headers
32 
33 namespace {
34 
35  struct RunBasedHistograms {
36  public:
37  RunBasedHistograms() :
38  // L1T and HLT configuration
39  hltConfig(),
40  // L1T and HLT results
41  tcds_bx_all(nullptr),
42  l1t_bx_all(nullptr),
43  hlt_bx_all(nullptr),
44  tcds_bx(),
45  l1t_bx(),
46  hlt_bx(),
47  tcds_bx_2d(),
48  l1t_bx_2d(),
49  hlt_bx_2d()
50  {
51  }
52 
53  public:
54  // HLT configuration
55  HLTConfigProvider hltConfig;
56 
57  // L1T and HLT results
58  ConcurrentMonitorElement tcds_bx_all;
59  ConcurrentMonitorElement l1t_bx_all;
60  ConcurrentMonitorElement hlt_bx_all;
61  std::vector<ConcurrentMonitorElement> tcds_bx;
62  std::vector<ConcurrentMonitorElement> l1t_bx;
63  std::vector<ConcurrentMonitorElement> hlt_bx;
64  std::vector<ConcurrentMonitorElement> tcds_bx_2d;
65  std::vector<ConcurrentMonitorElement> l1t_bx_2d;
66  std::vector<ConcurrentMonitorElement> hlt_bx_2d;
67  };
68 
69 }
70 
71 class TriggerBxMonitor : public DQMGlobalEDAnalyzer<RunBasedHistograms> {
72 public:
73  explicit TriggerBxMonitor(edm::ParameterSet const&);
74  ~TriggerBxMonitor() override = default;
75 
76  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
77 
78 private:
79  void dqmBeginRun(edm::Run const&, edm::EventSetup const&, RunBasedHistograms &) const override;
80  void bookHistograms(DQMStore::ConcurrentBooker &, edm::Run const&, edm::EventSetup const&, RunBasedHistograms &) const override;
81  void dqmAnalyze(edm::Event const&, edm::EventSetup const&, RunBasedHistograms const&) const override;
82 
83  // number of bunch crossings
84  static const unsigned int s_bx_range = 3564;
85 
86  // TCDS trigger types
87  // see https://twiki.cern.ch/twiki/bin/viewauth/CMS/TcdsEventRecord
88  static constexpr const char * s_tcds_trigger_types[] = {
89  "Empty", // 0 - No trigger
90  "Physics", // 1 - GT trigger
91  "Calibration", // 2 - Sequence trigger (calibration)
92  "Random", // 3 - Random trigger
93  "Auxiliary", // 4 - Auxiliary (CPM front panel NIM input) trigger
94  nullptr, // 5 - reserved
95  nullptr, // 6 - reserved
96  nullptr, // 7 - reserved
97  "Cyclic", // 8 - Cyclic trigger
98  "Bunch-pattern", // 9 - Bunch-pattern trigger
99  "Software", // 10 - Software trigger
100  "TTS", // 11 - TTS-sourced trigger
101  nullptr, // 12 - reserved
102  nullptr, // 13 - reserved
103  nullptr, // 14 - reserved
104  nullptr // 15 - reserved
105  };
106 
107  // module configuration
111  const bool m_make_1d_plots;
112  const bool m_make_2d_plots;
113  const uint32_t m_ls_range;
114 };
115 
116 // definition
118 
120 {
122  desc.addUntracked<edm::InputTag>( "l1tResults", edm::InputTag("gtStage2Digis"));
123  desc.addUntracked<edm::InputTag>( "hltResults", edm::InputTag("TriggerResults"));
124  desc.addUntracked<std::string>( "dqmPath", "HLT/TriggerBx" );
125  desc.addUntracked<bool>( "make1DPlots", true);
126  desc.addUntracked<bool>( "make2DPlots", false);
127  desc.addUntracked<uint32_t>( "lsRange", 4000);
128  descriptions.add("triggerBxMonitor", desc);
129 }
130 
132  // module configuration
133  m_l1t_results( consumes<GlobalAlgBlkBxCollection>( config.getUntrackedParameter<edm::InputTag>( "l1tResults" ) ) ),
134  m_hlt_results( consumes<edm::TriggerResults>( config.getUntrackedParameter<edm::InputTag>( "hltResults" ) ) ),
135  m_dqm_path( config.getUntrackedParameter<std::string>( "dqmPath" ) ),
136  m_make_1d_plots( config.getUntrackedParameter<bool>( "make1DPlots" ) ),
137  m_make_2d_plots( config.getUntrackedParameter<bool>( "make2DPlots" ) ),
138  m_ls_range( config.getUntrackedParameter<uint32_t>( "lsRange" ) )
139 {
140 }
141 
142 void TriggerBxMonitor::dqmBeginRun(edm::Run const& run, edm::EventSetup const& setup, RunBasedHistograms& histograms) const
143 {
144  // initialise the TCDS vector
145  if (m_make_1d_plots) {
146  histograms.tcds_bx.clear();
147  histograms.tcds_bx.resize(std::size(s_tcds_trigger_types));
148  }
149  if (m_make_2d_plots) {
150  histograms.tcds_bx_2d.clear();
151  histograms.tcds_bx_2d.resize(std::size(s_tcds_trigger_types));
152  }
153 
154  // cache the L1 trigger menu
155  if (m_make_1d_plots) {
156  histograms.l1t_bx.clear();
157  histograms.l1t_bx.resize(GlobalAlgBlk::maxPhysicsTriggers);
158  }
159  if (m_make_2d_plots) {
160  histograms.l1t_bx_2d.clear();
161  histograms.l1t_bx_2d.resize(GlobalAlgBlk::maxPhysicsTriggers);
162  }
163 
164  // initialise the HLTConfigProvider
165  bool changed = true;
167  labelsForToken(m_hlt_results, labels);
168  if (histograms.hltConfig.init(run, setup, labels.process, changed)) {
169  if (m_make_1d_plots) {
170  histograms.hlt_bx.clear();
171  histograms.hlt_bx.resize( histograms.hltConfig.size() );
172  }
173  if (m_make_2d_plots) {
174  histograms.hlt_bx_2d.clear();
175  histograms.hlt_bx_2d.resize( histograms.hltConfig.size() );
176  }
177  } else {
178  // HLTConfigProvider not initialised, skip the the HLT monitoring
179  edm::LogError("TriggerBxMonitor") << "failed to initialise HLTConfigProvider, the HLT bx distribution will not be monitored";
180  }
181 }
182 
184 {
185  // TCDS trigger type plots
186  {
188 
189  // book 2D histogram to monitor all TCDS trigger types in a single plot
191  histograms.tcds_bx_all = booker.book2D("TCDS Trigger Types", "TCDS Trigger Types vs. bunch crossing", s_bx_range + 1, -0.5, s_bx_range + 0.5, size, -0.5, size - 0.5);
192 
193  // book the individual histograms for the known TCDS trigger types
194  booker.setCurrentFolder( m_dqm_path + "/TCDS" );
195  for (unsigned int i = 0; i < size; ++i) {
196  if (s_tcds_trigger_types[i]) {
197  if (m_make_1d_plots) {
198  histograms.tcds_bx.at(i) = booker.book1D(s_tcds_trigger_types[i], s_tcds_trigger_types[i], s_bx_range + 1, -0.5, s_bx_range + 0.5);
199  }
200  if (m_make_2d_plots) {
201  std::string const& name_ls = std::string(s_tcds_trigger_types[i]) + " vs LS";
202  histograms.tcds_bx_2d.at(i) = 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);
203  }
204  histograms.tcds_bx_all.setBinLabel(i+1, s_tcds_trigger_types[i], 2); // Y axis
205  }
206  }
207  }
208 
209  // L1T plots
210  {
211  // book 2D histogram to monitor all L1 triggers in a single plot
212  booker.setCurrentFolder( m_dqm_path );
213  histograms.l1t_bx_all = booker.book2D("Level 1 Triggers", "Level 1 Triggers vs. bunch crossing", s_bx_range + 1, -0.5, s_bx_range + 0.5, GlobalAlgBlk::maxPhysicsTriggers, -0.5, GlobalAlgBlk::maxPhysicsTriggers - 0.5);
214 
215  // book the individual histograms for the L1 triggers that are included in the L1 menu
216  booker.setCurrentFolder( m_dqm_path + "/L1T" );
217  auto const& l1tMenu = edm::get<L1TUtmTriggerMenu, L1TUtmTriggerMenuRcd>(setup);
218  for (auto const& keyval: l1tMenu.getAlgorithmMap()) {
219  unsigned int bit = keyval.second.getIndex();
220  std::string const& name = (boost::format("%s (bit %d)") % keyval.first % bit).str();
221  if (m_make_1d_plots) {
222  histograms.l1t_bx.at(bit) = booker.book1D(name, name, s_bx_range + 1, -0.5, s_bx_range + 0.5);
223  }
224  if (m_make_2d_plots) {
225  std::string const& name_ls = name + " vs LS";
226  histograms.l1t_bx_2d.at(bit) = 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);
227  }
228  histograms.l1t_bx_all.setBinLabel(bit+1, keyval.first, 2); // Y axis
229  }
230  }
231 
232  // HLT plots
233  if (histograms.hltConfig.inited()) {
234  // book 2D histogram to monitor all HLT paths in a single plot
235  booker.setCurrentFolder( m_dqm_path );
236  histograms.hlt_bx_all = booker.book2D("High Level Triggers", "High Level Triggers vs. bunch crossing", s_bx_range + 1, -0.5, s_bx_range + 0.5, histograms.hltConfig.size(), -0.5, histograms.hltConfig.size() - 0.5);
237 
238  // book the individual HLT triggers histograms
239  booker.setCurrentFolder( m_dqm_path + "/HLT" );
240  for (unsigned int i = 0; i < histograms.hltConfig.size(); ++i) {
241  std::string const& name = histograms.hltConfig.triggerName(i);
242  if (m_make_1d_plots) {
243  histograms.hlt_bx[i] = booker.book1D(name, name, s_bx_range + 1, -0.5, s_bx_range + 0.5);
244  }
245  if (m_make_2d_plots) {
246  std::string const& name_ls = name + " vs LS";
247  histograms.hlt_bx_2d[i] = 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);
248  }
249  histograms.hlt_bx_all.setBinLabel(i+1, name, 2); // Y axis
250  }
251  }
252 }
253 
254 
255 void TriggerBxMonitor::dqmAnalyze(edm::Event const& event, edm::EventSetup const& setup, RunBasedHistograms const& histograms) const
256 {
257  unsigned int bx = event.bunchCrossing();
258  unsigned int ls = event.luminosityBlock();
259 
260  // monitor the bx distribution for the TCDS trigger types
261  {
263  unsigned int type = event.experimentType();
264  if (type < size) {
265  if (m_make_1d_plots and histograms.tcds_bx.at(type))
266  histograms.tcds_bx[type].fill(bx);
267  if (m_make_2d_plots and histograms.tcds_bx_2d.at(type))
268  histograms.tcds_bx_2d[type].fill(bx, ls);
269  }
270  histograms.tcds_bx_all.fill(bx, type);
271  }
272 
273  // monitor the bx distribution for the L1 triggers
274  {
275  auto const& bxvector = edm::get(event, m_l1t_results);
276  if (not bxvector.isEmpty(0)) {
277  auto const& results = bxvector.at(0, 0);
278  for (unsigned int i = 0; i < GlobalAlgBlk::maxPhysicsTriggers; ++i)
279  if (results.getAlgoDecisionFinal(i)) {
280  if (m_make_1d_plots and histograms.l1t_bx.at(i))
281  histograms.l1t_bx[i].fill(bx);
282  if (m_make_2d_plots and histograms.l1t_bx_2d.at(i))
283  histograms.l1t_bx_2d[i].fill(bx, ls);
284  histograms.l1t_bx_all.fill(bx, i);
285  }
286  }
287  }
288 
289  // monitor the bx distribution for the HLT triggers
290  if (histograms.hltConfig.inited()) {
291  auto const& hltResults = edm::get(event, m_hlt_results);
292  for (unsigned int i = 0; i < hltResults.size(); ++i) {
293  if (hltResults.at(i).accept()) {
294  if (m_make_1d_plots and histograms.hlt_bx.at(i))
295  histograms.hlt_bx[i].fill(bx);
296  if (m_make_2d_plots and histograms.hlt_bx_2d.at(i))
297  histograms.hlt_bx_2d[i].fill(bx, ls);
298  histograms.hlt_bx_all.fill(bx, i);
299  }
300  }
301  }
302 }
303 
304 
305 //define this as a plug-in
size
Write out results.
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
type
Definition: HCALResponse.h:21
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
void dqmBeginRun(edm::Run const &, edm::EventSetup const &, RunBasedHistograms &) const override
const edm::EDGetTokenT< GlobalAlgBlkBxCollection > m_l1t_results
#define nullptr
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:2
Definition: config.py:1
const edm::EDGetTokenT< edm::TriggerResults > m_hlt_results
char const * process
Definition: ProductLabels.h:7
static const char * s_tcds_trigger_types[]
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
ConcurrentMonitorElement book2D(Args &&...args)
Definition: DQMStore.h:163
example_global void dqmAnalyze(edm::Event const &,@example_global edm::EventSetup const &,@example_global Histograms___class__ const &) const override
ConcurrentMonitorElement book1D(Args &&...args)
Definition: DQMStore.h:160
format
Some error handling for the usage.
const uint32_t m_ls_range
const std::string m_dqm_path
TriggerBxMonitor(edm::ParameterSet const &)
def ls(path, rec=False)
Definition: eostools.py:349
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
static unsigned int maxPhysicsTriggers
Definition: GlobalAlgBlk.h:54
T const & get(Event const &event, InputTag const &tag)(false)
Definition: Event.h:658
const bool m_make_1d_plots
example_stream void bookHistograms(DQMStore::IBooker &,@example_stream edm::Run const &,@example_stream edm::EventSetup const &) override
static const unsigned int s_bx_range
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void labelsForToken(EDGetToken iToken, Labels &oLabels) const
void bookHistograms(DQMStore::ConcurrentBooker &, edm::Run const &, edm::EventSetup const &, RunBasedHistograms &) const override
HLT enums.
#define str(s)
const bool m_make_2d_plots
void dqmAnalyze(edm::Event const &, edm::EventSetup const &, RunBasedHistograms const &) const override
#define constexpr
Definition: event.py:1
Definition: Run.h:45