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