Belle II Software  release-08-01-10
ARICHChannelHist.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 #include <arich/utility/ARICHChannelHist.h>
10 #include <framework/logging/Logger.h>
11 #include <math.h>
12 #include <iostream>
13 #include <algorithm>
14 #include <TGraph.h>
15 #include <TVector2.h>
16 #include <TList.h>
17 
18 using namespace std;
19 using namespace Belle2;
20 
21 ARICHChannelHist::ARICHChannelHist(const char* name, const char* title, int type,
22  const std::vector<unsigned>& moduleIDs) : TH2Poly()
23 {
24 
25  m_type = type;
26  SetName(name);
27  SetTitle(title);
28  m_hapd2binMap.assign(420, 0);
29 
30  // positions of HAPDs and channel mapping (avoid using DB classes...)
31  const double rs[7] = {57.35, 65.81, 74.37, 82.868, 91.305, 99.794, 108.185};
32  const unsigned nhapds[7] = {42, 48, 54, 60, 66, 72, 78};
33  const unsigned chmap[144] = {88, 86, 96, 87, 75, 72, 97, 108, 73, 74, 98, 109, 84, 85, 120, 110, 76, 77, 132, 121, 136, 124, 99, 133, 125, 113, 122, 111, 101, 137, 134, 123, 89, 112, 100, 135, 52, 28, 3, 40, 41, 5, 15, 2, 17, 29, 27, 14, 4, 16, 1, 26, 53, 65, 0, 13, 48, 49, 39, 12, 61, 62, 25, 38, 63, 60, 24, 37, 64, 50, 51, 36, 91, 115, 140, 103, 102, 138, 128, 141, 126, 114, 116, 129, 139, 127, 142, 117, 90, 78, 143, 130, 95, 94, 104, 131, 82, 81, 118, 105, 80, 83, 119, 106, 79, 93, 92, 107, 55, 57, 47, 56, 68, 71, 46, 35, 70, 69, 45, 34, 59, 58, 23, 33, 67, 66, 11, 22, 7, 19, 44, 10, 18, 30, 21, 32, 42, 6, 9, 20, 54, 31, 43, 8};
34  const double chns[12] = { -2.88, -2.37, -1.86, -1.35, -0.84, -0.33, 0.33, 0.84, 1.35, 1.86, 2.37, 2.88};
35 
36  float size = 0.5 / 2. - 0.01;
37  if (m_type == 1) size = 7.0 / 2. - 0.5;
38  if (m_type == 2) size = 3.3 / 2.;
39 
40  float X[5], Y[5], globX[5], globY[5];
41  X[0] = -size; Y[0] = -size;
42  X[1] = size; Y[1] = -size;
43  X[2] = size; Y[2] = size;
44  X[3] = -size; Y[3] = size;
45  X[4] = -size; Y[4] = -size;
46 
47  int nhapd = 1;
48  unsigned iring = 0;
49  unsigned ihapd = 0;
50 
51  std::vector<unsigned> ids;
52  if (moduleIDs.size() > 0) ids = moduleIDs;
53  else {
54  for (int hapdID = 1; hapdID < 421; hapdID++) {
55  ids.push_back(hapdID);
56  }
57  }
58 
59  // HAPD bins
60  if (m_type == 1) {
61  for (int hapdID = 1; hapdID < 421; hapdID++) {
62  //for (unsigned hapdID : ids) {
63  //m_hapd2binMap[hapdID - 1] = nhapd;
64  //nhapd++;
65  float r = rs[iring];
66  float dphi = 2.*M_PI / nhapds[iring];
67  float fi = dphi / 2. + ihapd * dphi;
68  TVector2 centerPos(r * cos(fi), r * sin(fi));
69  for (int i = 0; i < 5; i++) {
70  float rotX = X[i] * cos(fi) - Y[i] * sin(fi);
71  float rotY = X[i] * sin(fi) + Y[i] * cos(fi);
72  globX[i] = rotX + centerPos.X();
73  globY[i] = rotY + centerPos.Y();
74  }
75  if (std::find(ids.begin(), ids.end(), hapdID) != ids.end()) {
76  m_hapd2binMap[hapdID - 1] = nhapd;
77  nhapd++;
78  TGraph* mybox = new TGraph(5, globX, globY);
79  mybox->SetName((to_string(hapdID)).c_str());
80  AddBin(mybox);
81  }
82  ihapd++;
83  if (ihapd == nhapds[iring]) { iring++; ihapd = 0;}
84  }
85 
86  } else if (m_type == 0) {
87  for (int hapdID = 1; hapdID < 421; hapdID++) {
88  //for (unsigned hapdID : ids) {
89  // m_hapd2binMap[hapdID - 1] = nhapd;
90  // nhapd++;
91  float dphi = 2.*M_PI / nhapds[iring];
92  float fi = dphi / 2. + ihapd * dphi;
93  float r = rs[iring];
94  for (int chID = 0; chID < 144; chID++) {
95 
96  unsigned chX = chmap[chID] % 12;
97  unsigned chY = chmap[chID] / 12;
98  TVector2 hapdPos(r * cos(fi), r * sin(fi));
99  TVector2 locPos(chns[chX], chns[chY]);
100  TVector2 centerPos = hapdPos + locPos.Rotate(fi);
101 
102  for (int i = 0; i < 5; i++) {
103  float rotX = X[i] * cos(fi) - Y[i] * sin(fi);
104  float rotY = X[i] * sin(fi) + Y[i] * cos(fi);
105  globX[i] = rotX + centerPos.X();
106  globY[i] = rotY + centerPos.Y();
107  }
108  if (std::find(ids.begin(), ids.end(), hapdID) != ids.end()) {
109  m_hapd2binMap[hapdID - 1] = nhapd;
110  if (chID == 143) nhapd++;
111  TGraph* mybox = new TGraph(5, globX, globY);
112  mybox->SetName((to_string(hapdID)).c_str());
113  AddBin(mybox);
114  }
115  }
116  ihapd++;
117  if (ihapd == nhapds[iring]) { iring++; ihapd = 0;}
118  }
119  } else if (m_type == 2) {
120  size += 0.2;
121  for (int hapdID = 1; hapdID < 421; hapdID++) {
122  float dphi = 2.*M_PI / nhapds[iring];
123  float fi = dphi / 2. + ihapd * dphi;
124  float r = rs[iring];
125  TVector2 hapdPos(r * cos(fi), r * sin(fi));
126  for (int chipID = 0; chipID < 4; chipID++) {
127  TVector2 locPos(-size + (chipID / 2)*size * 2, size - (chipID % 2)*size * 2);
128  TVector2 centerPos = hapdPos + locPos.Rotate(fi);
129 
130  for (int i = 0; i < 5; i++) {
131  float rotX = X[i] * cos(fi) - Y[i] * sin(fi);
132  float rotY = X[i] * sin(fi) + Y[i] * cos(fi);
133  globX[i] = rotX + centerPos.X();
134  globY[i] = rotY + centerPos.Y();
135  }
136 
137  if (std::find(ids.begin(), ids.end(), hapdID) != ids.end()) {
138  m_hapd2binMap[hapdID - 1] = nhapd;
139  if (chipID == 3) nhapd++;
140  TGraph* mybox = new TGraph(5, globX, globY);
141  mybox->SetName((to_string(hapdID)).c_str());
142  AddBin(mybox);
143  }
144  }
145  ihapd++;
146  if (ihapd == nhapds[iring]) { iring++; ihapd = 0;}
147  }
148  } else std::cout << "Invalid histogram type! use 0 for channel bins or 1 for HAPD bins" << std::endl;
149  SetOption("colz");
150  SetStats(0);
151  GetXaxis()->SetLimits(-115., 115.);
152  GetYaxis()->SetLimits(-115., 115.);
153 }
154 
155 void ARICHChannelHist::Draw(Option_t* option)
156 {
157  TH2Poly::Draw(option);
158  double rlin = 40;
159  double rlout = 113;
160  for (int isec = 0; isec < 6; isec++) {
161  double x1 = rlin * cos(M_PI / 3.*isec);
162  double x2 = rlout * cos(M_PI / 3.*isec);
163  double y1 = rlin * sin(M_PI / 3.*isec);
164  double y2 = rlout * sin(M_PI / 3.*isec);
165  lines[isec] = TLine(x1, y1, x2, y2);
166  lines[isec].Draw();
167  x1 = rlin * cos(M_PI / 3.*isec + M_PI / 6.);
168  y1 = rlin * sin(M_PI / 3.*isec + M_PI / 6.);
169  labels[isec] = TText(x1, y1, TString::Format("S-%d", isec + 1));
170  labels[isec].SetTextAlign(22);
171  labels[isec].SetTextSize(0.03);
172  labels[isec].Draw();
173  }
174 }
175 
176 
177 void ARICHChannelHist::fillBin(unsigned hapdID, unsigned chID, double weight)
178 {
179  unsigned chIndex = 0;
180  if (m_type == 0) chIndex = (m_hapd2binMap[hapdID - 1] - 1) * 144 + chID + 1;
181  if (m_type == 2) chIndex = (m_hapd2binMap[hapdID - 1] - 1) * 4 + chID + 1;
182  SetBinContent(chIndex, GetBinContent(chIndex) + weight);
183 }
184 
185 void ARICHChannelHist::setBinContent(unsigned hapdID, unsigned chID, double value)
186 {
187  unsigned chIndex = 0;
188  if (m_type == 0) chIndex = (m_hapd2binMap[hapdID - 1] - 1) * 144 + chID + 1;
189  if (m_type == 2) chIndex = (m_hapd2binMap[hapdID - 1] - 1) * 4 + chID + 1;
190  SetBinContent(chIndex, value);
191 }
192 
193 void ARICHChannelHist::setBinContent(unsigned hapdID, double value)
194 {
195  SetBinContent(m_hapd2binMap[hapdID - 1], value);
196 }
197 
198 void ARICHChannelHist::fillBin(unsigned hapdID, double weight)
199 {
200  SetBinContent(m_hapd2binMap[hapdID - 1], GetBinContent(hapdID) + weight);
201 }
202 
204 {
205 
206  int nbins = hist->GetNbinsX();
207  if (m_type == 1) {
208  if (nbins < 420) { B2ERROR("Number of bins in histogram small than number of ChannelHist bins!"); return;}
209  if (nbins == 420) for (int i = 0; i < nbins; i++) setBinContent(i + 1, hist->GetBinContent(i + 1));
210  if (nbins == 420 * 4) {
211  for (int i = 0; i < 420; i++) {
212  for (int j = 0; j < 4; j++) fillBin(i + 1, hist->GetBinContent(i * 4 + j + 1));
213  }
214  }
215  if (nbins == 420 * 144) {
216  for (int i = 0; i < 420; i++) {
217  for (int j = 0; j < 144; j++) fillBin(i + 1, hist->GetBinContent(i * 144 + j + 1));
218  }
219  }
220  } else if (m_type == 0) {
221  if (nbins < 420 * 144) { B2ERROR("Number of bins in histogram small than number of ChannelHist bins!"); return;}
222  for (int i = 0; i < 420; i++) {
223  for (int j = 0; j < 144; j++) setBinContent(i + 1, j, hist->GetBinContent(i * 144 + j + 1));
224  }
225  } else if (m_type == 2) {
226  if (nbins < 420 * 4) { B2ERROR("Number of bins in histogram small than number of ChannelHist bins!"); return;}
227  if (nbins == 420 * 4) for (int i = 0; i < nbins; i++) setBinContent(i / 4 + 1, i % 4, hist->GetBinContent(i + 1));
228  if (nbins == 420 * 144) {
229  for (int i = 0; i < 420; i++) {
230  for (int j = 0; j < 144; j++) fillBin(i + 1, j / 36, hist->GetBinContent(i * 144 + j + 1));
231  }
232  }
233  } else return;
234 }
235 
236 void ARICHChannelHist::setPoly(TH2Poly* poly)
237 {
238 
239  if (poly->GetNumberOfBins() == 0) {
240  for (const auto&& bin : *fBins) {
241  poly->AddBin((TGraph*)((TH2PolyBin*)bin)->GetPolygon());
242  }
243  }
244  if (poly->GetNumberOfBins() != GetNumberOfBins()) {std::cout << "Mismatch between number of bins in TH2Poly and ARICHChannelHist" << std::endl; return;}
245 
246  double max = poly->GetMaximum();
247  for (int i = 1; i < GetNumberOfBins() + 1; i++) {
248  poly->SetBinContent(i, GetBinContent(i) > max ? max : GetBinContent(i));
249  }
250  return;
251 }
252 
void fillFromTH1(TH1 *hist)
Fill the channelHist from the histogram Type 0 channelHist has to be filled with 420*144bin TH1 (each...
void Draw(Option_t *option="") override
Draw the histogram.
void setBinContent(unsigned hapdID, unsigned chID, double value)
Set content of bin corresponding to hapd hapdID and channel chID.
std::vector< unsigned > m_hapd2binMap
map of bins
void fillBin(unsigned hapdID, unsigned chID, double weight=1.)
Add entry to bin corresponding to hapd hapdID and channel chID.
int m_type
histogram type
void setPoly(TH2Poly *poly)
Fill pure TH2Poly from ARICHChannelHist, makes bins and fills content.
Abstract base class for different kinds of events.