Belle II Software  release-08-01-10
ARICHAerogelHist.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 //my
10 #include <arich/utility/ARICHAerogelHist.h>
11 
12 //c, c++
13 #include <iostream>
14 #include <string>
15 #include <iomanip>
16 #include <vector>
17 
18 //root
19 #include <TVector2.h>
20 #include <TCanvas.h>
21 #include <TMath.h>
22 #include <TPad.h>
23 #include <TString.h>
24 #include <TFile.h>
25 #include <TAxis.h>
26 
27 using namespace std;
28 using namespace Belle2;
29 
30 ARICHAerogelHist::ARICHAerogelHist(const char* name, const char* title) : TH2Poly(), m_histName(name), m_histTitle(title)
31 {
32  SetName(m_histName.Data());
33  SetTitle(m_histTitle.Data());
36 
37  if (m_verboseLevel > 0) {
38  std::cout << setw(10) << "ring" << setw(10) << "nTiles" << setw(10) << "tileRmin" << setw(10) << "tileRmax" << std::endl;
39  for (unsigned int i = 0; i < m_nTiles.size(); i++)
40  std::cout << setw(10) << i + 1 << setw(10) << m_nTiles[i] << setw(10) << m_tileRmin[i] << setw(10) << m_tileRmax[i] << std::endl;
42  }
43 
44  unsigned int n;
45  double* x;
46  double* y;
47 
48  double xold;
49  double yold;
50  double xnew;
51  double ynew;
52  double phi;
53 
54  //Add "poly bins" to the histogram
55  //Loop over rings with different radiuses
56  unsigned int iR = 0;
57  for (auto& m : m_verticesMap) {
58  //Loop aerogel tiles (bins) within one ring.
59  double phi0 = m_tileDeltaPhiCenter[iR] / 2.0;
60  double deltaPhi = m_tileDeltaPhiCenter[iR] + m_aerogelAriGapDeltaPhiCenter[iR];
61  for (Int_t j = 0; j < m_nTiles[iR]; j++) {
62  phi = phi0 + deltaPhi * j;
63  n = m.second.size();
64  x = new double [n];
65  y = new double [n];
66  //Loop over polygonal points which defines bins.
67  for (unsigned int i = 0; i < n; i++) {
68  xold = m.second[i].X();
69  yold = m.second[i].Y();
70  makeRotation(xold, yold, xnew, ynew, phi);
71  x[i] = xnew;
72  y[i] = ynew;
73  }
74  AddBin(n, x, y);
75  delete []x;
76  delete []y;
77  }
78  iR++;
79  }// for (auto& m: m_verticesMap) {
80 
81 }
82 
83 Int_t ARICHAerogelHist::GetBinIDFromRingColumn(Int_t ring, Int_t column)
84 {
85 
86  Int_t binID = 0;
87  unsigned int ringu = ring;
88  //cout<<" ewewe "<<m_nTiles[ring-1]<<endl;
89  if (ringu > m_nTiles.size() || ringu < 1)
90  return -999;
91  if (column > m_nTiles[ring - 1] || column < 1)
92  return -999;
93  for (Int_t i = 1; i < ring; i++)
94  binID += m_nTiles[i - 1];
95  return binID + column;
96 
97 }
98 
99 void ARICHAerogelHist::makeRotation(double xold, double yold, double& xnew, double& ynew, double phi)
100 {
101 
102  TVector2 v(xold, yold);
103  TVector2 vrot;
104  vrot = v.Rotate(phi);
105  xnew = vrot.X();
106  ynew = vrot.Y();
107 
108 }
109 
110 void ARICHAerogelHist::DrawHisto(TString opt = "ZCOLOT text same", TString outDirName = "./")
111 {
112 
113  TCanvas* c1 = new TCanvas("c1", "c1", 1000, 1000);
114  c1->SetTitle(m_histTitle.Data());
115  c1->SetName(m_histName.Data());
116  c1->SetRightMargin(0.17);
117  c1->SetLeftMargin(0.12);
118  c1->SetTopMargin(0.15);
119  c1->SetBottomMargin(0.15);
120 
121  TH2F* frame = new TH2F("h2", "h2", 40, -1200, 1200, 40, -1200, 1200);
122  //TH2F *frame = new TH2F("h2","h2",40, 400, 650, 40, -50, 200);
123  frame->GetXaxis()->SetTitle("x (mm)");
124  frame->GetYaxis()->SetTitle("y (mm)");
125  frame->GetXaxis()->CenterTitle();
126  frame->GetYaxis()->CenterTitle();
127  frame->GetYaxis()->SetTitleOffset(1.5);
128  frame->SetStats(kFALSE);
129  frame->SetTitle("");
130  frame->Draw();
131 
132  Draw(opt.Data());
133  c1->Modified();
134  c1->Update();
135  if (outDirName.Length() > 0) {
136  std::cout << "outDirName.Length() " << outDirName.Length() << std::endl;
137  TString outnamePDF = outDirName; outnamePDF += m_histName; outnamePDF += ".pdf";
138  TString outnameEPS = outDirName; outnameEPS += m_histName; outnameEPS += ".eps";
139  c1->SaveAs(outnamePDF.Data());
140  c1->SaveAs(outnameEPS.Data());
141  }
142 
143 }
144 
146 {
147 
148  unsigned int i;
149  for (auto& m : m_verticesMap) {
150  std::cout << " --> Aerogel ring : " << m.first << '\n';
151  for (i = 0; i < m.second.size(); i++) {
152  std::cout << " " << setw(15) << m.second[i].X() << setw(15) << m.second[i].Y() << std::endl;
153  }
154  }
155 
156 }
157 
159 {
160 
161  for (unsigned int i = 0; i < m_nTiles.size(); i++) {
162  std::vector<TVector2> vecTvec;
163  int nTiles = m_nTiles[i];
164  double rmin = m_tileRmin[i];
165  double rmax = m_tileRmax[i];
166  double lmin = 2 * TMath::Pi() * rmin / nTiles - m_aerogelTileGap;
167  double lmax = 2 * TMath::Pi() * rmax / nTiles - m_aerogelTileGap;
168  double phimin = lmin / rmin;
169  double phimax = lmax / rmax;
170  double x1 = rmin * TMath::Cos(phimin / 2.0);
171  double y1 = rmin * TMath::Sin(phimin / 2.0);
172  TVector2 v1(x1, y1);
173  vecTvec.push_back(v1);
174  double x2 = rmax * TMath::Cos(phimax / 2.0);
175  double y2 = rmax * TMath::Sin(phimax / 2.0);
176  TVector2 v2(x2, y2);
177  vecTvec.push_back(v2);
178  //Add circular points from outer radious (clockwise added)
179  if (m_nCircularPoints > 0) {
180  double dPhi = phimax / (m_nCircularPoints + 1);
181  for (int j = 0; j < m_nCircularPoints; j++) {
182  TVector2 v;
183  v.SetMagPhi(rmax, (phimax / 2.0 - dPhi * (j + 1)));
184  vecTvec.push_back(v);
185  }
186  }
187 
188  double x3 = x2;
189  double y3 = -y2;
190  TVector2 v3(x3, y3);
191  vecTvec.push_back(v3);
192  double x4 = x1;
193  double y4 = -y1;
194  TVector2 v4(x4, y4);
195  vecTvec.push_back(v4);
196  //Add circular points inner radious
197  if (m_nCircularPoints > 0) {
198  double dPhi = phimin / (m_nCircularPoints + 1);
199  for (int j = 0; j < m_nCircularPoints; j++) {
200  TVector2 v;
201  v.SetMagPhi(rmin, (-phimax / 2.0 + dPhi * (j + 1)));
202  vecTvec.push_back(v);
203  }
204  }
205 
206  vecTvec.push_back(v1);
207  m_verticesMap[i] = vecTvec;
208  }
209 
210 }
211 
213 {
214 
215  //Verbose level
216  m_verboseLevel = 0;
217  //Distance between aerogel tiles
218  m_aerogelTileGap = 1;
219  //Number of circular points
220  m_nCircularPoints = 10;
221  //Number of tiles per ring.
222  m_nTiles.push_back(22);
223  m_nTiles.push_back(28);
224  m_nTiles.push_back(34);
225  m_nTiles.push_back(40);
226  //Minimum radius of aerogel ring
227  m_tileRmin.push_back(441.0);
228  m_tileRmin.push_back(614.0);
229  m_tileRmin.push_back(787.0);
230  m_tileRmin.push_back(960.0);
231  //Maximum radius of aerogel ring
232  m_tileRmax.push_back(613.0);
233  m_tileRmax.push_back(786.0);
234  m_tileRmax.push_back(959.0);
235  m_tileRmax.push_back(1117.0);
236  for (unsigned int i = 0; i < m_nTiles.size(); i++) {
237  m_tileRcenter.push_back((m_tileRmax[i] + m_tileRmin[i]) / 2.0);
238  double r = m_tileRcenter[i];
239  double l = 2 * TMath::Pi() * r / m_nTiles[i] - m_aerogelTileGap;
240  double phi = l / r;
241  m_tileDeltaPhiCenter.push_back(phi);
243  }
244 
245 }
TString m_histTitle
Histogram title.
TString m_histName
Histogram name.
std::map< Int_t, std::vector< TVector2 > > m_verticesMap
Aerogel vertices map.
Int_t m_verboseLevel
Verbose level.
std::vector< double > m_tileRmax
Maximum radius of aerogel ring.
void SetUpVerticesMap()
Function for calculation vertices for one aerogel tile.
std::vector< Int_t > m_nTiles
Number of tiles per ring.
std::vector< double > m_tileRcenter
Center radius of aerogel ring.
void SetInitialParametersByDefault()
Function which set initial values of input parameters.
Int_t m_nCircularPoints
Number of circular points.
std::vector< double > m_tileDeltaPhiCenter
Angle opening (phi) of the aerogel tile measured between two rays (0.0,0.0 : and centre of the ring f...
std::vector< double > m_tileRmin
Minimum radius of aerogel ring.
void dumpVerticesMap()
Function to print vertices for one aerogel tile.
std::vector< double > m_aerogelAriGapDeltaPhiCenter
Angle opening (phi) of the air gap between aerogel tiles.
Int_t GetBinIDFromRingColumn(Int_t ring, Int_t column)
Function which return histogram bin id from ring and column id's.
double m_aerogelTileGap
Distance between aerogel tiles.
void makeRotation(double xold, double yold, double &xnew, double &ynew, double phi)
Function to rotate 2D point (x and y) around z axis by angle phi.
void DrawHisto(TString opt, TString outDirName)
Function to draw the histogram.
Abstract base class for different kinds of events.