Belle II Software  release-05-02-19
DQMHistAnalysisPXDCharge.cc
1 //+
2 // File : DQMHistAnalysisPXDCharge.cc
3 // Description : Analysis of PXD Cluster Charge
4 //
5 // Author : Bjoern Spruck, University Mainz
6 // Date : 2018
7 //-
8 
9 
10 #include <dqm/analysis/modules/DQMHistAnalysisPXDCharge.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(DQMHistAnalysisPXDCharge)
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("PXDER"));
34  addParam("RangeLow", m_rangeLow, "Lower boarder for fit", 30.);
35  addParam("RangeHigh", m_rangeHigh, "High border for fit", 85.);
36  addParam("PVPrefix", m_pvPrefix, "PV Prefix", std::string("DQM:PXD:Charge:"));
37  addParam("useEpics", m_useEpics, "useEpics", true);
38  B2DEBUG(99, "DQMHistAnalysisPXDCharge: Constructor done.");
39 }
40 
41 DQMHistAnalysisPXDChargeModule::~DQMHistAnalysisPXDChargeModule()
42 {
43 #ifdef _BELLE2_EPICS
44  if (m_useEpics) {
45  if (ca_current_context()) ca_context_destroy();
46  }
47 #endif
48 }
49 
50 void DQMHistAnalysisPXDChargeModule::initialize()
51 {
52  B2DEBUG(99, "DQMHistAnalysisPXDCharge: initialized.");
53 
54  m_monObj = getMonitoringObject("pxd");
55  const VXD::GeoCache& geo = VXD::GeoCache::getInstance();
56 
57  //collect the list of all PXD Modules in the geometry here
58  std::vector<VxdID> sensors = geo.getListOfSensors();
59  for (VxdID& aVxdID : sensors) {
60  VXD::SensorInfoBase info = geo.getSensorInfo(aVxdID);
61  if (info.getType() != VXD::SensorInfoBase::PXD) continue;
62  m_PXDModules.push_back(aVxdID); // reorder, sort would be better
63  }
64  std::sort(m_PXDModules.begin(), m_PXDModules.end()); // back to natural order
65 
66  gROOT->cd(); // this seems to be important, or strange things happen
67 
68  m_cCharge = new TCanvas((m_histogramDirectoryName + "/c_Charge").data());
69  m_hCharge = new TH1F("hPXDClusterCharge", "PXD Cluster Charge; Module; Track Cluster Charge", m_PXDModules.size(), 0,
70  m_PXDModules.size());
71  m_hCharge->SetDirectory(0);// dont mess with it, this is MY histogram
72  m_hCharge->SetStats(false);
73  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
74  TString ModuleName = (std::string)m_PXDModules[i];
75  m_hCharge->GetXaxis()->SetBinLabel(i + 1, ModuleName);
76  }
77  //Unfortunately this only changes the labels, but can't fill the bins by the VxdIDs
78  m_hCharge->Draw("");
79 
80  m_monObj->addCanvas(m_cCharge);
81 
83 // m_line1 = new TLine(0, 10, m_PXDModules.size(), 10);
84 // m_line2 = new TLine(0, 16, m_PXDModules.size(), 16);
85 // m_line3 = new TLine(0, 3, m_PXDModules.size(), 3);
86 // m_line1->SetHorizontal(true);
87 // m_line1->SetLineColor(3);// Green
88 // m_line1->SetLineWidth(3);
89 // m_line2->SetHorizontal(true);
90 // m_line2->SetLineColor(1);// Black
91 // m_line2->SetLineWidth(3);
92 // m_line3->SetHorizontal(true);
93 // m_line3->SetLineColor(1);
94 // m_line3->SetLineWidth(3);
95 
97 
98  m_fLandau = new TF1("f_Landau", "landau", m_rangeLow, m_rangeHigh);
99  m_fLandau->SetParameter(0, 1000);
100  m_fLandau->SetParameter(1, 50);
101  m_fLandau->SetParameter(2, 10);
102  m_fLandau->SetLineColor(4);
103  m_fLandau->SetNpx(60);
104  m_fLandau->SetNumberFitPoints(60);
105 
106  m_fMean = new TF1("f_Mean", "pol0", 0.5, m_PXDModules.size() - 0.5);
107  m_fMean->SetParameter(0, 50);
108  m_fMean->SetLineColor(5);
109  m_fMean->SetNpx(m_PXDModules.size());
110  m_fMean->SetNumberFitPoints(m_PXDModules.size());
111 
112 #ifdef _BELLE2_EPICS
113  if (m_useEpics) {
114  if (!ca_current_context()) SEVCHK(ca_context_create(ca_disable_preemptive_callback), "ca_context_create");
115  mychid.resize(3);
116  SEVCHK(ca_create_channel((m_pvPrefix + "Mean").data(), NULL, NULL, 10, &mychid[0]), "ca_create_channel failure");
117  SEVCHK(ca_create_channel((m_pvPrefix + "Diff").data(), NULL, NULL, 10, &mychid[1]), "ca_create_channel failure");
118  SEVCHK(ca_create_channel((m_pvPrefix + "Status").data(), NULL, NULL, 10, &mychid[2]), "ca_create_channel failure");
119  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
120  }
121 #endif
122 }
123 
124 
125 void DQMHistAnalysisPXDChargeModule::beginRun()
126 {
127  B2DEBUG(99, "DQMHistAnalysisPXDCharge: beginRun called.");
128 
129  m_cCharge->Clear();
130 }
131 
132 void DQMHistAnalysisPXDChargeModule::event()
133 {
134  if (!m_cCharge) return;
135  m_hCharge->Reset(); // dont sum up!!!
136 
137  bool enough = false;
138 
139  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
140  std::string name = "DQMER_PXD_" + (std::string)m_PXDModules[i] + "_ClusterCharge";
141  std::replace(name.begin(), name.end(), '.', '_');
142 
143  TH1* hh1 = findHist(name);
144  if (hh1 == NULL) {
145  hh1 = findHist(m_histogramDirectoryName, name);
146  }
147  if (hh1) {
148 // B2INFO("Histo " << name << " found in mem");
150  m_fLandau->SetParameter(0, 1000);
151  m_fLandau->SetParameter(1, 50);
152  m_fLandau->SetParLimits(1, 10, 80);
153  m_fLandau->SetParameter(2, 10);
154  m_fLandau->SetParLimits(2, 1, 50);
155  if (hh1->GetEntries() > 100) {
156  hh1->Fit(m_fLandau, "R");
157  m_hCharge->SetBinContent(i + 1, m_fLandau->GetParameter(1));
158  m_hCharge->SetBinError(i + 1, m_fLandau->GetParError(1));
159  } else {
160  m_hCharge->SetBinContent(i + 1, 50);
161  m_hCharge->SetBinError(i + 1, 100);
162  }
163  // m_hCharge->SetBinError(i + 1, m_fLandau->GetParameter(2) * 0.25); // arbitrary scaling
164  // cout << m_fLandau->GetParameter(0) << " " << m_fLandau->GetParameter(1) << " " << m_fLandau->GetParameter(2) << endl;
165 
166  // m_hCharge->Fill(i, hh1->GetMean());
167  m_cCharge->cd();
168  hh1->Draw();
169  if (hh1->GetEntries() > 100) m_fLandau->Draw("same");
170 // m_cCharge->Print(str(format("cc_%d.pdf") % i).data());
171 
172  if (hh1->GetEntries() > 1000) enough = true;
173  }
174  }
175  m_cCharge->cd();
176 
177  double data = 0;
178  double diff = 0;
179  if (m_hCharge && enough) {
180  double currentMin, currentMax;
181  m_hCharge->Draw("");
182 // m_line1->Draw();
183 // m_line2->Draw();
184 // m_line3->Draw();
185  m_hCharge->Fit(m_fMean, "R");
186  data = m_fMean->GetParameter(0); // we are more interessted in the maximum deviation from mean
187  m_hCharge->GetMinimumAndMaximum(currentMin, currentMax);
188  diff = fabs(data - currentMin) > fabs(currentMax - data) ? fabs(data - currentMin) : fabs(currentMax - data);
189  if (0) B2INFO("Mean: " << data << " Max Diff: " << diff);
190  }
191 
192 #ifdef _BELLE2_EPICS
193  if (m_useEpics) {
194  SEVCHK(ca_put(DBR_DOUBLE, mychid[0], (void*)&data), "ca_set failure");
195  SEVCHK(ca_put(DBR_DOUBLE, mychid[1], (void*)&diff), "ca_set failure");
196  }
197 #endif
198 
199  int status = 0;
200 
201  if (!enough) {
202  // not enough Entries
203  m_cCharge->Pad()->SetFillColor(kGray);// Magenta or Gray
204  status = 0; // default
205  } else {
207  if (fabs(data - 50.) > 20. || diff > 30) {
208  m_cCharge->Pad()->SetFillColor(kRed);// Red
209  status = 4;
210  } else if (fabs(data - 50) > 15. || diff > 10) {
211  m_cCharge->Pad()->SetFillColor(kYellow);// Yellow
212  status = 3;
213  } else {
214  m_cCharge->Pad()->SetFillColor(kGreen);// Green
215  status = 2;
216  }
217 
218  // FIXME overwrite for now
219  m_cCharge->Pad()->SetFillColor(kGreen);// Green
220  }
221 
222 #ifdef _BELLE2_EPICS
223  if (m_useEpics) {
224  SEVCHK(ca_put(DBR_INT, mychid[2], (void*)&status), "ca_set failure");
225  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
226  }
227 #endif
228 
229  auto tt = new TLatex(5.5, 0, "1.3.2 Module is excluded, please ignore");
230  tt->SetTextAngle(90);// Rotated
231  tt->SetTextAlign(12);// Centered
232  tt->Draw();
233 
234  m_cCharge->Modified();
235  m_cCharge->Update();
236 }
237 
238 void DQMHistAnalysisPXDChargeModule::endRun()
239 {
240  B2DEBUG(99, "DQMHistAnalysisPXDCharge : endRun called");
241 }
242 
243 
244 void DQMHistAnalysisPXDChargeModule::terminate()
245 {
246  // should delete canvas here, maybe hist, too? Who owns it?
247 #ifdef _BELLE2_EPICS
248  if (m_useEpics) {
249  for (auto m : mychid) SEVCHK(ca_clear_channel(m), "ca_clear_channel failure");
250  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
251  }
252 #endif
253  B2DEBUG(99, "DQMHistAnalysisPXDCharge: terminate called");
254 }
255 
Belle2::VxdID
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:43
Belle2::DQMHistAnalysisPXDChargeModule
DQM Histogram Analysis for PXD Cluster Charge.
Definition: DQMHistAnalysisPXDCharge.h:30
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