Belle II Software  release-06-01-15
DQMHistAnalysisPXDReduction.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 : DQMHistAnalysisPXDReduction.cc
10 // Description : Analysis of PXD Reduction
11 //-
12 
13 
14 #include <dqm/analysis/modules/DQMHistAnalysisPXDReduction.h>
15 #include <TROOT.h>
16 #include <vxd/geometry/GeoCache.h>
17 
18 using namespace std;
19 using namespace Belle2;
20 
21 //-----------------------------------------------------------------
22 // Register the Module
23 //-----------------------------------------------------------------
24 REG_MODULE(DQMHistAnalysisPXDReduction)
25 
26 //-----------------------------------------------------------------
27 // Implementation
28 //-----------------------------------------------------------------
29 
32 {
33  // This module CAN NOT be run in parallel!
34 
35  //Parameter definition
36  addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of Histogram dir", std::string("PXDDAQ"));
37  addParam("PVPrefix", m_pvPrefix, "PV Prefix", std::string("DQM:PXD:Red:"));
38  addParam("useEpics", m_useEpics, "Whether to update EPICS PVs.", false);
39  B2DEBUG(1, "DQMHistAnalysisPXDReduction: Constructor done.");
40 }
41 
42 DQMHistAnalysisPXDReductionModule::~DQMHistAnalysisPXDReductionModule()
43 {
44 #ifdef _BELLE2_EPICS
45  if (m_useEpics) {
46  if (ca_current_context()) ca_context_destroy();
47  }
48 #endif
49 }
50 
51 void DQMHistAnalysisPXDReductionModule::initialize()
52 {
53  B2DEBUG(1, "DQMHistAnalysisPXDReduction: initialized.");
54 
55  m_monObj = getMonitoringObject("pxd");
56  const VXD::GeoCache& geo = VXD::GeoCache::getInstance();
57 
58  //collect the list of all PXD Modules in the geometry here
59  std::vector<VxdID> sensors = geo.getListOfSensors();
60  for (VxdID& aVxdID : sensors) {
61  VXD::SensorInfoBase info = geo.getSensorInfo(aVxdID);
62  // B2DEBUG(20,"VXD " << aVxdID);
63  if (info.getType() != VXD::SensorInfoBase::PXD) continue;
64  m_PXDModules.push_back(aVxdID); // reorder, sort would be better
65 
66  }
67  std::sort(m_PXDModules.begin(), m_PXDModules.end()); // back to natural order
68 
69  gROOT->cd(); // this seems to be important, or strange things happen
70 
71  m_cReduction = new TCanvas((m_histogramDirectoryName + "/c_Reduction").data());
72  m_hReduction = new TH1F("hPXDReduction", "PXD Reduction; Module; Reduction", m_PXDModules.size(), 0, m_PXDModules.size());
73  m_hReduction->SetDirectory(0);// dont mess with it, this is MY histogram
74  m_hReduction->SetStats(false);
75  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
76  TString ModuleName = (std::string)m_PXDModules[i];
77  m_hReduction->GetXaxis()->SetBinLabel(i + 1, ModuleName);
78  }
79  //Unfortunately this only changes the labels, but can't fill the bins by the VxdIDs
80  m_hReduction->Draw("");
81  m_monObj->addCanvas(m_cReduction);
82 
84  m_line1 = new TLine(0, 10, m_PXDModules.size(), 10);
85 // m_line2 = new TLine(0, 16, m_PXDModules.size(), 16);
86 // m_line3 = new TLine(0, 3, m_PXDModules.size(), 3);
87  m_line1->SetHorizontal(true);
88  m_line1->SetLineColor(3);// Green
89  m_line1->SetLineWidth(3);
90 // m_line2->SetHorizontal(true);
91 // m_line2->SetLineColor(1);// Black
92 // m_line2->SetLineWidth(3);
93 // m_line3->SetHorizontal(true);
94 // m_line3->SetLineColor(1);
95 // m_line3->SetLineWidth(3);
96 
97 
98 #ifdef _BELLE2_EPICS
99  if (m_useEpics) {
100  if (!ca_current_context()) SEVCHK(ca_context_create(ca_disable_preemptive_callback), "ca_context_create");
101  mychid.resize(2);
102  SEVCHK(ca_create_channel((m_pvPrefix + "Status").data(), NULL, NULL, 10, &mychid[0]), "ca_create_channel failure");
103  SEVCHK(ca_create_channel((m_pvPrefix + "Value").data(), NULL, NULL, 10, &mychid[1]), "ca_create_channel failure");
104  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
105  }
106 #endif
107 }
108 
109 
110 void DQMHistAnalysisPXDReductionModule::beginRun()
111 {
112  B2DEBUG(1, "DQMHistAnalysisPXDReduction: beginRun called.");
113 
114  m_cReduction->Clear();
115 }
116 
117 void DQMHistAnalysisPXDReductionModule::event()
118 {
119  if (!m_cReduction) return;
120  m_hReduction->Reset(); // dont sum up!!!
121 
122  bool enough = false;
123  double ireduction = 0.0;
124  int ireductioncnt = 0;
125 
126  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
127  std::string name = "PXDDAQDHEDataReduction_" + (std::string)m_PXDModules[i ];
128  // std::replace( name.begin(), name.end(), '.', '_');
129 
130  TH1* hh1 = findHist(name);
131  if (hh1 == NULL) {
132  hh1 = findHist(m_histogramDirectoryName, name);
133  }
134  if (hh1) {
135  auto mean = hh1->GetMean();
136  m_hReduction->Fill(i, mean);
137  if (hh1->GetEntries() > 100) enough = true;
138  if (mean > 0) {
139  ireduction += mean; // well fit would be better
140  ireductioncnt++;
141  }
142  }
143  }
144  m_cReduction->cd();
145 
146  int status = 0;
147  // not enough Entries
148  if (!enough) {
149  status = 0; // Grey
150  m_cReduction->Pad()->SetFillColor(kGray);// Magenta or Gray
151  } else {
152  status = 1; // White
154 // if (value > m_up_err_limit || value < m_low_err_limit ) {
155 // m_cReduction->Pad()->SetFillColor(kRed);// Red
156 // } else if (value > m_up_warn_limit || value < m_low_warn_limit ) {
157 // m_cReduction->Pad()->SetFillColor(kYellow);// Yellow
158 // } else {
159 // m_cReduction->Pad()->SetFillColor(kGreen);// Green
160 // } else {
161  m_cReduction->Pad()->SetFillColor(kWhite);// White
162 // }
163  }
164 
165  double value = ireductioncnt > 0 ? ireduction / ireductioncnt : 0;
166 
167  if (m_hReduction) {
168  m_hReduction->Draw("");
169  if (status != 0) {
170  m_line1->SetY1(value);
171  m_line1->SetY2(value); // aka SetHorizontal
172  m_line1->Draw();
173  }
174 // m_line2->Draw();
175 // m_line3->Draw();
176  }
177 
178  m_monObj->setVariable("reduction", value);
179 
180  m_cReduction->Modified();
181  m_cReduction->Update();
182 #ifdef _BELLE2_EPICS
183  if (m_useEpics) {
184  SEVCHK(ca_put(DBR_INT, mychid[0], (void*)&status), "ca_set failure");
185  SEVCHK(ca_put(DBR_DOUBLE, mychid[1], (void*)&value), "ca_set failure");
186  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
187  }
188 #endif
189 }
190 
191 void DQMHistAnalysisPXDReductionModule::terminate()
192 {
193  B2DEBUG(1, "DQMHistAnalysisPXDReduction: terminate called");
194  // m_cReduction->Print("c1.pdf");
195  // should delete canvas here, maybe hist, too? Who owns it?
196 #ifdef _BELLE2_EPICS
197  if (m_useEpics) {
198  for (auto m : mychid) SEVCHK(ca_clear_channel(m), "ca_clear_channel failure");
199  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
200  }
201 #endif
202 }
203 
The base class for the histogram analysis module.
DQM Histogram Analysis for PXD Reduction.
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.