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