Belle II Software  release-06-02-00
DQMHistInjection.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 : DQMHistInjection.cc
10 // Description : DQM module, which gives histograms showing occupancies after injection
11 //-
12 
13 
14 #include <dqm/analysis/modules/DQMHistInjection.h>
15 #include <klm/dataobjects/KLMElementNumbers.h>
16 #include <TROOT.h>
17 
18 using namespace std;
19 using namespace Belle2;
20 
21 //-----------------------------------------------------------------
22 // Register the Module
23 //-----------------------------------------------------------------
24 REG_MODULE(DQMHistInjection)
25 
26 //-----------------------------------------------------------------
27 // Implementation
28 //-----------------------------------------------------------------
29 
31 {
32  // This module CAN NOT be run in parallel!
33 
34 // addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms were placed", std::string("PXDINJ"));
35  addParam("PVPrefix", m_pvPrefix, "PV Prefix", std::string("DQM:INJ:"));
36  addParam("useEpics", m_useEpics, "Whether to update EPICS PVs.", false);
37  B2DEBUG(1, "DQMHistInjection: Constructor done.");
38 }
39 
40 DQMHistInjectionModule::~DQMHistInjectionModule()
41 {
42 #ifdef _BELLE2_EPICS
43  if (m_useEpics) {
44  if (ca_current_context()) ca_context_destroy();
45  }
46 #endif
47 }
48 
49 void DQMHistInjectionModule::initialize()
50 {
51 
52  gROOT->cd(); // this seems to be important, or strange things happen
53 
54  m_cInjectionLERPXD = new TCanvas("PXDINJ/c_InjectionLERPXD");
55  m_cInjectionLERPXDOcc = new TCanvas("PXDINJ/c_InjectionLERPXDOcc");
56  m_cInjectionLERSVD = new TCanvas("SVDInjection/c_InjectionLERSVD");
57  m_cInjectionLERSVDOcc = new TCanvas("SVDInjection/c_InjectionLERSVDOcc");
58  m_cInjectionLERECL = new TCanvas("ECLINJ/c_InjectionLERECL");
59  m_cBurstLERECL = new TCanvas("ECLINJ/c_BurstInjectionLERECL");
60  m_cInjectionLERTOP = new TCanvas("TOP/c_InjectionLERTOP");
61  m_cInjectionLERARICH = new TCanvas("ARICH/c_InjectionLERARICH");
62  m_cInjectionLERKLM = new TCanvas("KLM/c_InjectionLERKLM");
63 
64  m_cInjectionHERPXD = new TCanvas("PXDINJ/c_InjectionHERPXD");
65  m_cInjectionHERPXDOcc = new TCanvas("PXDINJ/c_InjectionHERPXDOcc");
66  m_cInjectionHERSVD = new TCanvas("SVDInjection/c_InjectionHERSVD");
67  m_cInjectionHERSVDOcc = new TCanvas("SVDInjection/c_InjectionHERSVDOcc");
68  m_cInjectionHERECL = new TCanvas("ECLINJ/c_InjectionHERECL");
69  m_cBurstHERECL = new TCanvas("ECLINJ/c_BurstInjectionHERECL");
70  m_cInjectionHERTOP = new TCanvas("TOP/c_InjectionHERTOP");
71  m_cInjectionHERARICH = new TCanvas("ARICH/c_InjectionHERARICH");
72  m_cInjectionHERKLM = new TCanvas("KLM/c_InjectionHERKLM");
73 
74  m_hInjectionLERPXD = new TH1F("HitInjectionLERPXD", "PXD Hits after LER Injection;Time in #mus;Mean Hits/event", 4000, 0 , 20000);
75  m_hInjectionLERPXDOcc = new TH1F("HitInjectionPXDLEROcc", "PXD Occ after LER Injection;Time in #mus;Mean Occ in % per module", 4000,
76  0 , 20000);
77  m_hInjectionLERSVD = new TH1F("HitInjectionLERSVD", "SVD Hits after LER Injection;Time in #mus;Mean Hits/event", 4000, 0 , 20000);
78  m_hInjectionLERSVDOcc = new TH1F("HitInjectionSVDLEROcc", "SVD Occ after LER Injection;Time in #mus;Mean Occ in % per module", 4000,
79  0 , 20000);
80  m_hInjectionLERECL = new TH1F("HitInjectionLERECL", "ECL Hits after LER Injection;Time in #mus;Mean Hits/event", 4000, 0 , 20000);
81  m_hBurstLERECL = new TH1F("BurstInjectionLERECL", "ECL Bursts after LER Injection;Time in #mus;Suppressions/event (1 #mus bins)",
82  20000, 0 , 20000);
83  m_hInjectionLERTOP = new TH1F("HitInjectionLERTOP", "TOP Occ after LER Injection;Time in #mus;Mean Occ in % /event", 4000, 0 ,
84  20000);
85  m_hInjectionLERARICH = new TH1F("HitInjectionLERARICH", "ARICH Occ after LER Injection;Time in #mus;Mean Hits/event", 4000, 0 ,
86  20000);
87  m_hInjectionLERKLM = new TH1F("HitInjectionLERKLM",
88  "KLM occupancy after LER Injection;Time [#mus];Digits occupancy in % / (5 #mus)", 4000, 0 ,
89  20000);
90 
91  m_hInjectionHERPXD = new TH1F("HitInjectionHERPXD", "PXD Hits after HER Injection;Time in #mus;Mean Hits/event", 4000, 0 , 20000);
92  m_hInjectionHERPXDOcc = new TH1F("HitInjectionPXDHEROcc", "PXD Occ after HER Injection;Time in #mus;Mean Occ in % per modul", 4000,
93  0 , 20000);
94  m_hInjectionHERSVD = new TH1F("HitInjectionHERSVD", "SVD Hits after HER Injection;Time in #mus;Mean Hits/event", 4000, 0 , 20000);
95  m_hInjectionHERSVDOcc = new TH1F("HitInjectionSVDHEROcc", "SVD Occ after HER Injection;Time in #mus;Mean Occ in % per modul", 4000,
96  0 , 20000);
97  m_hInjectionHERECL = new TH1F("HitInjectionHERECL", "ECL Hits after HER Injection;Time in #mus;Mean Hits/event", 4000, 0 , 20000);
98  m_hBurstHERECL = new TH1F("BurstInjectionHERECL", "ECL Bursts after HER Injection;Time in #mus;Suppressions/event (1 #mus bins)",
99  20000, 0 , 20000);
100  m_hInjectionHERTOP = new TH1F("HitInjectionHERTOP", "TOP Occ after HER Injection;Time in #mus;Mean Occ in % /event", 4000, 0 ,
101  20000);
102  m_hInjectionHERARICH = new TH1F("HitInjectionHERARICH", "ARICH Occ after HER Injection;Time in #mus;Mean Hits/event", 4000, 0 ,
103  20000);
104  m_hInjectionHERKLM = new TH1F("HitInjectionHERKLM",
105  "KLM occupancy after HER Injection;Time [#mus];Digits occupancy in % / (5 #mus)", 4000, 0 ,
106  20000);
107 
108 #ifdef _BELLE2_EPICS
109  if (m_useEpics) {
110  if (!ca_current_context()) SEVCHK(ca_context_create(ca_disable_preemptive_callback), "ca_context_create");
111  m_nodes.resize(14);
112  SEVCHK(ca_create_channel((m_pvPrefix + "LER:Triggers").data(), NULL, NULL, 10, &m_nodes[0].mychid), "ca_create_channel failure");
113  m_nodes[0].histo = nullptr;
114  SEVCHK(ca_create_channel((m_pvPrefix + "LER:PXD").data(), NULL, NULL, 10, &m_nodes[1].mychid), "ca_create_channel failure");
115  m_nodes[1].histo = m_hInjectionLERPXD;
116  SEVCHK(ca_create_channel((m_pvPrefix + "LER:ECL").data(), NULL, NULL, 10, &m_nodes[2].mychid), "ca_create_channel failure");
117  m_nodes[2].histo = m_hInjectionLERECL;
118  SEVCHK(ca_create_channel((m_pvPrefix + "HER:Triggers").data(), NULL, NULL, 10, &m_nodes[3].mychid), "ca_create_channel failure");
119  m_nodes[3].histo = nullptr;
120  SEVCHK(ca_create_channel((m_pvPrefix + "HER:PXD").data(), NULL, NULL, 10, &m_nodes[4].mychid), "ca_create_channel failure");
121  m_nodes[4].histo = m_hInjectionHERPXD;
122  SEVCHK(ca_create_channel((m_pvPrefix + "HER:ECL").data(), NULL, NULL, 10, &m_nodes[5].mychid), "ca_create_channel failure");
123  m_nodes[5].histo = m_hInjectionHERECL;
124  SEVCHK(ca_create_channel((m_pvPrefix + "LER:TOP").data(), NULL, NULL, 10, &m_nodes[6].mychid), "ca_create_channel failure");
125  m_nodes[6].histo = m_hInjectionLERTOP;
126  SEVCHK(ca_create_channel((m_pvPrefix + "HER:TOP").data(), NULL, NULL, 10, &m_nodes[7].mychid), "ca_create_channel failure");
127  m_nodes[7].histo = m_hInjectionHERTOP;
128  SEVCHK(ca_create_channel((m_pvPrefix + "LER:SVD").data(), NULL, NULL, 10, &m_nodes[8].mychid), "ca_create_channel failure");
129  m_nodes[8].histo = m_hInjectionLERSVD;
130  SEVCHK(ca_create_channel((m_pvPrefix + "HER:SVD").data(), NULL, NULL, 10, &m_nodes[9].mychid), "ca_create_channel failure");
131  m_nodes[9].histo = m_hInjectionHERSVD;
132  SEVCHK(ca_create_channel((m_pvPrefix + "LER:ARICH").data(), NULL, NULL, 10, &m_nodes[10].mychid), "ca_create_channel failure");
133  m_nodes[10].histo = m_hInjectionLERARICH;
134  SEVCHK(ca_create_channel((m_pvPrefix + "HER:ARICH").data(), NULL, NULL, 10, &m_nodes[11].mychid), "ca_create_channel failure");
135  m_nodes[11].histo = m_hInjectionHERARICH;
136  SEVCHK(ca_create_channel((m_pvPrefix + "LER:KLM").data(), NULL, NULL, 10, &m_nodes[10].mychid), "ca_create_channel failure");
137  m_nodes[12].histo = m_hInjectionLERKLM;
138  SEVCHK(ca_create_channel((m_pvPrefix + "HER:KLM").data(), NULL, NULL, 10, &m_nodes[11].mychid), "ca_create_channel failure");
139  m_nodes[13].histo = m_hInjectionHERKLM;
140 
141  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
142  cleanPVs();
143  }
144 #endif
145  B2DEBUG(1, "DQMHistInjection: initialized.");
146 }
147 
148 
149 void DQMHistInjectionModule::beginRun()
150 {
151  B2DEBUG(1, "DQMHistInjection: beginRun called.");
152 
153 // m_cInjectionLERPXD->Clear(); // FIXME, unclear if this lets to crashes on new run?
154 // m_cInjectionLERPXDOcc->Clear();
155 // m_cInjectionLERECL->Clear();
156 // m_cInjectionHERPXD->Clear();
157 // m_cInjectionHERPXDOcc->Clear();
158 // m_cInjectionHERECL->Clear();
159 
160 // cleanPVs();
161 }
162 
163 
164 void DQMHistInjectionModule::event()
165 {
166  TH1* Hits = nullptr, *Triggers = nullptr;
167  TString locationHits = "";
168  TString locationTriggers = "";
169  //PXD
170  m_histogramDirectoryName = "PXDINJ";
171 
172  locationHits = "PXDOccInjLER";
173  if (m_histogramDirectoryName != "") {
174  locationHits = m_histogramDirectoryName + "/" + locationHits;
175  }
176  Hits = (TH1*)findHist(locationHits.Data());
177  locationTriggers = "PXDEOccInjLER";
178  if (m_histogramDirectoryName != "") {
179  locationTriggers = m_histogramDirectoryName + "/" + locationTriggers;
180  }
181  Triggers = (TH1*)findHist(locationTriggers.Data());
182 
183  //Finding only one of them should only happen in very strange situations...
184  //m_nodes[0].histo = Triggers;
185  if (Hits && Triggers) {
186  m_hInjectionLERPXD->Divide(Hits, Triggers);
187  m_hInjectionLERPXDOcc->Divide(Hits, Triggers, 100, 768 * 250); // to percent
188  }
189 
190  m_cInjectionLERPXD->Clear();
191  m_cInjectionLERPXD->cd(0);
192  m_hInjectionLERPXD->Draw("hist");
193 
194  m_cInjectionLERPXDOcc->Clear();
195  m_cInjectionLERPXDOcc->cd(0);
196  m_hInjectionLERPXDOcc->Draw("hist");
197 
198  locationHits = "PXDOccInjHER";
199  if (m_histogramDirectoryName != "") {
200  locationHits = m_histogramDirectoryName + "/" + locationHits;
201  }
202  Hits = (TH1*)findHist(locationHits.Data());
203  locationTriggers = "PXDEOccInjHER";
204  if (m_histogramDirectoryName != "") {
205  locationTriggers = m_histogramDirectoryName + "/" + locationTriggers;
206  }
207  Triggers = (TH1*)findHist(locationTriggers.Data());
208 
209  //Finding only one of them should only happen in very strange situations...
210  //m_nodes[3].histo = Triggers;
211  if (Hits && Triggers) {
212  m_hInjectionHERPXD->Divide(Hits, Triggers);
213  m_hInjectionHERPXDOcc->Divide(Hits, Triggers, 100, 768 * 250); // to percent
214  }
215 
216  m_cInjectionHERPXD->Clear();
217  m_cInjectionHERPXD->cd(0);
218  m_hInjectionHERPXD->Draw("hist");
219 
220  m_cInjectionHERPXDOcc->Clear();
221  m_cInjectionHERPXDOcc->cd(0);
222  m_hInjectionHERPXDOcc->Draw("hist");
223 
224  //SVD
225  m_histogramDirectoryName = "SVDInjection";
226 
227  locationHits = "SVDOccInjLER";
228  if (m_histogramDirectoryName != "") {
229  locationHits = m_histogramDirectoryName + "/" + locationHits;
230  }
231  Hits = (TH1*)findHist(locationHits.Data());
232  locationTriggers = "SVDTrgOccInjLER";
233  if (m_histogramDirectoryName != "") {
234  locationTriggers = m_histogramDirectoryName + "/" + locationTriggers;
235  }
236  Triggers = (TH1*)findHist(locationTriggers.Data());
237 
238  //Finding only one of them should only happen in very strange situations...
239  //m_nodes[0].histo = Triggers;
240  if (Hits && Triggers) {
241  m_hInjectionLERSVD->Divide(Hits, Triggers);
242  m_hInjectionLERSVDOcc->Divide(Hits, Triggers, 100, 768 * 7 * 2); // to percent (L3V has 768 strips * 2 * 7 sides)
243  }
244 
245  m_cInjectionLERSVD->Clear();
246  m_cInjectionLERSVD->cd(0);
247  m_hInjectionLERSVD->Draw("hist");
248 
249  m_cInjectionLERSVDOcc->Clear();
250  m_cInjectionLERSVDOcc->cd(0);
251  m_hInjectionLERSVDOcc->Draw("hist");
252 
253  locationHits = "SVDOccInjHER";
254  if (m_histogramDirectoryName != "") {
255  locationHits = m_histogramDirectoryName + "/" + locationHits;
256  }
257  Hits = (TH1*)findHist(locationHits.Data());
258  locationTriggers = "SVDTrgOccInjHER";
259  if (m_histogramDirectoryName != "") {
260  locationTriggers = m_histogramDirectoryName + "/" + locationTriggers;
261  }
262  Triggers = (TH1*)findHist(locationTriggers.Data());
263 
264  //Finding only one of them should only happen in very strange situations...
265  //m_nodes[3].histo = Triggers;
266  if (Hits && Triggers) {
267  m_hInjectionHERSVD->Divide(Hits, Triggers);
268  m_hInjectionHERSVDOcc->Divide(Hits, Triggers, 100, 768 * 2 * 7); // to percent (L3V has 768 strips * 2 * 7 sides)
269  }
270 
271  m_cInjectionHERSVD->Clear();
272  m_cInjectionHERSVD->cd(0);
273  m_hInjectionHERSVD->Draw("hist");
274 
275  m_cInjectionHERSVDOcc->Clear();
276  m_cInjectionHERSVDOcc->cd(0);
277  m_hInjectionHERSVDOcc->Draw("hist");
278 
279 
280  //ECL
281  m_histogramDirectoryName = "ECLINJ";
282 
283  locationHits = "ECLHitsInjLER";
284  if (m_histogramDirectoryName != "") {
285  locationHits = m_histogramDirectoryName + "/" + locationHits;
286  }
287  Hits = (TH1*)findHist(locationHits.Data());
288  locationTriggers = "ECLEHitsInjLER";
289  if (m_histogramDirectoryName != "") {
290  locationTriggers = m_histogramDirectoryName + "/" + locationTriggers;
291  }
292  Triggers = (TH1*)findHist(locationTriggers.Data());
293 
294  //Finding only one of them should only happen in very strange situations...
295 #ifdef _BELLE2_EPICS
296  if (m_useEpics) {
297  m_nodes[0].histo = Triggers;
298  }
299 #endif
300  if (Hits && Triggers) {
301  m_hInjectionLERECL->Divide(Hits, Triggers);
302  }
303 
304  m_cInjectionLERECL->Clear();
305  m_cInjectionLERECL->cd(0);
306  m_hInjectionLERECL->Draw("hist");
307 
308  locationHits = "ECLHitsInjHER";
309  if (m_histogramDirectoryName != "") {
310  locationHits = m_histogramDirectoryName + "/" + locationHits;
311  }
312  Hits = (TH1*)findHist(locationHits.Data());
313  locationTriggers = "ECLEHitsInjHER";
314  if (m_histogramDirectoryName != "") {
315  locationTriggers = m_histogramDirectoryName + "/" + locationTriggers;
316  }
317  Triggers = (TH1*)findHist(locationTriggers.Data());
318 
319  //Finding only one of them should only happen in very strange situations...
320 #ifdef _BELLE2_EPICS
321  if (m_useEpics) {
322  m_nodes[3].histo = Triggers;
323  }
324 #endif
325  if (Hits && Triggers) {
326  m_hInjectionHERECL->Divide(Hits, Triggers);
327  }
328 
329  m_cInjectionHERECL->Clear();
330  m_cInjectionHERECL->cd(0);
331  m_hInjectionHERECL->Draw("hist");
332 // =====
333  locationHits = "ECLBurstsInjLER";
334  if (m_histogramDirectoryName != "") {
335  locationHits = m_histogramDirectoryName + "/" + locationHits;
336  }
337  Hits = (TH1*)findHist(locationHits.Data());
338  locationTriggers = "ECLEBurstsInjLER";
339  if (m_histogramDirectoryName != "") {
340  locationTriggers = m_histogramDirectoryName + "/" + locationTriggers;
341  }
342  Triggers = (TH1*)findHist(locationTriggers.Data());
343 
344  if (Hits && Triggers) {
345  m_hBurstLERECL->Divide(Hits, Triggers);
346  }
347 
348  m_cBurstLERECL->Clear();
349  m_cBurstLERECL->cd(0);
350  m_hBurstLERECL->Draw("hist");
351 // =====
352 
353  locationHits = "ECLBurstsInjHER";
354  if (m_histogramDirectoryName != "") {
355  locationHits = m_histogramDirectoryName + "/" + locationHits;
356  }
357  Hits = (TH1*)findHist(locationHits.Data());
358  locationTriggers = "ECLEBurstsInjHER";
359  if (m_histogramDirectoryName != "") {
360  locationTriggers = m_histogramDirectoryName + "/" + locationTriggers;
361  }
362  Triggers = (TH1*)findHist(locationTriggers.Data());
363 
364  if (Hits && Triggers) {
365  m_hBurstHERECL->Divide(Hits, Triggers);
366  }
367 
368  m_cBurstHERECL->Clear();
369  m_cBurstHERECL->cd(0);
370  m_hBurstHERECL->Draw("hist");
371 // =====
372 
373 
374  //TOP
375  m_histogramDirectoryName = "TOP";
376 
377  locationHits = "TOPOccInjLER";
378  if (m_histogramDirectoryName != "") {
379  locationHits = m_histogramDirectoryName + "/" + locationHits;
380  }
381  Hits = (TH1*)findHist(locationHits.Data());
382  locationTriggers = "TOPEOccInjLER";
383  if (m_histogramDirectoryName != "") {
384  locationTriggers = m_histogramDirectoryName + "/" + locationTriggers;
385  }
386  Triggers = (TH1*)findHist(locationTriggers.Data());
387 
388  if (Hits && Triggers) {
389  m_hInjectionLERTOP->Divide(Hits, Triggers, 100, 8192);
390  }
391 
392  m_cInjectionLERTOP->Clear();
393  m_cInjectionLERTOP->cd(0);
394  m_hInjectionLERTOP->Draw("hist");
395 
396  locationHits = "TOPOccInjHER";
397  if (m_histogramDirectoryName != "") {
398  locationHits = m_histogramDirectoryName + "/" + locationHits;
399  }
400  Hits = (TH1*)findHist(locationHits.Data());
401  locationTriggers = "TOPEOccInjHER";
402  if (m_histogramDirectoryName != "") {
403  locationTriggers = m_histogramDirectoryName + "/" + locationTriggers;
404  }
405  Triggers = (TH1*)findHist(locationTriggers.Data());
406 
407  if (Hits && Triggers) {
408  m_hInjectionHERTOP->Divide(Hits, Triggers, 100, 8192);
409  }
410 
411  m_cInjectionHERTOP->Clear();
412  m_cInjectionHERTOP->cd(0);
413  m_hInjectionHERTOP->Draw("hist");
414 
415 
416 
417  //ARICH
418  m_histogramDirectoryName = "ARICH";
419 
420  locationHits = "ARICHOccInjLER";
421  if (m_histogramDirectoryName != "") {
422  locationHits = m_histogramDirectoryName + "/" + locationHits;
423  }
424  Hits = (TH1*)findHist(locationHits.Data());
425  locationTriggers = "ARICHEOccInjLER";
426  if (m_histogramDirectoryName != "") {
427  locationTriggers = m_histogramDirectoryName + "/" + locationTriggers;
428  }
429  Triggers = (TH1*)findHist(locationTriggers.Data());
430 
431  if (Hits && Triggers) {
432  m_hInjectionLERARICH->Divide(Hits, Triggers);
433  }
434 
435  m_cInjectionLERARICH->Clear();
436  m_cInjectionLERARICH->cd(0);
437  m_hInjectionLERARICH->Draw("hist");
438 
439  locationHits = "ARICHOccInjHER";
440  if (m_histogramDirectoryName != "") {
441  locationHits = m_histogramDirectoryName + "/" + locationHits;
442  }
443  Hits = (TH1*)findHist(locationHits.Data());
444  locationTriggers = "ARICHEOccInjHER";
445  if (m_histogramDirectoryName != "") {
446  locationTriggers = m_histogramDirectoryName + "/" + locationTriggers;
447  }
448  Triggers = (TH1*)findHist(locationTriggers.Data());
449 
450  if (Hits && Triggers) {
451  m_hInjectionHERARICH->Divide(Hits, Triggers);
452  }
453 
454  m_cInjectionHERARICH->Clear();
455  m_cInjectionHERARICH->cd(0);
456  m_hInjectionHERARICH->Draw("hist");
457 
458  // KLM
459  m_histogramDirectoryName = "KLM";
460 
461  locationHits = "KLMOccInjLER";
462  if (m_histogramDirectoryName != "") {
463  locationHits = m_histogramDirectoryName + "/" + locationHits;
464  }
465  Hits = (TH1*)findHist(locationHits.Data());
466  locationTriggers = "KLMTrigInjLER";
467  if (m_histogramDirectoryName != "") {
468  locationTriggers = m_histogramDirectoryName + "/" + locationTriggers;
469  }
470  Triggers = (TH1*)findHist(locationTriggers.Data());
471 
472  if (Hits && Triggers) {
473  m_hInjectionLERKLM->Divide(Hits, Triggers, 100, KLMElementNumbers::getTotalChannelNumber());
474  }
475 
476  m_cInjectionLERKLM->Clear();
477  m_cInjectionLERKLM->cd(0);
478  m_hInjectionLERKLM->Draw("hist");
479 
480  locationHits = "KLMOccInjHER";
481  if (m_histogramDirectoryName != "") {
482  locationHits = m_histogramDirectoryName + "/" + locationHits;
483  }
484  Hits = (TH1*)findHist(locationHits.Data());
485  locationTriggers = "KLMTrigInjHER";
486  if (m_histogramDirectoryName != "") {
487  locationTriggers = m_histogramDirectoryName + "/" + locationTriggers;
488  }
489  Triggers = (TH1*)findHist(locationTriggers.Data());
490 
491  if (Hits && Triggers) {
492  m_hInjectionHERKLM->Divide(Hits, Triggers, 100, KLMElementNumbers::getTotalChannelNumber());
493  }
494 
495  m_cInjectionHERKLM->Clear();
496  m_cInjectionHERKLM->cd(0);
497  m_hInjectionHERKLM->Draw("hist");
498 
499 #ifdef _BELLE2_EPICS
500  if (m_useEpics) {
501  for (auto& m : m_nodes) {
502  if (!m.mychid) continue;
503  int length = m.data.size();
504  if (length != int(ca_element_count(m.mychid)) && int(ca_element_count(m.mychid)) > 0) {
505  // FIXME, unclear why this is needed to prevent crashes on new run?
506  m.data.resize(int(ca_element_count(m.mychid)), 0.0);
507  length = m.data.size();
508  }
509  if (m.histo && m.histo->GetNcells() > 2 && length > 0 && length == int(ca_element_count(m.mychid))) {
510  // If bin count doesnt match, we loose bins but otherwise ca_array_put will complain
511  // We fill up the array with ZEROs otherwise
512  if (m.histo->GetDimension() == 1) {
513  int i = 0;
514  int nx = m.histo->GetNbinsX() + 1;
515  for (int x = 1; x < nx && i < length ; x++) {
516  m.data[i++] = m.histo->GetBinContent(x);
517  }
518 
519  } else if (m.histo->GetDimension() == 2) {
520  int i = 0;
521  int nx = m.histo->GetNbinsX() + 1;
522  int ny = m.histo->GetNbinsY() + 1;
523  for (int y = 1; y < ny && i < length; y++) {
524  for (int x = 1; x < nx && i < length ; x++) {
525  m.data[i++] = m.histo->GetBinContent(x, y);
526  }
527  }
528  }
529  SEVCHK(ca_array_put(DBR_DOUBLE, length, m.mychid, (void*)m.data.data()), "ca_put failure");
530  } else {
531  B2DEBUG(99, "Inj " << ca_name(m.mychid) << " , " << m.histo << " , " << (m.histo ? m.histo->GetNcells() : 0) << " , " << length <<
532  " , "
533  <<
534  ca_element_count(m.mychid));
535  }
536  }
537  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
538  }
539 #endif
540 }
541 
542 void DQMHistInjectionModule::cleanPVs(void)
543 {
544 #ifdef _BELLE2_EPICS
545  if (m_useEpics) {
546  for (auto m : m_nodes) {
547  if (m.mychid) {
548  int length = int(ca_element_count(m.mychid));
549  if (length > 0) {
550  m.data.resize(length, 0.0);
551  SEVCHK(ca_array_put(DBR_DOUBLE, length, m.mychid, (void*)m.data.data()), "ca_put failure");
552  } else {
553  B2DEBUG(99, "clean: lenght " << ca_name(m.mychid));
554  }
555  } else {
556  B2DEBUG(99, "clean: chid " << ca_name(m.mychid));
557  }
558  }
559  }
560 #endif
561 }
562 
563 void DQMHistInjectionModule::terminate()
564 {
565  B2DEBUG(1, "DQMHistInjection: terminate called");
566 #ifdef _BELLE2_EPICS
567  if (m_useEpics) {
568  for (auto m : m_nodes) {
569  SEVCHK(ca_clear_channel(m.mychid), "ca_clear_channel failure");
570  }
571  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
572  }
573 #endif
574 }
575 
The base class for the histogram analysis module.
DQM Histogram Analysis for PXD Efficiency.
#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.