Belle II Software development
DQMHistAnalysisPXDFits.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 : DQMHistAnalysisPXDFits.cc
10// Description : DQM module, which fits many PXD histograms and writes out fit parameters in new histograms
11//-
12
13
14#include <dqm/analysis/modules/DQMHistAnalysisPXDFits.h>
15#include <TROOT.h>
16
17#include <boost/format.hpp>
18
19using namespace std;
20using namespace Belle2;
21
22using boost::format;
23
24//-----------------------------------------------------------------
25// Register the Module
26//-----------------------------------------------------------------
27REG_MODULE(DQMHistAnalysisPXDFits);
28
29//-----------------------------------------------------------------
30// Implementation
31//-----------------------------------------------------------------
32
35{
36 // This module CAN NOT be run in parallel!
37 setDescription("DQM Analysis for PXD Raw");
38
39 // Parameter definition
40 addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms were placed",
41 std::string("PXDRAW"));
42// addParam("HistoName", m_histoname, "Name of Histogram (incl dir)", std::string(""));
43 for (auto i = 0, j = 0; i < 64; i++) {
44 auto layer = (((i >> 5) & 0x1) + 1);
45 auto ladder = ((i >> 1) & 0xF);
46 if (ladder == 0) continue; // numbering starts at 1
47 if (layer == 1 && ladder > 8) continue; // 8 inner ladders
48 if (layer == 2 && ladder > 12) continue; // 12 outer ladders
49 m_id_to_inx[i] = j; // i = id , j - index
50 m_inx_to_id[j] = i;
51 j++;
52 if (j == NUM_MODULES) break;
53 }
54 for (auto i = 0; i < NUM_MODULES; i++) {
55 m_hSignal[i] = nullptr;
56 m_hCommon[i] = nullptr;
57 m_hCounts[i] = nullptr;
58 m_cSignal[i] = nullptr;
59 m_cCommon[i] = nullptr;
60 m_cCounts[i] = nullptr;
61 }
62 B2DEBUG(1, "DQMHistAnalysisPXDFits: Constructor done.");
63}
64
66{
67 B2DEBUG(1, "DQMHistAnalysisPXDFits: initialized.");
68
69 gROOT->cd(); // this seems to be important, or strange things happen
70
71 TString a;
72
73 a = "pxdraw/hSignalAll";
74 a.ReplaceAll("/", "_");
75 m_cSignalAll = new TCanvas((m_histogramDirectoryName + "/c_").data() + a);
76 m_hSignalAll = new TH1F(a, a, NUM_MODULES, 0, NUM_MODULES);
77 m_hSignalAll->SetDirectory(0);// dont mess with it, this is MY histogram
78 m_hSignalAll->SetStats(false);
79
80 a = "pxdraw/hCommonAll";
81 a.ReplaceAll("/", "_");
82 m_cCommonAll = new TCanvas((m_histogramDirectoryName + "/c_").data() + a);
83 m_hCommonAll = new TH1F(a, a, NUM_MODULES, 0, NUM_MODULES);
84 m_hCommonAll->SetDirectory(0);// dont mess with it, this is MY histogram
85 m_hCommonAll->SetStats(false);
86
87 a = "pxdraw/hCountsAll";
88 a.ReplaceAll("/", "_");
89 m_cCountsAll = new TCanvas((m_histogramDirectoryName + "/c_").data() + a);
90 m_hCountsAll = new TH1F(a, a, NUM_MODULES, 0, NUM_MODULES);
91 m_hCountsAll->SetDirectory(0);// dont mess with it, this is MY histogram
92 m_hCountsAll->SetStats(false);
93
94 a = "pxdraw/hOccupancyAll";
95 a.ReplaceAll("/", "_");
96 m_cOccupancyAll = new TCanvas((m_histogramDirectoryName + "/c_").data() + a);
97 m_hOccupancyAll = new TH1F(a, a, NUM_MODULES, 0, NUM_MODULES);
98 m_hOccupancyAll->SetDirectory(0);// dont mess with it, this is MY histogram
99 m_hOccupancyAll->SetStats(false);
100
101 for (auto i = 0; i < NUM_MODULES; i++) {
102 auto id = m_inx_to_id[i];
103 auto layer = (((id >> 5) & 0x1) + 1);
104 auto ladder = ((id >> 1) & 0xF);
105 auto sensor = ((id & 0x1) + 1);
106 string s2 = str(format("_%d.%d.%d") % layer % ladder % sensor);
107
108 m_hSignalAll->GetXaxis()->SetBinLabel(i + 1, TString(s2));
109 m_hCommonAll->GetXaxis()->SetBinLabel(i + 1, TString(s2));
110 m_hCountsAll->GetXaxis()->SetBinLabel(i + 1, TString(s2));
111 m_hOccupancyAll->GetXaxis()->SetBinLabel(i + 1, TString(s2));
112
113 a = "pxdraw/hSignal";
114 a.ReplaceAll("/", "_");
115 a += s2;
116
117 m_cSignal[i] = new TCanvas((m_histogramDirectoryName + "/c_").data() + a);
118 m_hSignal[i] = new TH2F(a, a, 6, 0, 6, 4, 0, 4);
119 m_hSignal[i]->SetDirectory(0);// dont mess with it, this is MY histogram
120 m_hSignal[i]->SetStats(false);
121 m_hSignal[i]->SetMinimum(0);
122 m_hSignal[i]->SetMaximum(64);
123
124 a = "pxdraw/hCommon";
125 a.ReplaceAll("/", "_");
126 a += s2;
127 m_cCommon[i] = new TCanvas((m_histogramDirectoryName + "/c_").data() + a);
128 m_hCommon[i] = new TH2F(a, a, 6, 0, 6, 4, 0, 4);
129 m_hCommon[i]->SetDirectory(0);// dont mess with it, this is MY histogram
130 m_hCommon[i]->SetStats(false);
131 m_hCommon[i]->SetMinimum(0);
132 m_hCommon[i]->SetMaximum(256);
133
134 a = "pxdraw/hCounts";
135 a.ReplaceAll("/", "_");
136 a += s2;
137 m_cCounts[i] = new TCanvas((m_histogramDirectoryName + "/c_").data() + a);
138 m_hCounts[i] = new TH2F(a, a, 6, 0, 6, 4, 0, 4);
139 m_hCounts[i]->SetDirectory(0);// dont mess with it, this is MY histogram
140 m_hCounts[i]->SetStats(false);
141 }
142
143 m_fLandau = new TF1("f_Landau", "landau", 0, 256);
144 m_fLandau->SetParameter(0, 1000);
145 m_fLandau->SetParameter(1, 0);
146 m_fLandau->SetParameter(2, 10);
147 m_fLandau->SetLineColor(4);
148 m_fLandau->SetNpx(256);
149 m_fLandau->SetNumberFitPoints(256);
150
151 m_fGaus = new TF1("f_Gaus", "gaus", 0, 8096);
152 m_fGaus->SetParameter(0, 1000);
153 m_fGaus->SetParameter(1, 0);
154 m_fGaus->SetParameter(2, 10);
155 m_fGaus->SetLineColor(4);
156 m_fGaus->SetNpx(256);
157 m_fGaus->SetNumberFitPoints(256);
158}
159
160
162{
163 B2DEBUG(1, "DQMHistAnalysisPXDFits: beginRun called.");
164
165 // not much we can do here ... as we have created everything already
166 for (auto i = 0; i < NUM_MODULES; i++) {
167 m_cSignal[i]->Clear();
168 m_cCommon[i]->Clear();
169 m_cCounts[i]->Clear();
170 // no need to Cd and Draw yet
171 }
172}
173
175{
176// bool flag = false;
177
178 m_hSignalAll->Reset(); // dont sum up!!!
179 m_hCommonAll->Reset(); // dont sum up!!!
180 m_hCountsAll->Reset(); // dont sum up!!!
181 m_hOccupancyAll->Reset(); // dont sum up!!!
182
183 for (auto i = 0; i < NUM_MODULES; i++) {
184 auto id = m_inx_to_id[i];
185 auto layer = (((id >> 5) & 0x1) + 1);
186 auto ladder = ((id >> 1) & 0xF);
187 auto sensor = ((id & 0x1) + 1);
188
189 m_hSignal[i]->Reset(); // dont sum up!!!
190 m_hCommon[i]->Reset(); // dont sum up!!!
191 m_hCounts[i]->Reset(); // dont sum up!!!
192
193 for (auto j = 0; j < 6; j++) {
194 for (auto k = 0; k < 4; k++) {
195 //TH1* hh1 = NULL;
196 string s2 = str(format("_%d.%d.%d_%d_%d") % layer % ladder % sensor % j % k);
197
198 std::string name = "hrawPxdHitsCharge" + s2;
199 TH1* hh1 = findHist(name);
200 if (hh1 == NULL) {
202 }
203
204 if (hh1 != NULL) {
205// cout << "do da fit " << endl;
206// m_fLandau->SetParameter(0, 1000);
207// m_fLandau->SetParameter(1, 0);
208// m_fLandau->SetParameter(2, 10);
209// hh1->Fit(m_fLandau, "0");
210// m_hSignal[i]->Fill(j, k, m_fLandau->GetParameter(1));
211// cout << m_fLandau->GetParameter(0) << " " << m_fLandau->GetParameter(1) << " " << m_fLandau->GetParameter(2) << endl;
212 m_hSignal[i]->Fill(j, k, hh1->GetMean());
213 m_hSignalAll->Fill(i, hh1->GetMean());
214 } else {
215 B2INFO("Histo " << name << " not found");
216 }
217
218 name = "hrawPxdHitsCommonMode" + s2;
219 hh1 = findHist(name);
220 if (hh1 == NULL) {
222 }
223
224 if (hh1 != NULL) {
225// cout << "do da fit " << endl;
226// m_fGaus->SetParameter(0, 1000);
227// m_fGaus->SetParameter(1, 10);
228// m_fGaus->SetParameter(2, 10);
229// hh1->Fit(m_fGaus, "0");
230// m_hCommon[i]->Fill(j, k, m_fGaus->GetParameter(1));
231// cout << m_fGaus->GetParameter(0) << " " << m_fGaus->GetParameter(1) << " " << m_fGaus->GetParameter(2) << endl;
232 m_hCommon[i]->Fill(j, k, hh1->GetMean());
233 m_hCommonAll->Fill(i, hh1->GetMean());
234 } else {
235 B2INFO("Histo " << name << " not found");
236 }
237
238 name = "hrawPxdCount" + s2;
239 hh1 = findHist(name);
240 if (hh1 == NULL) {
242 }
243
244 if (hh1 != NULL) {
245// cout << "do da fit " << endl;
246// m_fGaus->SetParameter(0, 1000);
247// m_fGaus->SetParameter(1, 100);
248// m_fGaus->SetParameter(2, 10);
249// hh1->Fit(m_fGaus, "0");
250// m_hCounts[i]->Fill(j, k, m_fGaus->GetParameter(1));
251// cout << m_fGaus->GetParameter(0) << " " << m_fGaus->GetParameter(1) << " " << m_fGaus->GetParameter(2) << endl;
252 m_hCounts[i]->Fill(j, k, hh1->GetMean());
253 m_hCountsAll->Fill(i, hh1->GetMean());
254 m_hOccupancyAll->Fill(i, hh1->GetMean() / (250 * 768 / 24)); // Occupancy in percent
255 } else {
256 B2INFO("Histo " << name << " not found");
257 }
258 }
259 }
260 if (m_cSignal[i]) {
261 m_cSignal[i]->cd();
262 if (m_hSignal[i]) m_hSignal[i]->Draw("colz");
263 m_cSignal[i]->Modified();
264 m_cSignal[i]->Update();
265 }
266 if (m_cCommon[i]) {
267 m_cCommon[i]->cd();
268 if (m_hCommon[i]) m_hCommon[i]->Draw("colz");
269 m_cCommon[i]->Modified();
270 m_cCommon[i]->Update();
271 }
272 if (m_cCounts[i]) {
273 m_cCounts[i]->cd();
274 if (m_hCounts[i]) m_hCounts[i]->Draw("colz");
275 m_cCounts[i]->Modified();
276 m_cCounts[i]->Update();
277 }
278 }
279 if (m_cSignalAll) {
280 m_cSignalAll->cd();
281 if (m_hSignalAll) {
282 m_hSignalAll->Scale(1.0 / 24.0); // need to scale
283 m_hSignalAll->Draw();
284 }
285 m_cSignalAll->Modified();
286 m_cSignalAll->Update();
287 }
288 if (m_cCommonAll) {
289 m_cCommonAll->cd();
290 if (m_hCommonAll) {
291 m_hCommonAll->Scale(1.0 / 24.0); // need to scale
292 m_hCommonAll->Draw();
293 }
294 m_cCommonAll->Modified();
295 m_cCommonAll->Update();
296 }
297 if (m_cCountsAll) {
298 m_cCountsAll->cd();
299 if (m_hCountsAll) {
300 // m_hCountsAll->Scale(1.0/24.0); // dont scale counts, we want to have sum
301 // but we would need to scale if we directly calculate occupancy!
302 m_hCountsAll->Draw();
303 }
304 m_cCountsAll->Modified();
305 m_cCountsAll->Update();
306 }
307 if (m_cOccupancyAll) {
308 m_cOccupancyAll->cd();
309 if (m_hOccupancyAll) {
310 m_hOccupancyAll->Scale(1.0 / 24.0); // need to scale
311 m_hOccupancyAll->Draw();
312 }
313 m_cOccupancyAll->Modified();
314 m_cOccupancyAll->Update();
315 }
316}
317
319{
320 B2DEBUG(1, "DQMHistAnalysisPXDFits : endRun called");
321}
322
323
325{
326 B2DEBUG(1, "DQMHistAnalysisPXDFits: terminate called");
327}
328
The base class for the histogram analysis module.
static TH1 * findHist(const std::string &histname, bool onlyIfUpdated=false)
Get histogram from list (no other search).
TH2F * m_hSignal[NUM_MODULES]
2D Signal Histograms
TH2F * m_hCommon[NUM_MODULES]
2D Common Histograms
void initialize() override final
Initializer.
TCanvas * m_cSignal[NUM_MODULES]
2D Signal Canvases
TCanvas * m_cOccupancyAll
All Occupancy Canvas.
TCanvas * m_cCommonAll
All Common Canvas.
TCanvas * m_cCommon[NUM_MODULES]
2D Common Canvases
TH1F * m_hOccupancyAll
All Occupancy Histogram.
std::map< int, int > m_inx_to_id
maps from index to VXDid
TCanvas * m_cSignalAll
All Signal Canvas.
void terminate() override final
This method is called at the end of the event processing.
void event() override final
This method is called for each event.
TH1F * m_hSignalAll
All Signal Histogram.
std::string m_histogramDirectoryName
Histogram doirectory.
void endRun() override final
This method is called if the current run ends.
TH1F * m_hCountsAll
All Counts Histogram.
std::map< int, int > m_id_to_inx
maps from VXDid to index
void beginRun() override final
Called when entering a new run.
TH1F * m_hCommonAll
All Common Histogram.
TH2F * m_hCounts[NUM_MODULES]
2D Counts Histograms
TCanvas * m_cCounts[NUM_MODULES]
2D Counts Canvases
TCanvas * m_cCountsAll
All Counts Canvas.
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
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.
STL namespace.