14 #include <dqm/analysis/modules/DQMHistAnalysisPXDCharge.h>
17 #include <vxd/geometry/GeoCache.h>
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,
"Whether to update EPICS PVs.",
false);
42 B2DEBUG(99,
"DQMHistAnalysisPXDCharge: Constructor done.");
45 DQMHistAnalysisPXDChargeModule::~DQMHistAnalysisPXDChargeModule()
49 if (ca_current_context()) ca_context_destroy();
54 void DQMHistAnalysisPXDChargeModule::initialize()
56 B2DEBUG(99,
"DQMHistAnalysisPXDCharge: initialized.");
58 m_monObj = getMonitoringObject(
"pxd");
63 for (
VxdID& aVxdID : sensors) {
65 if (info.getType() != VXD::SensorInfoBase::PXD)
continue;
66 m_PXDModules.push_back(aVxdID);
68 std::sort(m_PXDModules.begin(), m_PXDModules.end());
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,
75 m_hCharge->SetDirectory(0);
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);
84 m_monObj->addCanvas(m_cCharge);
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);
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());
118 if (!ca_current_context()) SEVCHK(ca_context_create(ca_disable_preemptive_callback),
"ca_context_create");
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");
129 void DQMHistAnalysisPXDChargeModule::beginRun()
131 B2DEBUG(99,
"DQMHistAnalysisPXDCharge: beginRun called.");
136 void DQMHistAnalysisPXDChargeModule::event()
138 if (!m_cCharge)
return;
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(),
'.',
'_');
147 TH1* hh1 = findHist(name);
149 hh1 = findHist(m_histogramDirectoryName, name);
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));
164 m_hCharge->SetBinContent(i + 1, 50);
165 m_hCharge->SetBinError(i + 1, 100);
173 if (hh1->GetEntries() > 100) m_fLandau->Draw(
"same");
176 if (hh1->GetEntries() > 1000) enough =
true;
183 if (m_hCharge && enough) {
184 double currentMin, currentMax;
189 m_hCharge->Fit(m_fMean,
"R");
190 data = m_fMean->GetParameter(0);
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);
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");
207 m_cCharge->Pad()->SetFillColor(kGray);
211 if (fabs(data - 50.) > 20. || diff > 30) {
212 m_cCharge->Pad()->SetFillColor(kRed);
214 }
else if (fabs(data - 50) > 15. || diff > 10) {
215 m_cCharge->Pad()->SetFillColor(kYellow);
218 m_cCharge->Pad()->SetFillColor(kGreen);
223 m_cCharge->Pad()->SetFillColor(kGreen);
228 SEVCHK(ca_put(DBR_INT, mychid[2], (
void*)&status),
"ca_set failure");
229 SEVCHK(ca_pend_io(5.0),
"ca_pend_io failure");
233 auto tt =
new TLatex(5.5, 0,
"1.3.2 Module is excluded, please ignore");
234 tt->SetTextAngle(90);
235 tt->SetTextAlign(12);
238 m_cCharge->Modified();
242 void DQMHistAnalysisPXDChargeModule::endRun()
244 B2DEBUG(99,
"DQMHistAnalysisPXDCharge : endRun called");
248 void DQMHistAnalysisPXDChargeModule::terminate()
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");
257 B2DEBUG(99,
"DQMHistAnalysisPXDCharge: terminate called");
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...
const std::vector< VxdID > getListOfSensors() const
Get list of all sensors.
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Base class to provide Sensor Information for PXD and SVD.
Class to uniquely identify a any structure of the PXD and SVD.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Abstract base class for different kinds of events.