CMS 3D CMS Logo

TriggerBxMonitor.cc
Go to the documentation of this file.
1 // C++ headers
2 #include <string>
3 #include <cstring>
4 
5 // boost headers
6 #include <boost/regex.hpp>
7 #include <boost/format.hpp>
8 
9 // Root headers
10 #include <TH1F.h>
11 
12 // CMSSW headers
35 
36 // helper functions
37 template <typename T>
38 static
39 const T & get(const edm::Event & event, const edm::EDGetTokenT<T> & token) {
41  event.getByToken(token, handle);
42  if (not handle.isValid())
43  throw * handle.whyFailed();
44  return * handle.product();
45 }
46 
47 template <typename R, typename T>
48 static
49 const T & get(const edm::EventSetup & setup) {
51  setup.get<R>().get(handle);
52  return * handle.product();
53 }
54 
55 
57 public:
58  explicit TriggerBxMonitor(edm::ParameterSet const &);
59  ~TriggerBxMonitor() override;
60 
61  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
62 
63 private:
64  void dqmBeginRun(edm::Run const &, edm::EventSetup const &) override;
65  void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
66  void analyze(edm::Event const &, edm::EventSetup const &) override;
67 
68  // number of bunch crossings
69  static const unsigned int s_bx_range = 3564;
70 
71  // TCDS trigger types
72  // see https://twiki.cern.ch/twiki/bin/viewauth/CMS/TcdsEventRecord
73  static constexpr const char * const s_tcds_trigger_types[] = {
74  "Empty", // 0 - No trigger
75  "Physics", // 1 - GT trigger
76  "Calibration", // 2 - Sequence trigger (calibration)
77  "Random", // 3 - Random trigger
78  "Auxiliary", // 4 - Auxiliary (CPM front panel NIM input) trigger
79  nullptr, // 5 - reserved
80  nullptr, // 6 - reserved
81  nullptr, // 7 - reserved
82  "Cyclic", // 8 - Cyclic trigger
83  "Bunch-pattern", // 9 - Bunch-pattern trigger
84  "Software", // 10 - Software trigger
85  "TTS", // 11 - TTS-sourced trigger
86  nullptr, // 12 - reserved
87  nullptr, // 13 - reserved
88  nullptr, // 14 - reserved
89  nullptr // 15 - reserved
90  };
91 
92  // module configuration
96  const bool m_make_1d_plots;
97  const bool m_make_2d_plots;
98  const uint32_t m_ls_range;
99 
100  // L1T and HLT configuration
103 
104  // L1T and HLT results
106  TH2F * m_l1t_bx_all;
107  TH2F * m_hlt_bx_all;
108  std::vector<TH1F *> m_tcds_bx;
109  std::vector<TH1F *> m_l1t_bx;
110  std::vector<TH1F *> m_hlt_bx;
111  std::vector<TH2F *> m_tcds_bx_2d;
112  std::vector<TH2F *> m_l1t_bx_2d;
113  std::vector<TH2F *> m_hlt_bx_2d;
114 };
115 
116 // definition
118 
119 
121 {
123  desc.addUntracked<edm::InputTag>( "l1tResults", edm::InputTag("gtStage2Digis"));
124  desc.addUntracked<edm::InputTag>( "hltResults", edm::InputTag("TriggerResults"));
125  desc.addUntracked<std::string>( "dqmPath", "HLT/TriggerBx" );
126  desc.addUntracked<bool>( "make1DPlots", true);
127  desc.addUntracked<bool>( "make2DPlots", false);
128  desc.addUntracked<uint32_t>( "lsRange", 4000);
129  descriptions.add("triggerBxMonitor", desc);
130 }
131 
132 
134  // module configuration
135  m_l1t_results( consumes<GlobalAlgBlkBxCollection>( config.getUntrackedParameter<edm::InputTag>( "l1tResults" ) ) ),
136  m_hlt_results( consumes<edm::TriggerResults>( config.getUntrackedParameter<edm::InputTag>( "hltResults" ) ) ),
137  m_dqm_path( config.getUntrackedParameter<std::string>( "dqmPath" ) ),
138  m_make_1d_plots( config.getUntrackedParameter<bool>( "make1DPlots" ) ),
139  m_make_2d_plots( config.getUntrackedParameter<bool>( "make2DPlots" ) ),
140  m_ls_range( config.getUntrackedParameter<uint32_t>( "lsRange" ) ),
141  // L1T and HLT configuration
143  m_hltConfig(),
144  // L1T and HLT results
148  m_tcds_bx(),
149  m_l1t_bx(),
150  m_hlt_bx(),
151  m_tcds_bx_2d(),
152  m_l1t_bx_2d(),
153  m_hlt_bx_2d()
154 {
155 }
156 
158 {
159 }
160 
162 {
163  // initialise the TCDS vector
164  if (m_make_1d_plots) {
165  m_tcds_bx.clear();
166  m_tcds_bx.resize(sizeof(s_tcds_trigger_types) / sizeof(const char *), nullptr);
167  }
168  if (m_make_2d_plots) {
169  m_tcds_bx_2d.clear();
170  m_tcds_bx_2d.resize(sizeof(s_tcds_trigger_types) / sizeof(const char *), nullptr);
171  }
172 
173  // cache the L1 trigger menu
174  m_l1tMenu = & get<L1TUtmTriggerMenuRcd, L1TUtmTriggerMenu>(setup);
175  if (m_l1tMenu) {
176  if (m_make_1d_plots) {
177  m_l1t_bx.clear();
179  }
180  if (m_make_2d_plots) {
181  m_l1t_bx_2d.clear();
183  }
184  } else {
185  edm::LogError("TriggerBxMonitor") << "failed to read the L1 menu from the EventSetup, the L1 trigger bx distribution will not be monitored";
186  }
187 
188  // initialise the HLTConfigProvider
189  bool changed = true;
191  labelsForToken(m_hlt_results, labels);
192  if (m_hltConfig.init(run, setup, labels.process, changed)) {
193  if (m_make_1d_plots) {
194  m_hlt_bx.clear();
195  m_hlt_bx.resize( m_hltConfig.size(), nullptr );
196  }
197  if (m_make_2d_plots) {
198  m_hlt_bx_2d.clear();
199  m_hlt_bx_2d.resize( m_hltConfig.size(), nullptr );
200  }
201  } else {
202  // HLTConfigProvider not initialised, skip the the HLT monitoring
203  edm::LogError("TriggerBxMonitor") << "failed to initialise HLTConfigProvider, the HLT bx distribution will not be monitored";
204  }
205 }
206 
208 {
209  // TCDS trigger type plots
210  {
211  size_t size = sizeof(s_tcds_trigger_types) / sizeof(const char *);
212 
213  // book 2D histogram to monitor all TCDS trigger types in a single plot
214  booker.setCurrentFolder( m_dqm_path );
215  m_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)->getTH2F();
216 
217  // book the individual histograms for the known TCDS trigger types
218  booker.setCurrentFolder( m_dqm_path + "/TCDS" );
219  for (unsigned int i = 0; i < size; ++i) {
220  if (s_tcds_trigger_types[i]) {
221  if (m_make_1d_plots) {
222  m_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)->getTH1F();
223  }
224  if (m_make_2d_plots) {
225  std::string const & name_ls = std::string(s_tcds_trigger_types[i]) + " vs LS";
226  m_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)->getTH2F();
227  }
228  m_tcds_bx_all->GetYaxis()->SetBinLabel(i+1, s_tcds_trigger_types[i]);
229  }
230  }
231  }
232 
233  // L1T plots
234  if (m_l1tMenu) {
235  // book 2D histogram to monitor all L1 triggers in a single plot
236  booker.setCurrentFolder( m_dqm_path );
237  m_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)->getTH2F();
238 
239  // book the individual histograms for the L1 triggers that are included in the L1 menu
240  booker.setCurrentFolder( m_dqm_path + "/L1T" );
241  for (auto const & keyval: m_l1tMenu->getAlgorithmMap()) {
242  unsigned int bit = keyval.second.getIndex();
243  std::string const & name = (boost::format("%s (bit %d)") % keyval.first % bit).str();
244  if (m_make_1d_plots) {
245  m_l1t_bx.at(bit) = booker.book1D(name, name, s_bx_range + 1, -0.5, s_bx_range + 0.5)->getTH1F();
246  }
247  if (m_make_2d_plots) {
248  std::string const & name_ls = name + " vs LS";
249  m_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)->getTH2F();
250  }
251  m_l1t_bx_all->GetYaxis()->SetBinLabel(bit+1, keyval.first.c_str());
252  }
253  }
254 
255  // HLT plots
256  if (m_hltConfig.inited()) {
257  // book 2D histogram to monitor all HLT paths in a single plot
258  booker.setCurrentFolder( m_dqm_path );
259  m_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, m_hltConfig.size(), -0.5, m_hltConfig.size() - 0.5)->getTH2F();
260 
261  // book the individual HLT triggers histograms
262  booker.setCurrentFolder( m_dqm_path + "/HLT" );
263  for (unsigned int i = 0; i < m_hltConfig.size(); ++i) {
265  if (m_make_1d_plots) {
266  m_hlt_bx[i] = booker.book1D(name, name, s_bx_range + 1, -0.5, s_bx_range + 0.5)->getTH1F();
267  }
268  if (m_make_2d_plots) {
269  std::string const & name_ls = name + " vs LS";
270  m_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)->getTH2F();
271  }
272  m_hlt_bx_all->GetYaxis()->SetBinLabel(i+1, name.c_str());
273  }
274  }
275 }
276 
277 
279 {
280  unsigned int bx = event.bunchCrossing();
281  unsigned int ls = event.luminosityBlock();
282 
283  // monitor the bx distribution for the TCDS trigger types
284  {
285  size_t size = sizeof(s_tcds_trigger_types) / sizeof(const char *);
286  unsigned int type = event.experimentType();
287  if (type < size) {
288  if (m_make_1d_plots and m_tcds_bx.at(type))
289  m_tcds_bx[type]->Fill(bx);
290  if (m_make_2d_plots and m_tcds_bx_2d.at(type))
291  m_tcds_bx_2d[type]->Fill(bx, ls);
292  }
293  m_tcds_bx_all->Fill(bx, type);
294  }
295 
296  // monitor the bx distribution for the L1 triggers
297  if (m_l1tMenu) {
298  auto const & bxvector = get<GlobalAlgBlkBxCollection>(event, 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 m_l1t_bx.at(i))
304  m_l1t_bx[i]->Fill(bx);
305  if (m_make_2d_plots and m_l1t_bx_2d.at(i))
306  m_l1t_bx_2d[i]->Fill(bx, ls);
307  m_l1t_bx_all->Fill(bx, i);
308  }
309  }
310  }
311 
312  // monitor the bx distribution for the HLT triggers
313  if (m_hltConfig.inited()) {
314  auto const & hltResults = get<edm::TriggerResults>(event, 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 m_hlt_bx.at(i))
318  m_hlt_bx[i]->Fill(bx);
319  if (m_make_2d_plots and m_hlt_bx_2d.at(i))
320  m_hlt_bx_2d[i]->Fill(bx, ls);
321  m_hlt_bx_all->Fill(bx, i);
322  }
323  }
324  }
325 }
326 
327 
328 //define this as a plug-in
unsigned int size() const
number of trigger paths in trigger table
size
Write out results.
static const char *const s_tcds_trigger_types[]
HLTConfigProvider m_hltConfig
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
type
Definition: HCALResponse.h:21
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
const std::string & triggerName(unsigned int triggerIndex) const
std::vector< TH2F * > m_hlt_bx_2d
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
const edm::EDGetTokenT< GlobalAlgBlkBxCollection > m_l1t_results
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
Definition: config.py:1
const edm::EDGetTokenT< edm::TriggerResults > m_hlt_results
#define nullptr
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
#define constexpr
bool inited() const
Accessors (const methods)
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
char const * process
Definition: ProductLabels.h:7
~TriggerBxMonitor() override
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
std::vector< TH1F * > m_tcds_bx
std::vector< TH1F * > m_l1t_bx
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
format
Some error handling for the usage.
const uint32_t m_ls_range
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
bool isValid() const
Definition: HandleBase.h:74
const std::string m_dqm_path
TriggerBxMonitor(edm::ParameterSet const &)
std::vector< TH2F * > m_tcds_bx_2d
void dqmBeginRun(edm::Run const &, edm::EventSetup const &) override
L1TUtmTriggerMenu const * m_l1tMenu
def ls(path, rec=False)
Definition: eostools.py:348
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
static unsigned int maxPhysicsTriggers
Definition: GlobalAlgBlk.h:52
T const * product() const
Definition: Handle.h:81
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
bool init(const edm::Run &iRun, const edm::EventSetup &iSetup, const std::string &processName, bool &changed)
d&#39;tor
const bool m_make_1d_plots
static const unsigned int s_bx_range
TH1F * getTH1F(void) const
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< TH1F * > m_hlt_bx
const std::map< std::string, L1TUtmAlgorithm > & getAlgorithmMap() const
void labelsForToken(EDGetToken iToken, Labels &oLabels) const
HLT enums.
std::shared_ptr< cms::Exception > whyFailed() const
Definition: HandleBase.h:106
TH2F * getTH2F(void) const
long double T
const bool m_make_2d_plots
void analyze(edm::Event const &, edm::EventSetup const &) override
T const * product() const
Definition: ESHandle.h:86
Definition: event.py:1
Definition: Run.h:43
std::vector< TH2F * > m_l1t_bx_2d