Belle II Software development
ARICHAerogelHist Class Reference

Base class for geometry parameters. More...

#include <ARICHAerogelHist.h>

Inheritance diagram for ARICHAerogelHist:

Public Member Functions

 ARICHAerogelHist ()
 Default constructor.
 
 ~ARICHAerogelHist ()
 Default destructor.
 
 ARICHAerogelHist (const char *name, const char *title)
 Constructor with name, title.
 
Int_t GetBinIDFromRingColumn (Int_t ring, Int_t column)
 Function which return histogram bin id from ring and column id's.
 
void DrawHisto (TString opt, TString outDirName)
 Function to draw the histogram.
 

Protected Member Functions

void SetInitialParametersByDefault ()
 Function which set initial values of input parameters.
 
void SetUpVerticesMap ()
 Function for calculation vertices for one aerogel tile.
 
void dumpVerticesMap ()
 Function to print vertices for one aerogel tile.
 
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.
 
 ClassDef (ARICHAerogelHist, 2)
 ClassDef.
 

Protected Attributes

std::vector< Int_t > m_nTiles
 Number of tiles per ring.
 
std::vector< double > m_tileRmin
 Minimum radius of aerogel ring.
 
std::vector< double > m_tileRmax
 Maximum radius of aerogel ring.
 
std::vector< double > m_tileRcenter
 Center radius of aerogel ring.
 
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 from left/right).
 
std::vector< double > m_aerogelAriGapDeltaPhiCenter
 Angle opening (phi) of the air gap between aerogel tiles.
 
Int_t m_verboseLevel
 Verbose level.
 
Int_t m_nCircularPoints
 Number of circular points.
 
double m_aerogelTileGap
 Distance between aerogel tiles.
 
std::map< Int_t, std::vector< ROOT::Math::XYVector > > m_verticesMap
 Aerogel vertices map.
 
TString m_histName
 Histogram name.
 
TString m_histTitle
 Histogram title.
 

Detailed Description

Base class for geometry parameters.

Definition at line 28 of file ARICHAerogelHist.h.

Constructor & Destructor Documentation

◆ ARICHAerogelHist() [1/2]

ARICHAerogelHist ( )
inline

Default constructor.

Definition at line 35 of file ARICHAerogelHist.h.

35{};

◆ ~ARICHAerogelHist()

~ARICHAerogelHist ( )
inline

Default destructor.

Definition at line 40 of file ARICHAerogelHist.h.

40{};

◆ ARICHAerogelHist() [2/2]

ARICHAerogelHist ( const char *  name,
const char *  title 
)

Constructor with name, title.

Parameters
namename.
titletitle.

Definition at line 29 of file ARICHAerogelHist.cc.

29 : TH2Poly(), m_histName(name), m_histTitle(title)
30{
31 SetName(m_histName.Data());
32 SetTitle(m_histTitle.Data());
35
36 if (m_verboseLevel > 0) {
37 std::cout << setw(10) << "ring" << setw(10) << "nTiles" << setw(10) << "tileRmin" << setw(10) << "tileRmax" << std::endl;
38 for (unsigned int i = 0; i < m_nTiles.size(); i++)
39 std::cout << setw(10) << i + 1 << setw(10) << m_nTiles[i] << setw(10) << m_tileRmin[i] << setw(10) << m_tileRmax[i] << std::endl;
41 }
42
43 unsigned int n;
44 double* x;
45 double* y;
46
47 double xold;
48 double yold;
49 double xnew;
50 double ynew;
51 double phi;
52
53 //Add "poly bins" to the histogram
54 //Loop over rings with different radiuses
55 unsigned int iR = 0;
56 for (auto& m : m_verticesMap) {
57 //Loop aerogel tiles (bins) within one ring.
58 double phi0 = m_tileDeltaPhiCenter[iR] / 2.0;
59 double deltaPhi = m_tileDeltaPhiCenter[iR] + m_aerogelAriGapDeltaPhiCenter[iR];
60 for (Int_t j = 0; j < m_nTiles[iR]; j++) {
61 phi = phi0 + deltaPhi * j;
62 n = m.second.size();
63 x = new double [n];
64 y = new double [n];
65 //Loop over polygonal points which defines bins.
66 for (unsigned int i = 0; i < n; i++) {
67 xold = m.second[i].X();
68 yold = m.second[i].Y();
69 makeRotation(xold, yold, xnew, ynew, phi);
70 x[i] = xnew;
71 y[i] = ynew;
72 }
73 AddBin(n, x, y);
74 delete []x;
75 delete []y;
76 }
77 iR++;
78 }// for (auto& m: m_verticesMap) {
79
80}
TString m_histTitle
Histogram title.
TString m_histName
Histogram name.
Int_t m_verboseLevel
Verbose level.
std::vector< double > m_tileRmax
Maximum radius of aerogel ring.
std::map< Int_t, std::vector< ROOT::Math::XYVector > > m_verticesMap
Aerogel vertices map.
void SetUpVerticesMap()
Function for calculation vertices for one aerogel tile.
std::vector< Int_t > m_nTiles
Number of tiles per ring.
void SetInitialParametersByDefault()
Function which set initial values of input parameters.
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.
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.

Member Function Documentation

◆ DrawHisto()

void DrawHisto ( TString  opt = "ZCOLOT text same",
TString  outDirName = "./" 
)

Function to draw the histogram.

Parameters
optdraw option string default value : "ZCOLOT text same".
outDirNamename of epe and pdf to save the plots.

Definition at line 104 of file ARICHAerogelHist.cc.

105{
106
107 TCanvas* c1 = new TCanvas("c1", "c1", 1000, 1000);
108 c1->SetTitle(m_histTitle.Data());
109 c1->SetName(m_histName.Data());
110 c1->SetRightMargin(0.17);
111 c1->SetLeftMargin(0.12);
112 c1->SetTopMargin(0.15);
113 c1->SetBottomMargin(0.15);
114
115 TH2F* frame = new TH2F("h2", "h2", 40, -1200, 1200, 40, -1200, 1200);
116 //TH2F *frame = new TH2F("h2","h2",40, 400, 650, 40, -50, 200);
117 frame->GetXaxis()->SetTitle("x (mm)");
118 frame->GetYaxis()->SetTitle("y (mm)");
119 frame->GetXaxis()->CenterTitle();
120 frame->GetYaxis()->CenterTitle();
121 frame->GetYaxis()->SetTitleOffset(1.5);
122 frame->SetStats(kFALSE);
123 frame->SetTitle("");
124 frame->Draw();
125
126 Draw(opt.Data());
127 c1->Modified();
128 c1->Update();
129 if (outDirName.Length() > 0) {
130 std::cout << "outDirName.Length() " << outDirName.Length() << std::endl;
131 TString outnamePDF = outDirName; outnamePDF += m_histName; outnamePDF += ".pdf";
132 TString outnameEPS = outDirName; outnameEPS += m_histName; outnameEPS += ".eps";
133 c1->SaveAs(outnamePDF.Data());
134 c1->SaveAs(outnameEPS.Data());
135 }
136
137}

◆ dumpVerticesMap()

void dumpVerticesMap ( )
protected

Function to print vertices for one aerogel tile.

Definition at line 139 of file ARICHAerogelHist.cc.

140{
141
142 unsigned int i;
143 for (auto& m : m_verticesMap) {
144 std::cout << " --> Aerogel ring : " << m.first << '\n';
145 for (i = 0; i < m.second.size(); i++) {
146 std::cout << " " << setw(15) << m.second[i].X() << setw(15) << m.second[i].Y() << std::endl;
147 }
148 }
149
150}

◆ GetBinIDFromRingColumn()

Int_t GetBinIDFromRingColumn ( Int_t  ring,
Int_t  column 
)

Function which return histogram bin id from ring and column id's.

Parameters
ringringID number.
columncolumnID number.

Definition at line 82 of file ARICHAerogelHist.cc.

83{
84
85 Int_t binID = 0;
86 unsigned int ringu = ring;
87 //cout<<" ewewe "<<m_nTiles[ring-1]<<endl;
88 if (ringu > m_nTiles.size() || ringu < 1)
89 return -999;
90 if (column > m_nTiles[ring - 1] || column < 1)
91 return -999;
92 for (Int_t i = 1; i < ring; i++)
93 binID += m_nTiles[i - 1];
94 return binID + column;
95
96}

◆ makeRotation()

void makeRotation ( double  xold,
double  yold,
double &  xnew,
double &  ynew,
double  phi 
)
protected

Function to rotate 2D point (x and y) around z axis by angle phi.

Parameters
xoldold x coordinate.
yoldold y coordinate.
xnewnew x coordinate.
ynewnew y coordinate.
phiroration angle.

Definition at line 98 of file ARICHAerogelHist.cc.

99{
100 xnew = xold * std::cos(phi) - yold * std::sin(phi);
101 ynew = xold * std::sin(phi) + yold * std::cos(phi);
102}

◆ SetInitialParametersByDefault()

void SetInitialParametersByDefault ( )
protected

Function which set initial values of input parameters.

Definition at line 206 of file ARICHAerogelHist.cc.

207{
208
209 //Verbose level
210 m_verboseLevel = 0;
211 //Distance between aerogel tiles
213 //Number of circular points
215 //Number of tiles per ring.
216 m_nTiles.push_back(22);
217 m_nTiles.push_back(28);
218 m_nTiles.push_back(34);
219 m_nTiles.push_back(40);
220 //Minimum radius of aerogel ring
221 m_tileRmin.push_back(441.0);
222 m_tileRmin.push_back(614.0);
223 m_tileRmin.push_back(787.0);
224 m_tileRmin.push_back(960.0);
225 //Maximum radius of aerogel ring
226 m_tileRmax.push_back(613.0);
227 m_tileRmax.push_back(786.0);
228 m_tileRmax.push_back(959.0);
229 m_tileRmax.push_back(1117.0);
230 for (unsigned int i = 0; i < m_nTiles.size(); i++) {
231 m_tileRcenter.push_back((m_tileRmax[i] + m_tileRmin[i]) / 2.0);
232 double r = m_tileRcenter[i];
233 double l = 2 * TMath::Pi() * r / m_nTiles[i] - m_aerogelTileGap;
234 double phi = l / r;
235 m_tileDeltaPhiCenter.push_back(phi);
237 }
238
239}
std::vector< double > m_tileRcenter
Center radius of aerogel ring.
Int_t m_nCircularPoints
Number of circular points.
double m_aerogelTileGap
Distance between aerogel tiles.

◆ SetUpVerticesMap()

void SetUpVerticesMap ( )
protected

Function for calculation vertices for one aerogel tile.

Definition at line 152 of file ARICHAerogelHist.cc.

153{
154
155 for (unsigned int i = 0; i < m_nTiles.size(); i++) {
156 std::vector<ROOT::Math::XYVector> vecTvec;
157 int nTiles = m_nTiles[i];
158 double rmin = m_tileRmin[i];
159 double rmax = m_tileRmax[i];
160 double lmin = 2 * TMath::Pi() * rmin / nTiles - m_aerogelTileGap;
161 double lmax = 2 * TMath::Pi() * rmax / nTiles - m_aerogelTileGap;
162 double phimin = lmin / rmin;
163 double phimax = lmax / rmax;
164 double x1 = rmin * TMath::Cos(phimin / 2.0);
165 double y1 = rmin * TMath::Sin(phimin / 2.0);
166 ROOT::Math::XYVector v1(x1, y1);
167 vecTvec.push_back(v1);
168 double x2 = rmax * TMath::Cos(phimax / 2.0);
169 double y2 = rmax * TMath::Sin(phimax / 2.0);
170 ROOT::Math::XYVector v2(x2, y2);
171 vecTvec.push_back(v2);
172 //Add circular points from outer radious (clockwise added)
173 if (m_nCircularPoints > 0) {
174 double dPhi = phimax / (m_nCircularPoints + 1);
175 for (int j = 0; j < m_nCircularPoints; j++) {
176 const double phi = phimax / 2.0 - dPhi * (j + 1);
177 ROOT::Math::XYVector v(rmax * std::cos(phi), rmax * std::sin(phi));
178 vecTvec.push_back(v);
179 }
180 }
181
182 double x3 = x2;
183 double y3 = -y2;
184 ROOT::Math::XYVector v3(x3, y3);
185 vecTvec.push_back(v3);
186 double x4 = x1;
187 double y4 = -y1;
188 ROOT::Math::XYVector v4(x4, y4);
189 vecTvec.push_back(v4);
190 //Add circular points inner radious
191 if (m_nCircularPoints > 0) {
192 double dPhi = phimin / (m_nCircularPoints + 1);
193 for (int j = 0; j < m_nCircularPoints; j++) {
194 const double phi = -phimax / 2.0 + dPhi * (j + 1);
195 ROOT::Math::XYVector v(rmax * std::cos(phi), rmax * std::sin(phi));
196 vecTvec.push_back(v);
197 }
198 }
199
200 vecTvec.push_back(v1);
201 m_verticesMap[i] = vecTvec;
202 }
203
204}
const std::vector< double > v3
MATLAB generated random vector.
const std::vector< double > v2
MATLAB generated random vector.
const std::vector< double > v1
MATLAB generated random vector.

Member Data Documentation

◆ m_aerogelAriGapDeltaPhiCenter

std::vector<double> m_aerogelAriGapDeltaPhiCenter
protected

Angle opening (phi) of the air gap between aerogel tiles.

Measured between ray (0.0,0.0 : and centre of the ring from left/right).

Definition at line 97 of file ARICHAerogelHist.h.

◆ m_aerogelTileGap

double m_aerogelTileGap
protected

Distance between aerogel tiles.

Definition at line 100 of file ARICHAerogelHist.h.

◆ m_histName

TString m_histName
protected

Histogram name.

Definition at line 105 of file ARICHAerogelHist.h.

◆ m_histTitle

TString m_histTitle
protected

Histogram title.

Definition at line 106 of file ARICHAerogelHist.h.

◆ m_nCircularPoints

Int_t m_nCircularPoints
protected

Number of circular points.

Definition at line 99 of file ARICHAerogelHist.h.

◆ m_nTiles

std::vector<Int_t> m_nTiles
protected

Number of tiles per ring.

Definition at line 90 of file ARICHAerogelHist.h.

◆ m_tileDeltaPhiCenter

std::vector<double> m_tileDeltaPhiCenter
protected

Angle opening (phi) of the aerogel tile measured between two rays (0.0,0.0 : and centre of the ring from left/right).

Definition at line 95 of file ARICHAerogelHist.h.

◆ m_tileRcenter

std::vector<double> m_tileRcenter
protected

Center radius of aerogel ring.

Definition at line 93 of file ARICHAerogelHist.h.

◆ m_tileRmax

std::vector<double> m_tileRmax
protected

Maximum radius of aerogel ring.

Definition at line 92 of file ARICHAerogelHist.h.

◆ m_tileRmin

std::vector<double> m_tileRmin
protected

Minimum radius of aerogel ring.

Definition at line 91 of file ARICHAerogelHist.h.

◆ m_verboseLevel

Int_t m_verboseLevel
protected

Verbose level.

Definition at line 98 of file ARICHAerogelHist.h.

◆ m_verticesMap

std::map<Int_t, std::vector<ROOT::Math::XYVector> > m_verticesMap
protected

Aerogel vertices map.

Definition at line 103 of file ARICHAerogelHist.h.


The documentation for this class was generated from the following files: