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