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