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