CMS 3D CMS Logo

TriggerBxVsOrbitMonitor.cc
Go to the documentation of this file.
1 // C++ headers
2 #include <cstring>
3 #include <string>
4 
5 // CMSSW headers
27 
28 namespace {
31 
32  struct RunBasedHistograms {
33  dqm::reco::MonitorElement* orbit_bx_all;
34  std::vector<dqm::reco::MonitorElement*> orbit_bx;
35  std::vector<dqm::reco::MonitorElement*> orbit_bx_all_byLS;
36  };
37 } // namespace
38 
39 class TriggerBxVsOrbitMonitor : public DQMGlobalEDAnalyzer<RunBasedHistograms> {
40 public:
42  ~TriggerBxVsOrbitMonitor() override = default;
43 
44  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
45 
46 private:
47  void dqmBeginRun(edm::Run const&, edm::EventSetup const&, RunBasedHistograms&) const override;
48  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&, RunBasedHistograms&) const override;
49  void dqmAnalyze(edm::Event const&, edm::EventSetup const&, RunBasedHistograms const&) const override;
50 
51  // number of bunch crossings
52  static const unsigned int s_bx_range = 3564;
53  static const unsigned int s_orbit_range = 262144; // 2**18 orbits in 1 LS
54 
55  // TCDS trigger types
56  // see https://twiki.cern.ch/twiki/bin/viewauth/CMS/TcdsEventRecord
57  static constexpr const char* const s_tcds_trigger_types[] = {
58  "Empty", // 0 - No trigger
59  "Physics", // 1 - GT trigger
60  "Calibration", // 2 - Sequence trigger (calibration)
61  "Random", // 3 - Random trigger
62  "Auxiliary", // 4 - Auxiliary (CPM front panel NIM input) trigger
63  nullptr, // 5 - reserved
64  nullptr, // 6 - reserved
65  nullptr, // 7 - reserved
66  "Cyclic", // 8 - Cyclic trigger
67  "Bunch-pattern", // 9 - Bunch-pattern trigger
68  "Software", // 10 - Software trigger
69  "TTS", // 11 - TTS-sourced trigger
70  nullptr, // 12 - reserved
71  nullptr, // 13 - reserved
72  nullptr, // 14 - reserved
73  nullptr // 15 - reserved
74  };
75 
76  // module configuration
80  const int m_minLS;
81  const int m_maxLS;
82  const int m_minBX;
83  const int m_maxBX;
84 };
85 
86 // definition
87 constexpr const char* const TriggerBxVsOrbitMonitor::s_tcds_trigger_types[];
88 
91  desc.addUntracked<edm::InputTag>("l1tResults", edm::InputTag("gtStage2Digis"));
92  desc.addUntracked<edm::InputTag>("hltResults", edm::InputTag("TriggerResults"));
93  desc.addUntracked<std::string>("dqmPath", "HLT/TriggerBx");
94  desc.addUntracked<int>("minLS", 134);
95  desc.addUntracked<int>("maxLS", 136);
96  desc.addUntracked<int>("minBX", 894);
97  desc.addUntracked<int>("maxBX", 912);
98  descriptions.add("triggerBxVsOrbitMonitor", desc);
99 }
100 
102  : // module configuration
103  m_l1t_results(consumes<GlobalAlgBlkBxCollection>(config.getUntrackedParameter<edm::InputTag>("l1tResults"))),
104  m_hlt_results(consumes<edm::TriggerResults>(config.getUntrackedParameter<edm::InputTag>("hltResults"))),
105  m_dqm_path(config.getUntrackedParameter<std::string>("dqmPath")),
106  m_minLS(config.getUntrackedParameter<int>("minLS")),
107  m_maxLS(config.getUntrackedParameter<int>("maxLS")),
108  m_minBX(config.getUntrackedParameter<int>("minBX")),
109  m_maxBX(config.getUntrackedParameter<int>("maxBX")) {}
110 
112  edm::EventSetup const& setup,
113  RunBasedHistograms& histograms) const {
114  size_t nLS = m_maxLS - m_minLS + 1;
115 
116  histograms.orbit_bx_all_byLS.clear();
117  histograms.orbit_bx_all_byLS.resize(nLS);
118 
119  histograms.orbit_bx.clear();
120  histograms.orbit_bx.resize(std::size(s_tcds_trigger_types));
121 }
122 
124  edm::Run const& run,
125  edm::EventSetup const& setup,
126  RunBasedHistograms& histograms) const {
127  // TCDS trigger type plots
128  size_t size = std::size(s_tcds_trigger_types);
129  size_t nLS = m_maxLS - m_minLS + 1;
130  size_t nBX = m_maxBX - m_minBX + 1;
131 
132  // book 2D histogram to monitor all TCDS trigger types in a single plot
133  booker.setCurrentFolder(m_dqm_path + "/orbitVsBX");
134  histograms.orbit_bx_all = booker.book2D("OrbitVsBX",
135  "Event orbits vs. bunch crossing",
136  nBX,
137  m_minBX - 0.5,
138  m_maxBX + 0.5,
139  s_orbit_range + 1,
140  -0.5,
141  s_orbit_range + 0.5);
142  histograms.orbit_bx_all->setXTitle("BX");
143  histograms.orbit_bx_all->setYTitle("orbit");
144 
145  for (unsigned int i = 0; i < nLS; ++i) {
146  std::string iname = std::to_string(i);
147  histograms.orbit_bx_all_byLS[i] = booker.book2D("OrbitVsBX_LS" + iname,
148  "Event orbits vs. bunch crossing, for lumisection " + iname,
149  nBX,
150  m_minBX - 0.5,
151  m_maxBX + 0.5,
152  s_orbit_range + 1,
153  -0.5,
154  s_orbit_range + 0.5);
155  histograms.orbit_bx_all_byLS[i]->setXTitle("BX");
156  histograms.orbit_bx_all_byLS[i]->setYTitle("orbit");
157  }
158 
159  booker.setCurrentFolder(m_dqm_path + "/orbitVsBX/TCDS");
160  for (unsigned int i = 0; i < size; ++i) {
161  if (s_tcds_trigger_types[i]) {
162  histograms.orbit_bx[i] = booker.book2D("OrbitVsBX_" + std::string(s_tcds_trigger_types[i]),
163  "Event orbits vs. bunch crossing " + std::string(s_tcds_trigger_types[i]),
164  nBX,
165  m_minBX - 0.5,
166  m_maxBX + 0.5,
167  s_orbit_range + 1,
168  -0.5,
169  s_orbit_range + 0.5);
170  histograms.orbit_bx[i]->setXTitle("BX");
171  histograms.orbit_bx[i]->setYTitle("orbit");
172  }
173  }
174 }
175 
177  edm::EventSetup const& setup,
178  RunBasedHistograms const& histograms) const {
179  unsigned int type = event.experimentType();
180  unsigned int ls = event.id().luminosityBlock();
181  unsigned int orbit = event.orbitNumber() % s_orbit_range;
182  unsigned int bx = event.bunchCrossing();
183  histograms.orbit_bx_all->Fill(bx, orbit);
184 
185  int iLS = ls - m_minLS;
186  if (iLS >= 0 and iLS < int(histograms.orbit_bx_all_byLS.size()))
187  histograms.orbit_bx_all_byLS[iLS]->Fill(bx, orbit);
188 
189  // monitor the bx distribution for the TCDS trigger types
190  size_t size = std::size(s_tcds_trigger_types);
191  if (type < size and histograms.orbit_bx[type]) {
192  histograms.orbit_bx[type]->Fill(bx, orbit);
193  }
194 }
195 
196 //define this as a plug-in
size
Write out results.
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &, RunBasedHistograms &) const override
static const unsigned int s_orbit_range
static const unsigned int s_bx_range
Definition: config.py:1
static std::string to_string(const XMLCh *ch)
dqm::reco::DQMStore DQMStore
~TriggerBxVsOrbitMonitor() override=default
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const edm::EDGetTokenT< GlobalAlgBlkBxCollection > m_l1t_results
TriggerBxVsOrbitMonitor(edm::ParameterSet const &)
void dqmAnalyze(edm::Event const &, edm::EventSetup const &, RunBasedHistograms const &) const override
def ls(path, rec=False)
Definition: eostools.py:349
const edm::EDGetTokenT< edm::TriggerResults > m_hlt_results
dqm::legacy::MonitorElement MonitorElement
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:212
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void dqmBeginRun(edm::Run const &, edm::EventSetup const &, RunBasedHistograms &) const override
HLT enums.
static constexpr const char *const s_tcds_trigger_types[]
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: event.py:1
Definition: Run.h:45