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