Belle II Software  release-08-01-10
DQMHistAnalysisV0.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 : DQMHistAnalysisV0.cc
10 // Description : Overlay plotting for V0
11 //-
12 
13 
14 #include <dqm/analysis/modules/DQMHistAnalysisV0.h>
15 
16 #include <TFile.h>
17 #include <TROOT.h>
18 #include <TStyle.h>
19 #include <TString.h>
20 #include <TH2.h>
21 #include <TGraph.h>
22 
23 using namespace std;
24 using namespace Belle2;
25 
26 //-----------------------------------------------------------------
27 // Register the Module
28 //-----------------------------------------------------------------
29 REG_MODULE(DQMHistAnalysisV0);
30 
31 //-----------------------------------------------------------------
32 // Implementation
33 //-----------------------------------------------------------------
34 
35 DQMHistAnalysisV0Module::DQMHistAnalysisV0Module()
37 {
38  //Parameter definition
39  addParam("OverlayPath", m_OverlayPath, "Path to CAD drawings", std::string(""));
40 
41 }
42 
43 
45 {
46  B2INFO("DQMHistAnalysisV0: initialized.");
47 
48  gROOT->cd();
49  for (int i = 0; i < 32; i++) {
50  m_c_xvsy[i] = new TCanvas(Form("V0Object/c_xvsy[%i]", i), Form("c_xvsy[%i]", i), 800, 800);
51  }
52 
53  m_c_xvsz = new TCanvas("V0Object/c_xvsz", "c_xvsz", 1500, 400); //Stretch in x to enhance visibility
54 
55  auto* m_fh = new TFile(Form("%s/v0cad.root", m_OverlayPath.c_str()));
56 
57  contLevelXY.resize(32);
58  for (int i = 0; i < 32; i++) {
59  contLevelXY[i] = new TList();
60  m_fh->cd();
61  if (m_fh->cd(Form("h_%dc", i))) {
62  for (int j = 0; j < 500; j++) {
63  auto curv = gDirectory->Get(Form("Graph_%d", j));
64  if (!curv) break;
65  contLevelXY[i]->AddLast(curv);
66  }
67  }
68  }
69  {
70  contLevelXZ = new TList();
71  m_fh->cd();
72  if (m_fh->cd("h_xzc")) {
73  for (int j = 0; true; j++) {
74  auto curv = gDirectory->Get(Form("Graph_%d", j));
75  if (!curv) break;
76  contLevelXZ->AddLast(curv);
77  }
78  }
79  }
81  gROOT->cd();
82 }
83 
84 
86 {
87 
88  //gStyle requirements
89  gStyle->SetOptStat(0);
90  gStyle->SetPalette(kViridis, 0, 0.7);
91  for (int i = 0; i < 32; i++) {
92  TH2* h = (TH2*) findHist(Form("V0Objects/xvsy[%i]", i));
93  m_c_xvsy[i]->cd();
94  if (h) h->Draw("COLZ");
95 
96  TList* c = contLevelXY[i];
97  if (c && c->GetSize() > 0) {
98  auto* curv = (TGraph*)c->First();
99  for (int j = 0; j < c->GetSize(); j++) {
100  //auto* gc = (TGraph*)curv->Clone();
101  //gc->Draw("C");
102  curv->Draw("L");
103  curv = (TGraph*)c->After(curv); // Get Next graph
104  }
105  }
106 
107  m_c_xvsy[i]->Modified();
108  m_c_xvsy[i]->Update();
109  }
110 
111  TH2* hxz = (TH2*) findHist("V0Objects/xvsz");
112  m_c_xvsz->cd();
113 
114  // create a new pad for every event? -> mem leak!
115 
116  m_c_xvsz->cd();
117  if (hxz) hxz->Draw("COLZ");
118  {
119  TList* c = contLevelXZ;
120  if (c && c->GetSize() > 0) {
121  auto* curv = (TGraph*)c->First();
122  for (int j = 0; j < c->GetSize(); j++) {
123  //auto* gc = (TGraph*)curv->Clone();
124  //gc->Draw("C");
125  curv->Draw("L");
126  curv = (TGraph*)c->After(curv); // Get Next graph
127  }
128  }
129  }
130 
131  m_c_xvsz->Modified();
132  m_c_xvsz->Update();
133 
134 }
135 
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).
void initialize() override final
Module functions to be called from main process.
std::vector< TList * > contLevelXY
Vector of TList objects for contour in XY plane.
TCanvas * m_c_xvsz
TCanvas object for x vs z plot.
TList * contLevelXZ
TList objects for contour in XZ plane.
void event() override final
Module functions to be called from event process.
std::string m_OverlayPath
Path to overlay file.
TCanvas * m_c_xvsy[32]
TCanvas objects for x vs y plots.
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.