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
28 
29 namespace {
30  struct RunBasedHistograms {
31  ConcurrentMonitorElement orbit_bx_all;
32  std::vector<ConcurrentMonitorElement> orbit_bx;
33  std::vector<ConcurrentMonitorElement> orbit_bx_all_byLS;
34  };
35 }
36 
37 class TriggerBxVsOrbitMonitor : public DQMGlobalEDAnalyzer<RunBasedHistograms> {
38 public:
40  ~TriggerBxVsOrbitMonitor() override = default;
41 
42  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
43 
44 private:
45  void dqmBeginRun(edm::Run const&, edm::EventSetup const&, RunBasedHistograms &) const override;
46  void bookHistograms(DQMStore::ConcurrentBooker &, edm::Run const&, edm::EventSetup const&, RunBasedHistograms &) const override;
47  void dqmAnalyze(edm::Event const&, edm::EventSetup const&, RunBasedHistograms const&) const override;
48 
49  // number of bunch crossings
50  static const unsigned int s_bx_range = 3564;
51  static const unsigned int s_orbit_range = 262144; // 2**18 orbits in 1 LS
52 
53  // TCDS trigger types
54  // see https://twiki.cern.ch/twiki/bin/viewauth/CMS/TcdsEventRecord
55  static constexpr const char * const s_tcds_trigger_types[] = {
56  "Empty", // 0 - No trigger
57  "Physics", // 1 - GT trigger
58  "Calibration", // 2 - Sequence trigger (calibration)
59  "Random", // 3 - Random trigger
60  "Auxiliary", // 4 - Auxiliary (CPM front panel NIM input) trigger
61  nullptr, // 5 - reserved
62  nullptr, // 6 - reserved
63  nullptr, // 7 - reserved
64  "Cyclic", // 8 - Cyclic trigger
65  "Bunch-pattern", // 9 - Bunch-pattern trigger
66  "Software", // 10 - Software trigger
67  "TTS", // 11 - TTS-sourced trigger
68  nullptr, // 12 - reserved
69  nullptr, // 13 - reserved
70  nullptr, // 14 - reserved
71  nullptr // 15 - reserved
72  };
73 
74  // module configuration
78  const int m_minLS;
79  const int m_maxLS;
80  const int m_minBX;
81  const int m_maxBX;
82 };
83 
84 // definition
86 
87 
89 {
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 {
111 }
112 
113 void TriggerBxVsOrbitMonitor::dqmBeginRun(edm::Run const& run, edm::EventSetup const& setup, RunBasedHistograms & histograms) const
114 {
115  size_t nLS = m_maxLS - m_minLS + 1;
116 
117  histograms.orbit_bx_all_byLS.clear();
118  histograms.orbit_bx_all_byLS.resize(nLS);
119 
120  histograms.orbit_bx.clear();
121  histograms.orbit_bx.resize(std::size(s_tcds_trigger_types));
122 }
123 
125 {
126 // TCDS trigger type plots
128  size_t nLS = m_maxLS - m_minLS + 1;
129  size_t nBX = m_maxBX - m_minBX + 1;
130 
131  // book 2D histogram to monitor all TCDS trigger types in a single plot
132  booker.setCurrentFolder( m_dqm_path + "/orbitVsBX" );
133  histograms.orbit_bx_all = booker.book2D(
134  "OrbitVsBX",
135  "Event orbits vs. bunch crossing",
136  nBX, m_minBX - 0.5, m_maxBX + 0.5,
137  s_orbit_range + 1, -0.5, s_orbit_range + 0.5);
138  histograms.orbit_bx_all.setXTitle("BX");
139  histograms.orbit_bx_all.setYTitle("orbit");
140 
141  for (unsigned int i = 0; i < nLS; ++i) {
142  std::string iname = std::to_string(i);
143  histograms.orbit_bx_all_byLS[i] = booker.book2D(
144  "OrbitVsBX_LS" + iname,
145  "Event orbits vs. bunch crossing, for lumisection " + iname,
146  nBX, m_minBX - 0.5, m_maxBX + 0.5,
147  s_orbit_range + 1, -0.5, s_orbit_range + 0.5);
148  histograms.orbit_bx_all_byLS[i].setXTitle("BX");
149  histograms.orbit_bx_all_byLS[i].setYTitle("orbit");
150  }
151 
152  booker.setCurrentFolder( m_dqm_path + "/orbitVsBX/TCDS" );
153  for (unsigned int i = 0; i < size; ++i) {
154  if (s_tcds_trigger_types[i]) {
155  histograms.orbit_bx[i] = booker.book2D(
156  "OrbitVsBX_" + std::string(s_tcds_trigger_types[i]),
157  "Event orbits vs. bunch crossing " + std::string(s_tcds_trigger_types[i]),
158  nBX, m_minBX - 0.5, m_maxBX + 0.5,
159  s_orbit_range + 1, -0.5, s_orbit_range + 0.5);
160  histograms.orbit_bx[i].setXTitle("BX");
161  histograms.orbit_bx[i].setYTitle("orbit");
162  }
163  }
164 }
165 
166 
167 void TriggerBxVsOrbitMonitor::dqmAnalyze(edm::Event const& event, edm::EventSetup const& setup, RunBasedHistograms const& histograms) const
168 {
169  unsigned int type = event.experimentType();
170  unsigned int ls = event.id().luminosityBlock();
171  unsigned int orbit = event.orbitNumber() % s_orbit_range;
172  unsigned int bx = event.bunchCrossing();
173  histograms.orbit_bx_all.fill(bx, orbit);
174 
175  int iLS = ls - m_minLS;
176  if (iLS >= 0 and iLS < int(histograms.orbit_bx_all_byLS.size()))
177  histograms.orbit_bx_all_byLS[iLS].fill(bx, orbit);
178 
179  // monitor the bx distribution for the TCDS trigger types
181  if (type < size and histograms.orbit_bx[type]) {
182  histograms.orbit_bx[type].fill(bx, orbit);
183  }
184 }
185 
186 
187 //define this as a plug-in
size
Write out results.
type
Definition: HCALResponse.h:21
static const char *const s_tcds_trigger_types[]
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
void bookHistograms(DQMStore::ConcurrentBooker &, edm::Run const &, edm::EventSetup const &, RunBasedHistograms &) const override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
static const unsigned int s_orbit_range
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
Definition: config.py:1
#define constexpr
void bookHistograms(fwlite::EventContainer &eventCont)
ConcurrentMonitorElement book2D(Args &&...args)
Definition: DQMStore.h:244
const edm::EDGetTokenT< GlobalAlgBlkBxCollection > m_l1t_results
TriggerBxVsOrbitMonitor(edm::ParameterSet const &)
def ls(path, rec=False)
Definition: eostools.py:348
const edm::EDGetTokenT< edm::TriggerResults > m_hlt_results
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:279
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void dqmAnalyze(edm::Event const &, edm::EventSetup const &, RunBasedHistograms const &) const override
HLT enums.
void dqmBeginRun(edm::Run const &, edm::EventSetup const &, RunBasedHistograms &) const override
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: event.py:1
Definition: Run.h:43