Belle II Software  release-05-02-19
DQMHistAnalysisPXDDAQ.cc
1 //+
2 // File : DQMHistAnalysisPXDDAQ.cc
3 // Description : Analysis of PXD Reduction
4 //
5 // Author : Bjoern Spruck, Univerisity Mainz
6 // Date : 2020
7 //-
8 
9 
10 #include <dqm/analysis/modules/DQMHistAnalysisPXDDAQ.h>
11 #include <TROOT.h>
12 
13 using namespace std;
14 using namespace Belle2;
15 
16 //-----------------------------------------------------------------
17 // Register the Module
18 //-----------------------------------------------------------------
19 REG_MODULE(DQMHistAnalysisPXDDAQ)
20 
21 //-----------------------------------------------------------------
22 // Implementation
23 //-----------------------------------------------------------------
24 
27 {
28  // This module CAN NOT be run in parallel!
29 
30  //Parameter definition
31  addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of Histogram dir", std::string("PXDDAQ"));
32  addParam("PVPrefix", m_pvPrefix, "PV Prefix", std::string("DQM:PXD:DAQ:"));
33  addParam("useEpics", m_useEpics, "useEpics", true);
34  B2DEBUG(1, "DQMHistAnalysisPXDDAQ: Constructor done.");
35 
36 }
37 
38 DQMHistAnalysisPXDDAQModule::~DQMHistAnalysisPXDDAQModule()
39 {
40 #ifdef _BELLE2_EPICS
41  if (m_useEpics) {
42  if (ca_current_context()) ca_context_destroy();
43  }
44 #endif
45 }
46 
47 void DQMHistAnalysisPXDDAQModule::initialize()
48 {
49  B2DEBUG(1, "DQMHistAnalysisPXDDAQ: initialized.");
50 
51  m_monObj = getMonitoringObject("pxd");
52 
53  gROOT->cd(); // this seems to be important, or strange things happen
54 
55  m_cDAQError = new TCanvas((m_histogramDirectoryName + "/c_DAQError").data());
56  m_cMissingDHC = new TCanvas((m_histogramDirectoryName + "/c_MissingDHC").data());
57  m_cMissingDHE = new TCanvas((m_histogramDirectoryName + "/c_MissingDHE").data());
58  m_cMissingDHP = new TCanvas((m_histogramDirectoryName + "/c_MissingDHP").data());
59  m_cStatistic = new TCanvas((m_histogramDirectoryName + "/c_Statistic").data());
60 
61  m_monObj->addCanvas(m_cDAQError);
62  m_monObj->addCanvas(m_cMissingDHC);
63  m_monObj->addCanvas(m_cMissingDHE);
64  m_monObj->addCanvas(m_cMissingDHP);
65  m_monObj->addCanvas(m_cStatistic);
66 
67  m_hMissingDHC = new TH2F("hPXDMissingDHC", "PXD Missing DHC", 16, 0, 16, 2, 0, 2);
68  m_hMissingDHE = new TH2F("hPXDMissingDHE", "PXD Missing DHE", 64, 0, 64, 2, 0, 2);
69 
70 #ifdef _BELLE2_EPICS
71  mychid.resize(20);
72  if (m_useEpics) {
73  if (!ca_current_context()) SEVCHK(ca_context_create(ca_disable_preemptive_callback), "ca_context_create");
74  SEVCHK(ca_create_channel((m_pvPrefix + "HLTRej").data(), NULL, NULL, 10, &mychid[0]), "ca_create_channel failure");
75  SEVCHK(ca_create_channel((m_pvPrefix + "Trunc").data(), NULL, NULL, 10, &mychid[1]), "ca_create_channel failure");
76  SEVCHK(ca_create_channel((m_pvPrefix + "HER_Trunc").data(), NULL, NULL, 10, &mychid[2]), "ca_create_channel failure");
77  SEVCHK(ca_create_channel((m_pvPrefix + "LER_Trunc").data(), NULL, NULL, 10, &mychid[3]), "ca_create_channel failure");
78  SEVCHK(ca_create_channel((m_pvPrefix + "CM63").data(), NULL, NULL, 10, &mychid[4]), "ca_create_channel failure");
79  SEVCHK(ca_create_channel((m_pvPrefix + "HER_CM63").data(), NULL, NULL, 10, &mychid[5]), "ca_create_channel failure");
80  SEVCHK(ca_create_channel((m_pvPrefix + "LER_CM63").data(), NULL, NULL, 10, &mychid[6]), "ca_create_channel failure");
81  SEVCHK(ca_create_channel((m_pvPrefix + "HER_CM63_1ms").data(), NULL, NULL, 10, &mychid[7]), "ca_create_channel failure");
82  SEVCHK(ca_create_channel((m_pvPrefix + "LER_CM63_1ms").data(), NULL, NULL, 10, &mychid[8]), "ca_create_channel failure");
83  SEVCHK(ca_create_channel((m_pvPrefix + "HER_Trunc_1ms").data(), NULL, NULL, 10, &mychid[9]), "ca_create_channel failure");
84  SEVCHK(ca_create_channel((m_pvPrefix + "LER_Trunc_1ms").data(), NULL, NULL, 10, &mychid[10]), "ca_create_channel failure");
85  SEVCHK(ca_create_channel((m_pvPrefix + "MissFrame").data(), NULL, NULL, 10, &mychid[11]), "ca_create_channel failure");
86  SEVCHK(ca_create_channel((m_pvPrefix + "Timeout").data(), NULL, NULL, 10, &mychid[12]), "ca_create_channel failure");
87  SEVCHK(ca_create_channel((m_pvPrefix + "LinkDown").data(), NULL, NULL, 10, &mychid[13]), "ca_create_channel failure");
88  SEVCHK(ca_create_channel((m_pvPrefix + "Mismatch").data(), NULL, NULL, 10, &mychid[14]), "ca_create_channel failure");
89  SEVCHK(ca_create_channel((m_pvPrefix + "HER_Miss").data(), NULL, NULL, 10, &mychid[15]), "ca_create_channel failure");
90  SEVCHK(ca_create_channel((m_pvPrefix + "LER_Miss").data(), NULL, NULL, 10, &mychid[16]), "ca_create_channel failure");
91  SEVCHK(ca_create_channel((m_pvPrefix + "HER_Miss_1ms").data(), NULL, NULL, 10, &mychid[17]), "ca_create_channel failure");
92  SEVCHK(ca_create_channel((m_pvPrefix + "LER_Miss_1ms").data(), NULL, NULL, 10, &mychid[18]), "ca_create_channel failure");
93  SEVCHK(ca_create_channel((m_pvPrefix + "unused").data(), NULL, NULL, 10, &mychid[19]), "ca_create_channel failure");
94  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
95  }
96 #endif
97 }
98 
99 void DQMHistAnalysisPXDDAQModule::beginRun()
100 {
101  B2DEBUG(1, "DQMHistAnalysisPXDDAQ: beginRun called.");
102 
103  m_cMissingDHP->Clear();
104 }
105 
106 void DQMHistAnalysisPXDDAQModule::event()
107 {
108 // double data = 0.0;
109  if (m_cMissingDHP == nullptr || m_cMissingDHE == nullptr || m_cMissingDHC == nullptr
110  || m_cStatistic == nullptr) return; // we could assume this
111 
112  m_cDAQError->Clear();
113  m_cMissingDHP->Clear();
114  m_cMissingDHE->Clear();
115  m_cMissingDHC->Clear();
116  m_cStatistic->Clear();
117 
118 
119  {
120  std::string name = "PXDDAQError";
121 
122  if (m_hDAQError) { delete m_hDAQError; m_hDAQError = nullptr;}
123 
124  TH1* hh1 = findHist(name);
125  if (hh1 == NULL) {
126  hh1 = findHist(m_histogramDirectoryName, name);
127  }
128  m_cDAQError->cd();
129  if (hh1) {
130  m_hDAQError = (TH1F*)hh1->DrawClone("text");
131  m_hDAQError->SetName("hPXDDAQError");
132  m_hDAQError->SetTitle("PXD Fraction of DAQ Errors");
133  if (m_hDAQError->GetBinContent(0)) {
134  m_hDAQError->Scale(1.0 / m_hDAQError->GetBinContent(0));
135  }
136  m_hDAQError->Draw("text,hist");
137  }
138  }
139  {
140  // DHC histogram
141  std::string name = "PXDDAQDHCError";
142 
143  TH1* hh1 = findHist(name);
144  if (hh1 == NULL) {
145  hh1 = findHist(m_histogramDirectoryName, name);
146  }
147  m_cMissingDHC->cd();
148  if (hh1) {
149  auto events = hh1->GetBinContent(hh1->GetBin(-1, -1));
150  // first, we have to relate the per-DHC overflow (DHC object count) to the overall overflow (event count)
151  // second, we have to relate the "fake data" DHC bin to the per-DHC overflow (DHC object count)
152  m_hMissingDHC->Reset();
153  for (int i = 0; i < 16; i++) {
154  auto dhecount = hh1->GetBinContent(hh1->GetBin(i, -1));
155  if (events > 0) m_hMissingDHC->Fill((double)i, 0.0, 1.0 - dhecount / events);
156  // c_FAKE_NO_DATA_TRIG = 1ull << 29,
157  if (dhecount > 0) m_hMissingDHC->Fill((double)i, 1.0, hh1->GetBinContent(hh1->GetBin(i, 29) / dhecount));
158  }
159  m_hMissingDHC->Draw("text");
160  }
161  }
162 
163  {
164  // DHE histogram
165  std::string name = "PXDDAQDHEError";
166 
167  TH1* hh1 = findHist(name);
168  if (hh1 == NULL) {
169  hh1 = findHist(m_histogramDirectoryName, name);
170  }
171  m_cMissingDHE->cd();
172  if (hh1) {
173  auto events = hh1->GetBinContent(hh1->GetBin(-1, -1));
174  // first, we have to relate the per-DHE overflow (DHE object count) to the overall overflow (event count)
175  // second, we have to relate the "fake data" DHE bin to the per-DHE overflow (DHE object count)
176  m_hMissingDHE->Reset();
177  for (int i = 0; i < 64; i++) {
178  auto dhecount = hh1->GetBinContent(hh1->GetBin(i, -1));
179  if (events > 0) m_hMissingDHE->Fill((double)i, 0.0, 1.0 - dhecount / events);
180  // c_FAKE_NO_DATA_TRIG = 1ull << 29,
181  if (dhecount > 0) m_hMissingDHE->Fill((double)i, 1.0, hh1->GetBinContent(hh1->GetBin(i, 29) / dhecount));
182  }
183  m_hMissingDHE->Draw("text");
184  }
185  }
186 
187  {
188  // DHP histogram
189  if (m_hMissingDHP) { delete m_hMissingDHP; m_hMissingDHP = nullptr;}
190 
191  std::string name = "PXDDAQDHPDataMissing";
192 
193  TH1* hh1 = findHist(name);
194  if (hh1 == NULL) {
195  hh1 = findHist(m_histogramDirectoryName, name);
196  }
197  m_cMissingDHP->cd();
198  if (hh1) {
199  m_hMissingDHP = (TH1F*)hh1->DrawClone("text");
200  if (m_hMissingDHP->GetBinContent(0)) {
201  m_hMissingDHP->Scale(1.0 / m_hMissingDHP->GetBinContent(0));
202  m_hMissingDHP->Draw("text");
203  }
204  }
205  // double data = m_hMissingDHP->Max???;
206  //
207  // m_monObj->setVariable("missingDHPFraction", data);
208  //
209  m_cMissingDHP->Modified();
210  m_cMissingDHP->Update();
211  }
212 
213  double data_HLTRej = 0.0;
214  double data_Trunc = 0.0;
215  double data_HER_Trunc = 0.0;
216  double data_LER_Trunc = 0.0;
217  double data_CM63 = 0.0;
218  double data_HER_CM63 = 0.0;
219  double data_LER_CM63 = 0.0;
220  double data_HER_CM63_1ms = 0.0;
221  double data_LER_CM63_1ms = 0.0;
222  double data_HER_Trunc_1ms = 0.0;
223  double data_LER_Trunc_1ms = 0.0;
224  double data_MissFrame = 0.0;
225  double data_Timeout = 0.0;
226  double data_LinkDown = 0.0;
227  double data_Mismatch = 0.0;
228  double data_HER_Miss = 0.0;
229  double data_LER_Miss = 0.0;
230  double data_HER_Miss_1ms = 0.0;
231  double data_LER_Miss_1ms = 0.0;
232  double data_unused = 0.0;
233 
234  // Stat histogram
235  if (m_hStatistic) { delete m_hStatistic; m_hStatistic = nullptr;}
236 
237  std::string name = "PXDDAQStat";
238 
239  TH1* hh1 = findHist(name);
240  if (hh1 == NULL) {
241  hh1 = findHist(m_histogramDirectoryName, name);
242  }
243  m_cStatistic->cd();
244  if (hh1) {
245  m_hStatistic = (TH1F*)hh1->DrawClone("text");
246  if (m_hStatistic->GetBinContent(0)) {
247  m_hStatistic->Scale(1.0 / m_hStatistic->GetBinContent(0));
248  m_hStatistic->Draw("text");
249  }
250  double scale = hh1->GetBinContent(0);// underflow is event counter
251  if (scale != 0.0) scale = 1.0 / scale; // just avoid dive by zero
252  data_HLTRej = hh1->GetBinContent(1 + 0) * scale;
253  data_Trunc = hh1->GetBinContent(1 + 1) * scale;
254  data_HER_Trunc = hh1->GetBinContent(1 + 2) * scale;
255  data_LER_Trunc = hh1->GetBinContent(1 + 3) * scale;
256  data_CM63 = hh1->GetBinContent(1 + 4) * scale;
257  data_HER_CM63 = hh1->GetBinContent(1 + 5) * scale;
258  data_LER_CM63 = hh1->GetBinContent(1 + 6) * scale;
259  data_HER_CM63_1ms = hh1->GetBinContent(1 + 7) * scale;
260  data_LER_CM63_1ms = hh1->GetBinContent(1 + 8) * scale;
261  data_HER_Trunc_1ms = hh1->GetBinContent(1 + 9) * scale;
262  data_LER_Trunc_1ms = hh1->GetBinContent(1 + 10) * scale;
263  data_MissFrame = hh1->GetBinContent(1 + 11) * scale;
264  data_Timeout = hh1->GetBinContent(1 + 12) * scale;
265  data_LinkDown = hh1->GetBinContent(1 + 13) * scale;
266  data_Mismatch = hh1->GetBinContent(1 + 14) * scale;
267  data_HER_Miss = hh1->GetBinContent(1 + 15) * scale;
268  data_LER_Miss = hh1->GetBinContent(1 + 16) * scale;
269  data_HER_Miss_1ms = hh1->GetBinContent(1 + 17) * scale;
270  data_LER_Miss_1ms = hh1->GetBinContent(1 + 18) * scale;
271  data_unused = hh1->GetBinContent(1 + 19) * scale;
272  }
273 
274  m_monObj->setVariable("HLTReject", data_HLTRej);
275  m_monObj->setVariable("Trunc", data_Trunc);
276  m_monObj->setVariable("HER_Trunc", data_HER_Trunc);
277  m_monObj->setVariable("LER_Trunc", data_LER_Trunc);
278  m_monObj->setVariable("CM63", data_CM63);
279  m_monObj->setVariable("HER_CM63", data_HER_CM63);
280  m_monObj->setVariable("LER_CM63", data_LER_CM63);
281  m_monObj->setVariable("HER_CM63_1ms", data_HER_CM63_1ms);
282  m_monObj->setVariable("LER_CM63_1ms", data_LER_CM63_1ms);
283  m_monObj->setVariable("HER_Trunc_1ms", data_HER_Trunc_1ms);
284  m_monObj->setVariable("LER_Trunc_1ms", data_LER_Trunc_1ms);
285  m_monObj->setVariable("MissFrame", data_MissFrame);
286  m_monObj->setVariable("Timeout", data_Timeout);
287  m_monObj->setVariable("LinkDown", data_LinkDown);
288  m_monObj->setVariable("Mismatch", data_Mismatch);
289  m_monObj->setVariable("HER_Miss", data_HER_Miss);
290  m_monObj->setVariable("LER_Miss", data_LER_Miss);
291  m_monObj->setVariable("HER_Miss_1ms", data_HER_Miss_1ms);
292  m_monObj->setVariable("LER_Miss_1ms", data_LER_Miss_1ms);
293 
294  m_cStatistic->Modified();
295  m_cStatistic->Update();
296 
297 #ifdef _BELLE2_EPICS
298  if (m_useEpics) {
299  SEVCHK(ca_put(DBR_DOUBLE, mychid[0], (void*)&data_HLTRej), "ca_set failure");
300  SEVCHK(ca_put(DBR_DOUBLE, mychid[1], (void*)&data_Trunc), "ca_set failure");
301  SEVCHK(ca_put(DBR_DOUBLE, mychid[2], (void*)&data_HER_Trunc), "ca_set failure");
302  SEVCHK(ca_put(DBR_DOUBLE, mychid[3], (void*)&data_LER_Trunc), "ca_set failure");
303  SEVCHK(ca_put(DBR_DOUBLE, mychid[4], (void*)&data_CM63), "ca_set failure");
304  SEVCHK(ca_put(DBR_DOUBLE, mychid[5], (void*)&data_HER_CM63), "ca_set failure");
305  SEVCHK(ca_put(DBR_DOUBLE, mychid[6], (void*)&data_LER_CM63), "ca_set failure");
306  SEVCHK(ca_put(DBR_DOUBLE, mychid[7], (void*)&data_HER_CM63_1ms), "ca_set failure");
307  SEVCHK(ca_put(DBR_DOUBLE, mychid[8], (void*)&data_LER_CM63_1ms), "ca_set failure");
308  SEVCHK(ca_put(DBR_DOUBLE, mychid[9], (void*)&data_HER_Trunc_1ms), "ca_set failure");
309  SEVCHK(ca_put(DBR_DOUBLE, mychid[10], (void*)&data_LER_Trunc_1ms), "ca_set failure");
310  SEVCHK(ca_put(DBR_DOUBLE, mychid[11], (void*)&data_MissFrame), "ca_set failure");
311  SEVCHK(ca_put(DBR_DOUBLE, mychid[12], (void*)&data_Timeout), "ca_set failure");
312  SEVCHK(ca_put(DBR_DOUBLE, mychid[13], (void*)&data_LinkDown), "ca_set failure");
313  SEVCHK(ca_put(DBR_DOUBLE, mychid[14], (void*)&data_Mismatch), "ca_set failure");
314  SEVCHK(ca_put(DBR_DOUBLE, mychid[15], (void*)&data_HER_Miss), "ca_set failure");
315  SEVCHK(ca_put(DBR_DOUBLE, mychid[16], (void*)&data_LER_Miss), "ca_set failure");
316  SEVCHK(ca_put(DBR_DOUBLE, mychid[17], (void*)&data_HER_Miss_1ms), "ca_set failure");
317  SEVCHK(ca_put(DBR_DOUBLE, mychid[18], (void*)&data_LER_Miss_1ms), "ca_set failure");
318  SEVCHK(ca_put(DBR_DOUBLE, mychid[19], (void*)&data_unused), "ca_set failure");
319  // write out
320  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
321  }
322 #endif
323 }
324 
325 void DQMHistAnalysisPXDDAQModule::terminate()
326 {
327  B2DEBUG(1, "DQMHistAnalysisPXDDAQ: terminate called");
328  // should delete canvas here, maybe hist, too? Who owns it?
329 }
330 
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::DQMHistAnalysisPXDDAQModule
DQM Histogram Analysis for PXD DAQ.
Definition: DQMHistAnalysisPXDDAQ.h:29
Belle2::DQMHistAnalysisModule
The base class for the histogram analysis module.
Definition: DQMHistAnalysis.h:27