Belle II Software development
DQMHistAnalysisPXDTrackCharge.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 : DQMHistAnalysisPXDTrackCharge.cc
10// Description : Analysis of PXD Cluster Charge
11//-
12
13#include <dqm/analysis/modules/DQMHistAnalysisPXDTrackCharge.h>
14#include <TROOT.h>
15#include <TStyle.h>
16#include <TLatex.h>
17#include <vxd/geometry/GeoCache.h>
18
19#include <RooDataHist.h>
20#include <RooAbsPdf.h>
21#include <RooPlot.h>
22#include <RooFitResult.h>
23
24
25using namespace std;
26using namespace Belle2;
27
28//-----------------------------------------------------------------
29// Register the Module
30//-----------------------------------------------------------------
31REG_MODULE(DQMHistAnalysisPXDTrackCharge);
32
33//-----------------------------------------------------------------
34// Implementation
35//-----------------------------------------------------------------
36
39{
40 // This module CAN NOT be run in parallel!
41 setDescription("DQM Analysis for PXD Track-Cluster Charge");
42
43 // Parameter definition
44 addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of Histogram dir", std::string("PXDER"));
45 addParam("RangeLow", m_rangeLow, "Lower border for fit", 20.);
46 addParam("RangeHigh", m_rangeHigh, "High border for fit", 80.);
47// addParam("PeakBefore", m_peakBefore, "Range for fit before peak (positive)", 5.);
48// addParam("PeakAfter", m_peakAfter, "Range for after peak", 40.);
49 addParam("ColorAlert", m_color, "Whether to show the color alert", true);
50 addParam("excluded", m_excluded, "excluded module (indizes starting from 0 to 39)");
51 B2DEBUG(99, "DQMHistAnalysisPXDTrackCharge: Constructor done.");
52}
53
55{
56}
57
59{
60 B2DEBUG(99, "DQMHistAnalysisPXDTrackCharge: initialized.");
61
64
65 m_rfws = new RooWorkspace("w");
66 m_rfws->factory("Landau::landau(x[0,100],ml[20,10,50],sl[5,1,30])");
67 m_rfws->factory("Gaussian::gauss(x,mg[0],sg[2,0.1,10])");
68 m_rfws->factory("FCONV::lxg(x,landau,gauss)");
69
70 m_x = m_rfws->var("x");
71 m_x->setRange("signal", m_rangeLow, m_rangeHigh);
72
73 // collect the list of all PXD Modules in the geometry here
74 std::vector<VxdID> sensors = geo.getListOfSensors();
75 for (VxdID& aVxdID : sensors) {
76 VXD::SensorInfoBase info = geo.getSensorInfo(aVxdID);
77 if (info.getType() != VXD::SensorInfoBase::PXD) continue;
78 m_PXDModules.push_back(aVxdID); // reorder, sort would be better
79
80 std::string name = "PXD_Track_Cluster_Charge_" + (std::string)aVxdID;
81 std::replace(name.begin(), name.end(), '.', '_');
82 m_cChargeMod[aVxdID] = new TCanvas((m_histogramDirectoryName + "/c_Fit_" + name).data());
83 if (aVxdID == VxdID("1.5.1")) {
84 for (int s = 0; s < 6; s++) {
85 for (int d = 0; d < 4; d++) {
86 m_cChargeModASIC[aVxdID][s][d] = new TCanvas((m_histogramDirectoryName + "/c_Fit_" + name + Form("_s%d_d%d", s + 1, d + 1)).data());
87 }
88 }
89 m_hChargeModASIC2d[aVxdID] = new TH2F(("hPXD_TCChargeMPV_" + name).data(),
90 ("PXD TCCharge MPV " + name + ";Switcher;DCD;MPV").data(),
91 6, 0.5, 6.5, 4, 0.5, 4.5);
92 m_cChargeModASIC2d[aVxdID] = new TCanvas((m_histogramDirectoryName + "/c_TCCharge_MPV_" + name).data());
93 }
94 }
95 if (m_PXDModules.size() == 0) {
96 // Backup if no geometry is present (testing...)
97 B2WARNING("No PXDModules in Geometry found! Use hard-coded setup.");
98 std::vector <string> mod = {
99 "1.1.1", "1.1.2", "1.2.1", "1.2.2", "1.3.1", "1.3.2", "1.4.1", "1.4.2",
100 "1.5.1", "1.5.2", "1.6.1", "1.6.2", "1.7.1", "1.7.2", "1.8.1", "1.8.2",
101 "2.1.1", "2.1.2", "2.2.1", "2.2.2", "2.3.1", "2.3.2", "2.4.1", "2.4.2",
102 "2.5.1", "2.5.2", "2.6.1", "2.6.2", "2.7.1", "2.7.2", "2.8.1", "2.8.2",
103 "2.9.1", "2.9.2", "2.10.1", "2.10.2", "2.11.1", "2.11.2", "2.12.1", "2.12.2"
104 };
105 for (auto& it : mod) m_PXDModules.push_back(VxdID(it));
106 }
107 std::sort(m_PXDModules.begin(), m_PXDModules.end()); // back to natural order
108
109 gROOT->cd(); // this seems to be important, or strange things happen
110
111 m_cTrackedClusters = new TCanvas((m_histogramDirectoryName + "/c_TrackedClusters").data());
112 m_hTrackedClusters = new TH1F("hPXDTrackedClusters", "PXD Tracked Clusters/Event;Module", 40, 0, 40);
113 m_hTrackedClusters->Draw();
114 auto ax = m_hTrackedClusters->GetXaxis();
115 if (ax) {
116 ax->Set(m_PXDModules.size(), 0, m_PXDModules.size());
117 for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
118 TString ModuleName = (std::string)m_PXDModules[i];
119 ax->SetBinLabel(i + 1, ModuleName);
120 }
121 } else B2ERROR("no axis");
122
123 m_cCharge = new TCanvas((m_histogramDirectoryName + "/c_TrackCharge").data());
125
126 m_gCharge = new TGraphErrors();
127 m_gCharge->SetName("Track_Cluster_Charge");
128 m_gCharge->SetTitle("Track Cluster Charge");
129
131 m_line_up = new TLine(0, 10, m_PXDModules.size(), 10);
132 m_line_mean = new TLine(0, 16, m_PXDModules.size(), 16);
133 m_line_low = new TLine(0, 3, m_PXDModules.size(), 3);
134 m_line_up->SetHorizontal(true);
135 m_line_up->SetLineColor(kMagenta);// Green
136 m_line_up->SetLineWidth(3);
137 m_line_up->SetLineStyle(7);
138 m_line_mean->SetHorizontal(true);
139 m_line_mean->SetLineColor(kGreen);// Black
140 m_line_mean->SetLineWidth(3);
141 m_line_mean->SetLineStyle(4);
142 m_line_low->SetHorizontal(true);
143 m_line_low->SetLineColor(kMagenta);
144 m_line_low->SetLineWidth(3);
145 m_line_low->SetLineStyle(7);
146
147 m_fMean = new TF1("f_Mean", "pol0", 0, m_PXDModules.size());
148 m_fMean->SetParameter(0, 50);
149 m_fMean->SetLineColor(kYellow);
150 m_fMean->SetLineWidth(3);
151 m_fMean->SetLineStyle(7);
152 m_fMean->SetNpx(m_PXDModules.size());
153 m_fMean->SetNumberFitPoints(m_PXDModules.size());
154
155 registerEpicsPV("PXD:TrackCharge:Mean", "Mean");
156 registerEpicsPV("PXD:TrackCharge:Diff", "Diff");
157 registerEpicsPV("PXD:TrackCharge:Status", "Status");
158}
159
160
162{
163 B2DEBUG(99, "DQMHistAnalysisPXDTrackCharge: beginRun called.");
164
165 m_cCharge->Clear();
166}
167
169{
170 if (!m_cCharge) return;
171
172 gStyle->SetOptStat(0);
173 gStyle->SetStatStyle(1);
174 gStyle->SetOptDate(22);// Date and Time in Bottom Right, does no work
175
176 bool enough = false;
177
178// auto landau = m_rfws->pdf("landau");
179// auto gauss = m_rfws->pdf("gauss");
180 auto model = m_rfws->pdf("lxg");
181
182 auto ml = m_rfws->var("ml");
183// auto sl = m_rfws->var("sl");
184// auto mg = m_rfws->var("mg");
185// auto sg = m_rfws->var("sg");
186
187 {
188 std::string name = "Tracked_Clusters"; // new name
189 TH1* hh2 = findHist(m_histogramDirectoryName, "PXD_Tracked_Clusters", true);
190 if (hh2) {// update only if histogram is updated
191 m_cTrackedClusters->Clear();
192 m_cTrackedClusters->cd();
193 m_hTrackedClusters->Reset();
194
195 auto scale = hh2->GetBinContent(0);// overflow misused as event counter!
196 if (scale > 0) {
197 int j = 1;
198 for (int i = 0; i < 64; i++) {
199 auto layer = (((i >> 5) & 0x1) + 1);
200 auto ladder = ((i >> 1) & 0xF);
201 auto sensor = ((i & 0x1) + 1);
202
203 auto id = Belle2::VxdID(layer, ladder, sensor);
204 // Check if sensor exist
205 if (Belle2::VXD::GeoCache::getInstance().validSensorID(id)) {
206 m_hTrackedClusters->SetBinContent(j, hh2->GetBinContent(i + 1) / scale);
207 j++;
208 }
209 }
210 }
211 m_hTrackedClusters->SetName(name.data());
212 m_hTrackedClusters->SetTitle("Tracked Clusters/Event");
213 m_hTrackedClusters->SetFillColor(kWhite);
214 m_hTrackedClusters->SetStats(kFALSE);
215 m_hTrackedClusters->SetLineStyle(1);// 2 or 3
216 m_hTrackedClusters->SetLineColor(kBlue);
217 m_hTrackedClusters->Draw("hist");
218
219 // get ref histogram, assumes no m_histogramDirectoryName directory in ref file
220// auto href2 = findHistInFile(m_refFile, name);
221//
222// if (href2) {
223// href2->SetLineStyle(3);// 2 or 3
224// href2->SetLineColor(kBlack);
225// href2->Draw("same,hist");
226// }
227
228 // keep this commented code as we may have excluded modules in phase4
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
235 }
236 }
237
238 m_gCharge->Set(0);
239
240 for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
241 TCanvas* canvas = m_cChargeMod[m_PXDModules[i]];
242 if (canvas == nullptr) continue;
243
244 std::string name = "PXD_Track_Cluster_Charge_" + (std::string)m_PXDModules[i];
245 std::replace(name.begin(), name.end(), '.', '_');
246
247 TH1* hh1 = findHist(m_histogramDirectoryName, name, true);
248 if (hh1) {// update only if histo was updated
249 canvas->cd();
250 canvas->Clear();
251 UpdateCanvas(canvas);
252
253 if (hh1->GetEntries() > 50) {
254
255 auto hdata = new RooDataHist(hh1->GetName(), hh1->GetTitle(), *m_x, (const TH1*) hh1);
256 auto plot = m_x->frame(RooFit::Title(hh1->GetTitle()));
257 /*auto r =*/ model->fitTo(*hdata, RooFit::Range("signal"));
258
259 model->paramOn(plot, RooFit::Format("NELU", RooFit::AutoPrecision(2)), RooFit::Layout(0.6, 0.9, 0.9));
260 hdata->plotOn(plot, RooFit::LineColor(kBlue)/*, RooFit::Range("plot"), RooFit::NormRange("signal")*/);
261 model->plotOn(plot, RooFit::LineColor(kRed), RooFit::Range("signal"), RooFit::NormRange("signal"));
262
263 plot->Draw("");
264
265// model->Print("");
266// ml->Print("");
267// sl->Print("");
268// mg->Print("");
269// sg->Print("");
270// cout << "ZZZ , " << Form("%d%02d%d ,", std::get<0>(t), std::get<1>(t), std::get<2>(t)) << ml->getValV() << "," << ml->getError() << "," << sl->getValV() << "," << sl->getError() << "," << sg->getValV() << "," << sg->getError() << endl;
271
272
273 int p = m_gCharge->GetN();
274 m_gCharge->SetPoint(p, i + 0.49, ml->getValV());
275 m_gCharge->SetPointError(p, 0.1, ml->getError()); // error in x is useless
276 m_monObj->setVariable(("trackcharge_" + (std::string)m_PXDModules[i]).c_str(), ml->getValV(), ml->getError());
277 } else {
278 hh1->Draw("hist"); // avoid to confuse people by showing nothing for low stat
279 }
280
281 // get ref histogram, assumes no m_histogramDirectoryName directory in ref file
282// auto hist2 = findHistInFile(m_refFile, name);
283// if (hist2) {
284// B2DEBUG(20, "Draw Normalized " << hist2->GetName());
285// hist2->SetLineStyle(3);// 2 or 3
286// hist2->SetLineColor(kBlack);
287//
288// canvas->cd();
289//
290// // if draw normalized
291// TH1* h = (TH1*)hist2->Clone(); // Annoying ... Maybe an memory leak? TODO
292// // would it work to scale it each time again?
293// if (abs(hist2->GetEntries()) > 0) h->Scale(hh1->GetEntries() / hist2->GetEntries());
294//
295// h->SetStats(kFALSE);
296// h->Draw("same,hist");
297// }
298
299 // add coloring, cuts? based on fit, compare with ref?
300 canvas->Pad()->SetFrameFillColor(10);
301 if (m_color) {
302 if (hh1->GetEntries() < 100) {
303 // not enough Entries
304 canvas->Pad()->SetFillColor(kGray);
305 } else {
306 canvas->Pad()->SetFillColor(kGreen);
307 }
308 } else {
309 canvas->Pad()->SetFillColor(kWhite);// White
310 }
311
312 canvas->Modified();
313 canvas->Update();
314 UpdateCanvas(canvas);
315
316 // means if ANY plot is > 100 entries, all plots are assumed to be o.k.
317 if (hh1->GetEntries() >= 1000) enough = true;
318 }
319 }
320
321 // now loop per module over asics pairs (1.5.1)
322 for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
323// TCanvas* canvas = m_cChargeMod[m_PXDModules[i]];
324 VxdID& aVxdID = m_PXDModules[i];
325
326 if (m_hChargeModASIC2d[aVxdID]) m_hChargeModASIC2d[aVxdID]->Reset();
327 if (m_cChargeModASIC2d[aVxdID]) m_cChargeModASIC2d[aVxdID]->Clear();
328
329 for (int s = 1; s <= 6; s++) {
330 for (int d = 1; d <= 4; d++) {
331 std::string name = "PXD_Track_Cluster_Charge_" + (std::string)m_PXDModules[i] + Form("_sw%d_dcd%d", s, d);
332 std::replace(name.begin(), name.end(), '.', '_');
333
334 if (m_cChargeModASIC[aVxdID][s - 1][d - 1]) {
335 m_cChargeModASIC[aVxdID][s - 1][d - 1]->Clear();
336 m_cChargeModASIC[aVxdID][s - 1][d - 1]->cd();
337 }
338
339 TH1* hh1 = findHist(m_histogramDirectoryName, name);
340 if (hh1) {
341 double mpv = 0.0;
342 if (hh1->GetEntries() > 50) {
343 auto hdata = new RooDataHist(hh1->GetName(), hh1->GetTitle(), *m_x, (const TH1*) hh1);
344 auto plot = m_x->frame(RooFit::Title(hh1->GetTitle()));
345 /*auto r =*/ model->fitTo(*hdata, RooFit::Range("signal"));
346
347 if (m_cChargeModASIC[aVxdID][s - 1][d - 1]) {
348 model->paramOn(plot, RooFit::Format("NELU", RooFit::AutoPrecision(2)), RooFit::Layout(0.6, 0.9, 0.9));
349 hdata->plotOn(plot, RooFit::LineColor(kBlue)/*, RooFit::Range("plot"), RooFit::NormRange("signal")*/);
350 model->plotOn(plot, RooFit::LineColor(kRed), RooFit::Range("signal"), RooFit::NormRange("signal"));
351 }
352 plot->Draw("");
353
354 mpv = ml->getValV();
355 }
356
357 if (m_hChargeModASIC2d[aVxdID]) {
358 if (mpv > 0.0) m_hChargeModASIC2d[aVxdID]->Fill(s, d, mpv); // TODO check what is s, d
359 }
360 }
361 }
362 }
363
364 // Overview map of ASCI combinations
365 if (m_hChargeModASIC2d[aVxdID] && m_cChargeModASIC2d[aVxdID]) {
366 m_cChargeModASIC2d[aVxdID]->cd();
367 m_hChargeModASIC2d[aVxdID]->Draw("colz");
369 }
370 }
371
372 m_cCharge->cd();
373 m_cCharge->Clear();
374 m_gCharge->SetMinimum(0);
375 m_gCharge->SetMaximum(70);
376 auto ax = m_gCharge->GetXaxis();
377 if (ax) {
378 ax->Set(m_PXDModules.size(), 0, m_PXDModules.size());
379 for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
380 TString ModuleName = (std::string)m_PXDModules[i];
381 ax->SetBinLabel(i + 1, ModuleName);
382 }
383 } else B2ERROR("no axis");
384
385 m_gCharge->SetLineColor(4);
386 m_gCharge->SetLineWidth(2);
387 m_gCharge->SetMarkerStyle(8);
388 m_gCharge->Draw("AP");
389
390 for (auto& it : m_excluded) {
391 static std::map <int, TLatex*> ltmap;
392 auto tt = ltmap[it];
393 if (!tt) {
394 tt = new TLatex(it + 0.5, 0, (" " + std::string(m_PXDModules[it]) + " Module is excluded, please ignore").c_str());
395 tt->SetTextSize(0.035);
396 tt->SetTextAngle(90);// Rotated
397 tt->SetTextAlign(12);// Centered
398 ltmap[it] = tt;
399 }
400 tt->Draw();
401 }
402 m_cCharge->cd(0);
403 m_cCharge->Modified();
404 m_cCharge->Update();
406
407 double data = 0;
408 double diff = 0;
409 if (m_gCharge && enough) {
410// double currentMin, currentMax;
411 m_gCharge->Fit(m_fMean, "R");
412 double mean = m_gCharge->GetMean(2);
413 double maxi = mean + 15;
414 double mini = mean - 15;
415 m_line_up->SetY1(maxi);
416 m_line_up->SetY2(maxi);
417 m_line_mean->SetY1(mean);
418 m_line_mean->SetY2(mean);
419 m_line_low->SetY1(mini);
420 m_line_low->SetY2(mini);
421 data = mean; // m_fMean->GetParameter(0); // we are more interested in the maximum deviation from mean
422 // m_gCharge->GetMinimumAndMaximum(currentMin, currentMax);
423 diff = m_gCharge->GetRMS(2);// RMS of Y
424 // better, max deviation as fabs(data - currentMin) > fabs(currentMax - data) ? fabs(data - currentMin) : fabs(currentMax - data);
425 m_line_up->Draw();
426 m_line_mean->Draw();
427 m_line_low->Draw();
428
429 m_monObj->setVariable("trackcharge", mean, diff);
430 }
431
432 setEpicsPV("Mean", data);
433 setEpicsPV("Diff", diff);
434
435 int status = 0;
436
437 if (!enough) {
438 // not enough Entries
439 m_cCharge->Pad()->SetFillColor(kGray);// Magenta or Gray
440 status = 0; // default
441 } else {
443 if (fabs(data - 30.) > 20. || diff > 12) {
444 m_cCharge->Pad()->SetFillColor(kRed);// Red
445 status = 4;
446 } else if (fabs(data - 30) > 15. || diff > 8) {
447 m_cCharge->Pad()->SetFillColor(kYellow);// Yellow
448 status = 3;
449 } else {
450 m_cCharge->Pad()->SetFillColor(kGreen);// Green
451 status = 2;
452 }
453
454 // FIXME overwrite for now
455 // m_cCharge->Pad()->SetFillColor(kGreen);// Green
456
457 }
458
459 setEpicsPV("Status", status);
460}
461
463{
464 B2DEBUG(99, "DQMHistAnalysisPXDTrackCharge : endRun called");
465}
466
467
469{
470 B2DEBUG(99, "DQMHistAnalysisPXDTrackCharge: terminate called");
471
472 if (m_rfws) delete m_rfws;
473
476 if (m_cCharge) delete m_cCharge;
477 if (m_gCharge) delete m_gCharge;
478 if (m_line_up) delete m_line_up;
479 if (m_line_mean) delete m_line_mean;
480 if (m_line_low) delete m_line_low;
481 if (m_fMean) delete m_fMean;
482}
The base class for the histogram analysis module.
static MonitoringObject * getMonitoringObject(const std::string &name)
Get MonitoringObject with given name (new object is created if non-existing)
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.
int registerEpicsPV(std::string pvname, std::string keyname="")
EPICS related Functions.
void UpdateCanvas(std::string name, bool updated=true)
Mark canvas as updated (or not)
void terminate(void) override final
This method is called at the end of the event processing.
bool m_color
Whether to use the color code for warnings and errors.
std::map< VxdID, std::array< std::array< TCanvas *, 4 >, 6 > > m_cChargeModASIC
Final Canvases for Fit and Ref per ASIC.
std::map< VxdID, TH2F * > m_hChargeModASIC2d
Final Canvas Fit and Ref per ASIC.
void initialize(void) override final
Initializer.
void endRun(void) override final
This method is called if the current run ends.
TLine * m_line_up
TLine object for upper limit of track cluster charge.
TGraphErrors * m_gCharge
Graph covering all modules.
std::vector< VxdID > m_PXDModules
IDs of all PXD Modules to iterate over.
std::string m_histogramDirectoryName
name of histogram directory
TLine * m_line_mean
TLine object for mean of track cluster charge.
TLine * m_line_low
TLine object for lower limit of track cluster charge.
std::vector< int > m_excluded
Indizes of excluded PXD Modules.
TH1F * m_hTrackedClusters
Histogram for TrackedClusters.
TCanvas * m_cTrackedClusters
Final Canvas for TrackedClusters.
std::map< VxdID, TCanvas * > m_cChargeMod
Final Canvases for Fit and Ref.
void beginRun(void) override final
Called when entering a new run.
void event(void) override final
This method is called for each event.
std::map< VxdID, TCanvas * > m_cChargeModASIC2d
Final Canvas Fit and Ref per ASIC.
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setVariable(const std::string &var, float val, float upErr=-1., float dwErr=-1)
set value to float variable (new variable is made if not yet existing)
void addCanvas(TCanvas *canv)
Add Canvas to monitoring object.
Class to facilitate easy access to sensor information of the VXD like coordinate transformations or p...
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 reference 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:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:649
Abstract base class for different kinds of events.
STL namespace.