Belle II Software  release-06-02-00
DQMHistAnalysisPXDCM.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 : DQMHistAnalysisPXDCM.cc
10 // Description : Analysis of PXD Common Modes
11 //-
12 
13 
14 #include <dqm/analysis/modules/DQMHistAnalysisPXDCM.h>
15 #include <TROOT.h>
16 #include <TLatex.h>
17 #include <TPaveText.h>
18 #include <vxd/geometry/GeoCache.h>
19 #include <framework/core/ModuleParam.templateDetails.h>
20 
21 using namespace std;
22 using namespace Belle2;
23 
24 //-----------------------------------------------------------------
25 // Register the Module
26 //-----------------------------------------------------------------
27 REG_MODULE(DQMHistAnalysisPXDCM)
28 
29 //-----------------------------------------------------------------
30 // Implementation
31 //-----------------------------------------------------------------
32 
35 {
36  // This module CAN NOT be run in parallel!
37 
38  //Parameter definition
39  addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of Histogram dir", std::string("PXDCM"));
40  addParam("PVPrefix", m_pvPrefix, "PV Prefix", std::string("DQM:PXD:CommonMode:"));
41  addParam("useEpics", m_useEpics, "Whether to update EPICS PVs.", false);
42  addParam("minEntries", m_minEntries, "minimum number of new entries for last time slot", 10000);
43 
44  addParam("warnMeanAdhoc", m_warnMeanAdhoc, "warn level for peak position", 2.0);
45  addParam("errorMeanAdhoc", m_errorMeanAdhoc, "error level for peak position", 3.0);
46  addParam("warnOutsideAdhoc", m_warnOutsideAdhoc, "warn level for outside fraction", 1e-5);
47  addParam("errorOutsideAdhoc", m_errorOutsideAdhoc, "error level for outside fraction", 1e-4);
48  addParam("upperLineAdhoc", m_upperLineAdhoc, "upper threshold and line for outside fraction", 17);
49 
50  addParam("warnMeanFull", m_warnMeanFull, "warn level for peak position", 2.0);
51  addParam("errorMeanFull", m_errorMeanFull, "error level for peak position", 3.0);
52  addParam("warnOutsideFull", m_warnOutsideFull, "warn level for outside fraction", 1e-5);
53  addParam("errorOutsideFull", m_errorOutsideFull, "error level for outside fraction", 1e-4);
54  addParam("upperLineFull", m_upperLineFull, "upper threshold and line for outside fraction", 17);
55 
56  addParam("gateMaskModuleList", m_par_module_list, "Module List for Gate Masking");
57  addParam("gateMaskGateList", m_par_gate_list, "Gate List for Gate Masking");
58 
59  B2DEBUG(99, "DQMHistAnalysisPXDCM: Constructor done.");
60 }
61 
62 DQMHistAnalysisPXDCMModule::~DQMHistAnalysisPXDCMModule()
63 {
64 #ifdef _BELLE2_EPICS
65  if (m_useEpics) {
66  if (ca_current_context()) ca_context_destroy();
67  }
68 #endif
69 }
70 
71 void DQMHistAnalysisPXDCMModule::initialize()
72 {
73  m_monObj = getMonitoringObject("pxd");
74  const VXD::GeoCache& geo = VXD::GeoCache::getInstance();
75 
76  // collect the list of all PXD Modules in the geometry here
77  std::vector<VxdID> sensors = geo.getListOfSensors();
78  for (VxdID& aVxdID : sensors) {
79  VXD::SensorInfoBase info = geo.getSensorInfo(aVxdID);
80  if (info.getType() != VXD::SensorInfoBase::PXD) continue;
81  m_PXDModules.push_back(aVxdID); // reorder, sort would be better
82  }
83  std::sort(m_PXDModules.begin(), m_PXDModules.end()); // back to natural order
84 
85  gROOT->cd(); // this seems to be important, or strange things happen
86 
87  m_cCommonMode = new TCanvas((m_histogramDirectoryName + "/c_CommonMode").data());
88  m_cCommonModeDelta = new TCanvas((m_histogramDirectoryName + "/c_CommonModeDelta").data());
89 
90  m_hCommonMode = new TH2D("hPXDCommonMode", "PXD CommonMode; Module; CommonMode", m_PXDModules.size(), 0, m_PXDModules.size(), 63, 0,
91  63);
92  m_hCommonMode->SetDirectory(0);// dont mess with it, this is MY histogram
93  m_hCommonMode->SetStats(false);
94  m_hCommonModeDelta = new TH2D("hPXDCommonModeAdhoc", "PXD CommonMode Adhoc; Module; CommonMode", m_PXDModules.size(), 0,
95  m_PXDModules.size(), 63, 0, 63);
96  m_hCommonModeDelta->SetDirectory(0);// dont mess with it, this is MY histogram
97  m_hCommonModeDelta->SetStats(false);
98  m_hCommonModeOld = new TH2D("hPXDCommonModeOld", "PXD CommonMode Old; Module; CommonMode", m_PXDModules.size(), 0,
99  m_PXDModules.size(),
100  63, 0, 63);
101  m_hCommonModeOld->SetDirectory(0);// dont mess with it, this is MY histogram
102  m_hCommonModeOld->SetStats(false);
103 
104  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
105  TString ModuleName = (std::string)m_PXDModules[i];
106  m_hCommonMode->GetXaxis()->SetBinLabel(i + 1, ModuleName);
107  m_hCommonModeDelta->GetXaxis()->SetBinLabel(i + 1, ModuleName);
108  m_hCommonModeOld->GetXaxis()->SetBinLabel(i + 1, ModuleName);
109  }
110  //Unfortunately this only changes the labels, but can't fill the bins by the VxdIDs
111  m_hCommonMode->Draw("colz");
112  m_hCommonModeDelta->Draw("colz");
113  m_hCommonModeOld->Draw("colz");
114 
115  m_monObj->addCanvas(m_cCommonMode);
116 
117  if (m_par_module_list.size() != m_par_gate_list.size()) {
118  B2FATAL("Parameter list need same length");
119  return;
120  }
121  for (size_t i = 0; i < m_par_module_list.size(); i++) {
122  for (auto n : m_par_gate_list[i]) {
123  m_masked_gates[VxdID(m_par_module_list[i])].push_back(n);
124  }
125  }
126 
128  m_line1 = new TLine(0, 10, m_PXDModules.size(), 10);
129  m_lineA = new TLine(0, m_upperLineAdhoc, m_PXDModules.size(), m_upperLineAdhoc);
130  m_lineF = new TLine(0, m_upperLineFull, m_PXDModules.size(), m_upperLineFull);
131  m_line1->SetHorizontal(true);
132  m_line1->SetLineColor(3);// Green
133  m_line1->SetLineWidth(3);
134  m_lineA->SetHorizontal(true);
135  m_lineA->SetLineColor(1);// Black
136  m_lineA->SetLineWidth(3);
137  m_lineF->SetHorizontal(true);
138  m_lineF->SetLineColor(1);// Black
139  m_lineF->SetLineWidth(3);
140 
141 
142 #ifdef _BELLE2_EPICS
143  if (m_useEpics) {
144  if (!ca_current_context()) SEVCHK(ca_context_create(ca_disable_preemptive_callback), "ca_context_create");
145  mychid.resize(4);
146  SEVCHK(ca_create_channel((m_pvPrefix + "Outside").data(), NULL, NULL, 10, &mychid[0]), "ca_create_channel failure");
147  SEVCHK(ca_create_channel((m_pvPrefix + "Status").data(), NULL, NULL, 10, &mychid[1]), "ca_create_channel failure");
148  SEVCHK(ca_create_channel((m_pvPrefix + "CM63").data(), NULL, NULL, 10, &mychid[2]), "ca_create_channel failure");
149  SEVCHK(ca_create_channel((m_pvPrefix + "Status_Adhoc").data(), NULL, NULL, 10, &mychid[3]), "ca_create_channel failure");
150 
151  for (VxdID& aPXDModule : m_PXDModules) {
152  TString buff = (std::string)aPXDModule;
153  buff.ReplaceAll(".", "_");
154  auto& my = mychid_mean[aPXDModule];
155  SEVCHK(ca_create_channel((m_pvPrefix + "Mean:" + buff).Data(), NULL, NULL, 10, &my), "ca_create_channel failure");
156  }
157  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
158  }
159 #endif
160  B2DEBUG(99, "DQMHistAnalysisPXDCM: initialized.");
161 }
162 
163 void DQMHistAnalysisPXDCMModule::beginRun()
164 {
165  B2DEBUG(99, "DQMHistAnalysisPXDCM: beginRun called.");
166 
167  m_cCommonMode->Clear();
168  m_cCommonModeDelta->Clear();
169  m_cCommonMode->SetLogz();
170  m_cCommonModeDelta->SetLogz();
171 
172  // this is needed at least for the "Old" and "Delta" one or update doesnt work
173  m_hCommonMode->Reset();
174  m_hCommonModeDelta->Reset();
175  m_hCommonModeOld->Reset();
176 }
177 
178 void DQMHistAnalysisPXDCMModule::event()
179 {
180  double all_outside = 0.0, all = 0.0;
181  double all_cm = 0.0;
182  bool error_full_flag = false;
183  bool warn_full_flag = false;
184  bool error_adhoc_flag = false;
185  bool warn_adhoc_flag = false;
186  bool anyupdate = false;
187  if (!m_cCommonMode) return;
188  m_hCommonMode->Reset(); // dont sum up!!!
189 
190  auto leg = new TPaveText(0.1, 0.6, 0.90, 0.95, "NDC");
191  leg->SetFillStyle(0);
192  leg->SetBorderSize(0);
193 
194  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
195  auto modname = (std::string)m_PXDModules[i];
196  std::string name = "PXDDAQCM_" + modname;
197  // std::replace( name.begin(), name.end(), '.', '_');
198 
199  TH1* hh1 = findHist(name);
200  if (hh1 == NULL) {
201  hh1 = findHist(m_histogramDirectoryName, name);
202  }
203  if (hh1) {
204  double current_full = 0.0;
205  double outside_full = 0.0;
206 
207  auto nevent = hh1->GetBinContent(0); // misuse underflow as event counter
208  double scale = nevent - m_hCommonModeOld->GetBinContent(i + 1, 0); // number of new events for delta
209  bool update = scale > m_minEntries ;
210  anyupdate |= update;
211  if (update) m_hCommonModeOld->SetBinContent(i + 1, 0, nevent);
212  if (scale > 0) scale = 1.0 / scale;
213  else scale = 1.; // worst case, no events at run start
214 
215  auto& gm = m_masked_gates[m_PXDModules[i]];
216  // We loop over a 2d histogram!
217  // loop CM values
218  double dhpc = 0.0;
219  for (int bin = 1; bin <= 64; bin++) { // including CM63!!!
220  // loop gates*asics
221  double v = 0;
222  for (int gate = 0; gate < 192; gate++) {
223  // attention, gate is not bin nr!
224  if (std::find(gm.begin(), gm.end(), gate) == gm.end()) {
225  v += hh1->GetBinContent(hh1->GetBin(gate + 1 + 192 * 0, bin)) +
226  hh1->GetBinContent(hh1->GetBin(gate + 1 + 192 * 1, bin)) +
227  hh1->GetBinContent(hh1->GetBin(gate + 1 + 192 * 2, bin)) +
228  hh1->GetBinContent(hh1->GetBin(gate + 1 + 192 * 3, bin));
229  }
230  }
231  m_hCommonMode->SetBinContent(i + 1, bin, v); // attention, mixing bin nr and index
232  // integration intervalls depend on CM default value, this seems to be agreed =10
233  // FIXME currently we have to much noise below the line ... thus excluding this to avoid false alarms
234  // outside_full += hh1->Integral(1 /*0*/, 5); /// FIXME we exclude bin 0 as we use it for debugging/timing pixels
235  // attention, n bins!
236  // we integrate up including value 62 (cm overflow), but not 63 (fifo full)
237  if (bin == 63 + 1) { // CM63
238  dhpc += v;
239  } else { // excluding CM63
240  current_full += v;
241  if (bin > m_upperLineFull + 1) outside_full += v;
242  }
243  if (nevent < m_minEntries) {
244  m_hCommonModeDelta->SetBinContent(i + 1, bin, v * scale); // attention, mixing bin nr and index
245  } else if (update) {
246  auto old = m_hCommonModeOld->GetBinContent(i + 1, bin); // attention, mixing bin nr and index
247  m_hCommonModeDelta->SetBinContent(i + 1, bin, (v - old)*scale); // attention, mixing bin nr and index
248  m_hCommonModeOld->SetBinContent(i + 1, bin, v); // attention, mixing bin nr and index
249  }
250  }
251 
252  all_outside += outside_full;
253  all += current_full;
254  all_cm += dhpc;
255  if (current_full > 1) {
256  error_full_flag |= (outside_full / current_full > m_errorOutsideFull);
257  warn_full_flag |= (outside_full / current_full > m_warnOutsideFull);
258  }
259 
260  if (update) {
261  Double_t mean_adhoc = 0.;
262  Double_t entries_adhoc = 0.;
263  Double_t outside_adhoc = 0.;
264 
265  // Attention, Bins
266  // we do not need to re-scale it as the scale is the same for all bins
267  for (int cm_y = 0; cm_y < m_upperLineAdhoc; cm_y++) {
268  auto v = m_hCommonModeDelta->GetBinContent(m_hCommonModeDelta->GetBin(i + 1, cm_y + 1));
269  entries_adhoc += v;
270  mean_adhoc += v * (cm_y + 1);
271  }
272  // Attention, Bins
273  // We ignore CM63 in outside and overall count
274  for (int cm_y = m_upperLineAdhoc; cm_y < 63; cm_y++) {
275  auto v = m_hCommonModeDelta->GetBinContent(m_hCommonModeDelta->GetBin(i + 1, cm_y + 1));
276  entries_adhoc += v;
277  outside_adhoc += v;
278  }
279  if (entries_adhoc > 0 && scale < 1e-3) { // ignore modules with minimum events
280  // scale <1e-3 == >1000 events
281  mean_adhoc /= entries_adhoc; // calculate mean
282  auto warn_tmp_m = fabs(10.0 - mean_adhoc) > m_warnMeanAdhoc;
283  auto err_tmp_m = fabs(10.0 - mean_adhoc) > m_errorMeanAdhoc;
284  auto warn_tmp_os = outside_adhoc / entries_adhoc > m_warnOutsideAdhoc;
285  auto err_tmp_os = outside_adhoc / entries_adhoc > m_errorOutsideAdhoc;
286  warn_adhoc_flag |= warn_tmp_m || warn_tmp_os;
287  error_adhoc_flag |= err_tmp_m || err_tmp_os;
288 
289  if (warn_tmp_m || err_tmp_m) {
290  TString tmp;
291  tmp.Form("%s: Mean %f", modname.c_str(), mean_adhoc);
292  leg->AddText(tmp);
293  B2INFO(name << " Mean " << mean_adhoc << " " << warn_tmp_m << err_tmp_m);
294  }
295  if (warn_tmp_os || err_tmp_os) {
296  TString tmp;
297  tmp.Form("%s: Outside %f %%", modname.c_str(), 100. * outside_adhoc / entries_adhoc);
298  leg->AddText(tmp);
299  B2INFO(name << " Outside " << outside_adhoc / entries_adhoc << " (" << outside_adhoc << "/" << entries_adhoc << ") " << warn_tmp_os
300  << err_tmp_os);
301  }
302  m_monObj->setVariable(("cm_" + modname).c_str(), mean_adhoc);
303 #ifdef _BELLE2_EPICS
304  if (m_useEpics) {
305  auto my = mychid_mean[m_PXDModules[i]];
306  if (my) SEVCHK(ca_put(DBR_DOUBLE, my, (void*)&mean_adhoc), "ca_set failure");
307  }
308 #endif
309  }
310  }
311  }
312  }
313 
314  int status = 0;
315  {
316  m_cCommonMode->cd();
317  // not enough Entries
318  if (all < 100.) {
319  m_cCommonMode->Pad()->SetFillColor(kGray);// Magenta or Gray
320  status = 0; // default
321  } else {
323  if (all_outside / all > m_errorOutsideFull || error_full_flag) {
324  m_cCommonMode->Pad()->SetFillColor(kRed);// Red
325  status = 4;
326  } else if (all_outside / all > m_warnOutsideFull || warn_full_flag) {
327  m_cCommonMode->Pad()->SetFillColor(kYellow);// Yellow
328  status = 3;
329  } else if (all_outside == 0. /*&& all_cm == 0.*/) {
330  // do not react on all_cm, we better monitor it elsewhere for clearity
331  m_cCommonMode->Pad()->SetFillColor(kGreen);// Green
332  status = 2;
333  } else { // between 0 and 50 ...
334  m_cCommonMode->Pad()->SetFillColor(kWhite);// White
335  status = 1;
336  }
337  }
338 
339  if (m_hCommonMode) {
340  m_hCommonMode->Draw("colz");
341  m_line1->Draw();
342  m_lineA->Draw();
343  }
344 
345  auto tt = new TLatex(5.5, 3, "1.3.2 Module is excluded, please ignore");
346  tt->SetTextAngle(90);// Rotated
347  tt->SetTextAlign(12);// Centered
348  tt->Draw();
349 
350  m_cCommonMode->Modified();
351  m_cCommonMode->Update();
352  }
353 
354  {
355  int status_adhoc = 0;
356  m_cCommonModeDelta->cd();
357  // not enough Entries
358 
359  if (all < 100.) { // delta cannot be more than all
360  m_cCommonModeDelta->Pad()->SetFillColor(kGray);// Magenta or Gray
361  status_adhoc = 0; // default
362  } else {
364  if (error_adhoc_flag) {
365  m_cCommonModeDelta->Pad()->SetFillColor(kRed);// Red
366  status_adhoc = 4;
367  } else if (warn_adhoc_flag) {
368  m_cCommonModeDelta->Pad()->SetFillColor(kYellow);// Yellow
369  status_adhoc = 3;
370  } else {
371  m_cCommonModeDelta->Pad()->SetFillColor(kGreen);// Green
372  status_adhoc = 2;
373  /* } else { // between 0 and 50 ...
374  m_cCommonModeDelta->Pad()->SetFillColor(kWhite);// White
375  status_adhoc = 1;*/
376  }
377  }
378 #ifdef _BELLE2_EPICS
379  if (m_useEpics && anyupdate) SEVCHK(ca_put(DBR_INT, mychid[3], (void*)&status_adhoc), "ca_set failure");
380 #endif
381  if (m_hCommonModeDelta) {
382  m_hCommonModeDelta->Draw("colz");
383  leg->Draw();
384  m_line1->Draw();
385  m_lineA->Draw();
386  }
387 
388  // keep this commented code as we may have excluded modules in phase4
389 // auto tt = new TLatex(5.5, 3, "1.3.2 Module is excluded, please ignore");
390 // tt->SetTextAngle(90);// Rotated
391 // tt->SetTextAlign(12);// Centered
392 // tt->Draw();
393 
394  m_cCommonModeDelta->Modified();
395  m_cCommonModeDelta->Update();
396  }
397 
398 #ifdef _BELLE2_EPICS
399  if (m_useEpics) {
400  double data = all > 0 ? (all_outside / all) : 0;
401  double data2 = all > 0 ? (all_cm / all) : 0;
402  SEVCHK(ca_put(DBR_DOUBLE, mychid[0], (void*)&data), "ca_set failure");
403  SEVCHK(ca_put(DBR_INT, mychid[1], (void*)&status), "ca_set failure");
404  SEVCHK(ca_put(DBR_DOUBLE, mychid[2], (void*)&data2), "ca_set failure");
405  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
406  }
407 #endif
408 }
409 
410 void DQMHistAnalysisPXDCMModule::terminate()
411 {
412  B2DEBUG(99, "DQMHistAnalysisPXDCM: terminate called");
413  // should delete canvas here, maybe hist, too? Who owns it?
414 #ifdef _BELLE2_EPICS
415  if (m_useEpics) {
416  for (auto m : mychid) SEVCHK(ca_clear_channel(m), "ca_clear_channel failure");
417  for (auto& m : mychid_mean) SEVCHK(ca_clear_channel(m.second), "ca_clear_channel failure");
418  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
419  }
420 #endif
421 }
422 
The base class for the histogram analysis module.
DQM Histogram Analysis for PXD Common Modes.
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:39
const std::vector< VxdID > getListOfSensors() const
Get list of all sensors.
Definition: GeoCache.cc:58
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:66
Base class to provide Sensor Information for PXD and SVD.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
#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.