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