Belle II Software  release-08-01-10
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 
31 DQMHistAnalysisPXDChargeModule::DQMHistAnalysisPXDChargeModule()
33 {
34  // This module CAN NOT be run in parallel!
35  setDescription("DQM Analysis for PXD Cluster Charge");
36 
37  // Parameter definition
38  addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of Histogram dir", std::string("PXDER"));
39  addParam("RangeLow", m_rangeLow, "Lower boarder for fit", 30.);
40  addParam("RangeHigh", m_rangeHigh, "High border for fit", 85.);
41  addParam("excluded", m_excluded, "excluded module (indizes starting from 0 to 39)");
42  B2DEBUG(99, "DQMHistAnalysisPXDCharge: Constructor done.");
43 }
44 
46 {
47 }
48 
50 {
51  B2DEBUG(99, "DQMHistAnalysisPXDCharge: initialized.");
52 
55 
56  //collect the list of all PXD Modules in the geometry here
57  std::vector<VxdID> sensors = geo.getListOfSensors();
58  for (VxdID& aVxdID : sensors) {
59  VXD::SensorInfoBase info = geo.getSensorInfo(aVxdID);
60  if (info.getType() != VXD::SensorInfoBase::PXD) continue;
61  m_PXDModules.push_back(aVxdID); // reorder, sort would be better
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_cCharge = new TCanvas((m_histogramDirectoryName + "/c_Charge").data());
68  m_hCharge = new TH1F("hPXDClusterCharge", "PXD Cluster Charge; Module; Track Cluster Charge", m_PXDModules.size(), 0,
69  m_PXDModules.size());
70  m_hCharge->SetDirectory(0);// dont mess with it, this is MY histogram
71  m_hCharge->SetStats(false);
72  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
73  TString ModuleName = (std::string)m_PXDModules[i];
74  m_hCharge->GetXaxis()->SetBinLabel(i + 1, ModuleName);
75  }
76  //Unfortunately this only changes the labels, but can't fill the bins by the VxdIDs
77  m_hCharge->Draw("");
78 
80 
82 // m_line1 = new TLine(0, 10, m_PXDModules.size(), 10);
83 // m_line2 = new TLine(0, 16, m_PXDModules.size(), 16);
84 // m_line3 = new TLine(0, 3, m_PXDModules.size(), 3);
85 // m_line1->SetHorizontal(true);
86 // m_line1->SetLineColor(3);// Green
87 // m_line1->SetLineWidth(3);
88 // m_line2->SetHorizontal(true);
89 // m_line2->SetLineColor(1);// Black
90 // m_line2->SetLineWidth(3);
91 // m_line3->SetHorizontal(true);
92 // m_line3->SetLineColor(1);
93 // m_line3->SetLineWidth(3);
94 
96 
97  m_fLandau = new TF1("f_Landau", "landau", m_rangeLow, m_rangeHigh);
98  m_fLandau->SetParameter(0, 1000);
99  m_fLandau->SetParameter(1, 50);
100  m_fLandau->SetParameter(2, 10);
101  m_fLandau->SetLineColor(4);
102  m_fLandau->SetNpx(60);
103  m_fLandau->SetNumberFitPoints(60);
104 
105  m_fMean = new TF1("f_Mean", "pol0", 0.5, m_PXDModules.size() - 0.5);
106  m_fMean->SetParameter(0, 50);
107  m_fMean->SetLineColor(5);
108  m_fMean->SetNpx(m_PXDModules.size());
109  m_fMean->SetNumberFitPoints(m_PXDModules.size());
110 
111  registerEpicsPV("PXD:Charge:Mean", "Mean");
112  registerEpicsPV("PXD:Charge:Diff", "Diff");
113  registerEpicsPV("PXD:Charge:Status", "Status");
114 }
115 
116 
118 {
119  B2DEBUG(99, "DQMHistAnalysisPXDCharge: beginRun called.");
120 
121  m_cCharge->Clear();
122 }
123 
125 {
126  if (!m_cCharge) return;
127  m_hCharge->Reset(); // dont sum up!!!
128 
129  bool enough = false;
130 
131  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
132  std::string name = "DQMER_PXD_" + (std::string)m_PXDModules[i] + "_ClusterCharge";
133  std::replace(name.begin(), name.end(), '.', '_');
134 
135  TH1* hh1 = findHist(name);
136  if (hh1 == NULL) {
137  hh1 = findHist(m_histogramDirectoryName, name);
138  }
139  if (hh1) {
140 // B2INFO("Histo " << name << " found in mem");
142  m_fLandau->SetParameter(0, 1000);
143  m_fLandau->SetParameter(1, 50);
144  m_fLandau->SetParLimits(1, 10, 80);
145  m_fLandau->SetParameter(2, 10);
146  m_fLandau->SetParLimits(2, 1, 50);
147  if (hh1->GetEntries() > 100) {
148  hh1->Fit(m_fLandau, "R");
149  m_hCharge->SetBinContent(i + 1, m_fLandau->GetParameter(1));
150  m_hCharge->SetBinError(i + 1, m_fLandau->GetParError(1));
151  } else {
152  m_hCharge->SetBinContent(i + 1, 50);
153  m_hCharge->SetBinError(i + 1, 100);
154  }
155  // m_hCharge->SetBinError(i + 1, m_fLandau->GetParameter(2) * 0.25); // arbitrary scaling
156  // cout << m_fLandau->GetParameter(0) << " " << m_fLandau->GetParameter(1) << " " << m_fLandau->GetParameter(2) << endl;
157 
158  // m_hCharge->Fill(i, hh1->GetMean());
159  m_cCharge->cd();
160  hh1->Draw();
161  if (hh1->GetEntries() > 100) m_fLandau->Draw("same");
162 // m_cCharge->Print(str(format("cc_%d.pdf") % i).data());
163 
164  if (hh1->GetEntries() > 1000) enough = true;
165  }
166  }
167  m_cCharge->cd();
168 
169  double data = 0;
170  double diff = 0;
171  if (m_hCharge && enough) {
172  double currentMin, currentMax;
173  m_hCharge->Draw("");
174 // m_line1->Draw();
175 // m_line2->Draw();
176 // m_line3->Draw();
177  m_hCharge->Fit(m_fMean, "R");
178  data = m_fMean->GetParameter(0); // we are more interessted in the maximum deviation from mean
179  m_hCharge->GetMinimumAndMaximum(currentMin, currentMax);
180  diff = fabs(data - currentMin) > fabs(currentMax - data) ? fabs(data - currentMin) : fabs(currentMax - data);
181  if (0) B2INFO("Mean: " << data << " Max Diff: " << diff);
182  }
183 
184  setEpicsPV("Mean", data);
185  setEpicsPV("Diff", diff);
186 
187  int status = 0;
188 
189  if (!enough) {
190  // not enough Entries
191  m_cCharge->Pad()->SetFillColor(kGray);// Magenta or Gray
192  status = 0; // default
193  } else {
195  if (fabs(data - 50.) > 20. || diff > 30) {
196  m_cCharge->Pad()->SetFillColor(kRed);// Red
197  status = 4;
198  } else if (fabs(data - 50) > 15. || diff > 10) {
199  m_cCharge->Pad()->SetFillColor(kYellow);// Yellow
200  status = 3;
201  } else {
202  m_cCharge->Pad()->SetFillColor(kGreen);// Green
203  status = 2;
204  }
205 
206  // FIXME overwrite for now
207  m_cCharge->Pad()->SetFillColor(kGreen);// Green
208  }
209 
210  setEpicsPV("Status", status);
211 
212  for (auto& it : m_excluded) {
213  auto tt = new TLatex(it + 0.5, 0, (" " + std::string(m_PXDModules[it]) + " Module is excluded, please ignore").c_str());
214  tt->SetTextSize(0.035);
215  tt->SetTextAngle(90);// Rotated
216  tt->SetTextAlign(12);// Centered
217  tt->Draw();
218  }
219 
220  m_cCharge->Modified();
221  m_cCharge->Update();
223 }
224 
226 {
227  B2DEBUG(99, "DQMHistAnalysisPXDCharge : endRun called");
228 }
229 
230 
232 {
233  B2DEBUG(99, "DQMHistAnalysisPXDCharge: terminate called");
234 }
235 
The base class for the histogram analysis module.
int registerEpicsPV(std::string pvname, std::string keyname="", bool update_pvs=true)
EPICS related Functions.
static TH1 * findHist(const std::string &histname, bool onlyIfUpdated=false)
Get histogram from list (no other search).
void setEpicsPV(std::string keyname, double value)
Write value to a EPICS PV.
void UpdateCanvas(std::string name, bool updated=true)
Mark canvas as updated (or not)
static MonitoringObject * getMonitoringObject(const std::string &histname)
Get MonitoringObject with given name (new object is created if non-existing)
void terminate(void) override final
This method is called at the end of the event processing.
double m_rangeLow
fit range lo edge for landau
void initialize(void) override final
Initializer.
TH1F * m_hCharge
Histogram covering all modules.
void endRun(void) override final
This method is called if the current run ends.
MonitoringObject * m_monObj
Monitoring Object.
TF1 * m_fMean
Fit the Mean for all modules.
std::vector< VxdID > m_PXDModules
IDs of all PXD Modules to iterate over.
std::string m_histogramDirectoryName
name of histogram directory
double m_rangeHigh
fit range hi edge for landau
std::vector< int > m_excluded
Indizes of excluded PXD Modules.
TF1 * m_fLandau
only one fit function for all Landaus
void beginRun(void) override final
Called when entering a new run.
void event(void) override final
This method is called for each event.
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void addCanvas(TCanvas *canv)
Add Canvas to monitoring object.
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:59
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:67
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
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
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#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.