Belle II Software  release-06-01-15
hitMapMaker.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 // Own include
10 #include <arich/modules/arichDQM/hitMapMaker.h>
11 
12 // ARICH
13 #include <arich/dbobjects/ARICHChannelMapping.h>
14 #include <arich/dbobjects/ARICHMergerMapping.h>
15 #include <arich/dataobjects/ARICHHit.h>
16 
17 #include <mdst/dataobjects/Track.h>
18 
19 // Dataobject classes
20 #include <framework/database/DBObjPtr.h>
21 
22 #include <TH1F.h>
23 #include <TH2F.h>
24 #include <TCanvas.h>
25 #include <TStyle.h>
26 #include <TColor.h>
27 #include <TExec.h>
28 
29 #include <vector>
30 #include <fstream>
31 
32 using namespace std;
33 
34 namespace Belle2 {
40  TH2* moduleHitMap(TH1* hitMap, int moduleID)
41  {
42 
43  int m_moduleID = moduleID;
44 
46 
47  TH2* m_moduleHitMap = new TH2D(Form("HAPDHitMapMod%d", moduleID), Form("HAPD hit map module %d;nChX;nChY", moduleID), 12, 0.5, 12.5,
48  12,
49  0.5, 12.5);
50  TH1* m_hitMap = hitMap;
51 
52  for (int i = 0; i < 144; i++) {
53  int hitsNum = m_hitMap->GetBinContent((m_moduleID - 1) * 144 + i);
54  int xChn, yChn;
55  arichChMap->getXYFromAsic(i, xChn, yChn);
56  if (hitsNum != 0) {
57  m_moduleHitMap->Fill(xChn + 1, yChn + 1, hitsNum);
58  }
59  }
60 
61  return m_moduleHitMap;
62  }
63 
64  TH2* moduleDeadMap(TH1* hitMap, int moduleID)
65  {
66  int m_moduleID = moduleID;
67 
69 
70  TH2* m_moduleDeadMap = new TH2D(Form("HAPDDeadMapMod%d", moduleID), Form("HAPD alive/dead module %d;nChX;nChY", moduleID), 2, 0.5,
71  2.5,
72  2, 0.5, 2.5);
73  TH1* m_hitMap = hitMap;
74 
75  int deadCh[2][2] = {};
76 
77  for (int i = 0; i < 144; i++) {
78  int hitsNum = m_hitMap->GetBinContent((m_moduleID - 1) * 144 + i);
79  int xChn, yChn;
80  arichChMap->getXYFromAsic(i, xChn, yChn);
81  if (hitsNum == 0) deadCh[(int)xChn / 6][(int)yChn / 6]++;
82  }
83 
84  for (int j = 0; j < 2; j++) {
85  for (int k = 0; k < 2; k++) {
86  if (deadCh[j][k] > 18) {
87  m_moduleDeadMap->Fill(j + 1, k + 1, 1);
88  } else {
89  m_moduleDeadMap->Fill(j + 1, k + 1, 10);
90  }
91  }
92  }
93 
94  return m_moduleDeadMap;
95  }
96 
97  TH1* mergerClusterHitMap1D(TH1* hitMap, int mergerID)
98  {
99 
100  int m_mergerID = mergerID;
101 
103  DBObjPtr<ARICHMergerMapping> arichMergerMap;
104 
105  TH1* m_hitMap = hitMap;
106 
107  std::vector<int> moduleIDs;
108  for (int i = 1; i < 7; i++) {
109  moduleIDs.push_back(arichMergerMap->getModuleID(m_mergerID, i));
110  }
111 
112  TH1D* m_mergerHitMap1D = new TH1D("MergerHitMap1D", Form("Hit map in Merger Board %d", m_mergerID), 144 * 6, -0.5, 144 * 6 - 0.5);
113  for (int i = 1; i < 7; i++) {
114  for (int j = 0; j < 144; j++) {
115  int hitsNum = m_hitMap->GetBinContent((moduleIDs[i] - 1) * 144 + i);
116  m_mergerHitMap1D->Fill(144 * (i - 1) + j, hitsNum);
117  }
118  }
119  return m_mergerHitMap1D;
120  }
121 
122  TCanvas* mergerClusterHitMap2D(TH1* hitMap, int mergerID)
123  {
124  int m_mergerID = mergerID;
125 
127  DBObjPtr<ARICHMergerMapping> arichMergerMap;
128 
129  TH1* m_hitMap = hitMap;
130 
131 
132  std::vector<int> moduleIDs;
133  for (int i = 1; i < 7; i++) {
134  moduleIDs.push_back(arichMergerMap->getModuleID(m_mergerID, i));
135  }
136 
137  TCanvas* m_mergerHitMap = new TCanvas("MergerHitMap", "Hit map in Merger Board", 600, 400);
138  m_mergerHitMap->Divide(3, 2);
139  for (int i = 1; i < 7; i++) {
140  m_mergerHitMap->cd(i);
141  moduleHitMap(m_hitMap, moduleIDs[i])->Draw("coloz");
142  }
143  return m_mergerHitMap;
144  }
145 
146  TCanvas* sectorHitMap(TH1* hitMap, int sector)
147  {
148  // TH1* m_hitMap = NULL;
149  TH2* m_moduleHitMap = NULL;
150 
151  //m_hitMap = hitMap;
152  TPad* p_hitMaps[421] = {};
153  TCanvas* m_sectorHitMap = new TCanvas(Form("c_hitMap%d", sector - 1), Form("Hit map of sector%d", sector - 1), 600, 400);
154  for (int i = 1; i < 421; i++) {
155  for (int iRing = 0; iRing < 7; iRing++) {
156  if (((iRing + 13)*iRing / 2) * 6 + (iRing + 7) * (sector - 1) < i && i <= ((iRing + 13)*iRing / 2) * 6 + (iRing + 7) * (sector)) {
157  m_sectorHitMap->cd();
158  p_hitMaps[i] = new TPad(Form("p_hitMap%d", i), "",
159  (double)((double)(6 - iRing) / 2 + ((i - ((iRing + 13)*iRing / 2) * 6 + (iRing + 7) * (sector - 1) - 1) % (iRing + 7))) / 13,
160  (double)iRing / 7, (double)((double)(8 - iRing) / 2 + ((i - ((iRing + 13)*iRing / 2) * 6 + (iRing + 7) * (sector - 1) - 1) %
161  (iRing + 7))) / 13, (double)(iRing + 1) / 7);
162  p_hitMaps[i]->Draw();
163  p_hitMaps[i]->SetNumber(i);
164  m_sectorHitMap->cd(i);
165  m_moduleHitMap = moduleHitMap(hitMap, i);
166  m_moduleHitMap->SetTitleSize(0, "xyz");
167  m_moduleHitMap->SetTitle(0);
168  m_moduleHitMap->SetLineWidth(1);
169  gStyle->SetLabelColor(0, "xyz");
170  m_moduleHitMap->SetStats(0);
171  m_moduleHitMap->Draw("col");
172  }
173  }
174  }
175  return m_sectorHitMap;
176  }
177 
178  TCanvas* sectorDeadMap(TH1* hitMap, int sector)
179  {
180  //TH1* m_hitMap = NULL;
181  TH2* m_moduleDeadMap = NULL;
182  TExec* ex1 = NULL;
183 
184  //m_hitMap = hitMap;
185  TPad* p_hitMaps[421] = {};
186  TCanvas* m_sectorDeadMap = new TCanvas(Form("c_deadMap%d", sector - 1), Form("Dead chip map of sector%d", sector - 1), 600, 400);
187  for (int i = 1; i < 421; i++) {
188  for (int iRing = 0; iRing < 7; iRing++) {
189  if (((iRing + 13)*iRing / 2) * 6 + (iRing + 7) * (sector - 1) < i && i <= ((iRing + 13)*iRing / 2) * 6 + (iRing + 7) * (sector)) {
190  m_sectorDeadMap->cd();
191  p_hitMaps[i] = new TPad(Form("p_deadMap%d", i), "",
192  (double)((double)(6 - iRing) / 2 + ((i - ((iRing + 13)*iRing / 2) * 6 + (iRing + 7) * (sector - 1) - 1) % (iRing + 7))) / 13,
193  (double)iRing / 7, (double)((double)(8 - iRing) / 2 + ((i - ((iRing + 13)*iRing / 2) * 6 + (iRing + 7) * (sector - 1) - 1) %
194  (iRing + 7))) / 13, (double)(iRing + 1) / 7);
195  p_hitMaps[i]->Draw();
196  p_hitMaps[i]->SetNumber(i);
197  m_sectorDeadMap->cd(i);
198  m_moduleDeadMap = moduleDeadMap(hitMap, i);
199  m_moduleDeadMap->SetTitleSize(0, "xyz");
200  m_moduleDeadMap->SetTitle(0);
201  m_moduleDeadMap->SetLineWidth(1);
202  m_moduleDeadMap->SetStats(0);
203  gStyle->SetLabelColor(0, "xyz");
204  ex1 = new TExec("ex1", "deadPalette();");
205  ex1->Draw();
206  m_moduleDeadMap->Draw("colz");
207  }
208  }
209  }
210  return m_sectorDeadMap;
211  }
212 
213  void deadPalette()
214  {
215  static Int_t colors[50];
216  static Bool_t initialized = kFALSE;
217  Double_t Red[3] = { 1.00, 0.00, 0.00};
218  Double_t Green[3] = { 0.00, 0.00, 0.00};
219  Double_t Blue[3] = { 0.00, 1.00, 1.00};
220  Double_t Length[3] = { 0.00, 0.20, 1.00 };
221  if (!initialized) {
222  Int_t FI = TColor::CreateGradientColorTable(3, Length, Red, Green, Blue, 50);
223  for (int i = 0; i < 50; i++) colors[i] = FI + i;
224  initialized = kTRUE;
225  return;
226  }
227  gStyle->SetPalette(50, colors);
228  }
229 
231 }
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
TH2 * moduleHitMap(TH1 *hitMap, int moduleID)
Make hit map in HAPD view (12*12 channels)
Definition: hitMapMaker.cc:40
TH2 * moduleDeadMap(TH1 *hitMap, int moduleID)
Make chip dead/alive map in HAPD view (2*2 chips)
Definition: hitMapMaker.cc:64
TCanvas * mergerClusterHitMap2D(TH1 *hitMap, int mergerID)
Make display of 6 HAPDs' 2D hit map of the Merger Board.
Definition: hitMapMaker.cc:122
TH1 * mergerClusterHitMap1D(TH1 *hitMap, int mergerID)
Make 1D hit map of specified Merger Board.
Definition: hitMapMaker.cc:97
TCanvas * sectorDeadMap(TH1 *hitMap, int sector)
Make display of 70 HAPDs' 2D dead/alive map of the sector.
Definition: hitMapMaker.cc:178
TCanvas * sectorHitMap(TH1 *hitMap, int sector)
Make display of 70 HAPDs' 2D hit map of the sector.
Definition: hitMapMaker.cc:146
Abstract base class for different kinds of events.