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 26 of file ARICHAerogelHist.cc.

26 : TH2Poly(), m_histName(name), m_histTitle(title)
27{
28 SetName(m_histName.Data());
29 SetTitle(m_histTitle.Data());
32
33 if (m_verboseLevel > 0) {
34 std::cout << setw(10) << "ring" << setw(10) << "nTiles" << setw(10) << "tileRmin" << setw(10) << "tileRmax" << std::endl;
35 for (unsigned int i = 0; i < m_nTiles.size(); i++)
36 std::cout << setw(10) << i + 1 << setw(10) << m_nTiles[i] << setw(10) << m_tileRmin[i] << setw(10) << m_tileRmax[i] << std::endl;
38 }
39
40 unsigned int n;
41 double* x;
42 double* y;
43
44 double xold;
45 double yold;
46 double xnew;
47 double ynew;
48 double phi;
49
50 //Add "poly bins" to the histogram
51 //Loop over rings with different radiuses
52 unsigned int iR = 0;
53 for (auto& m : m_verticesMap) {
54 //Loop aerogel tiles (bins) within one ring.
55 double phi0 = m_tileDeltaPhiCenter[iR] / 2.0;
56 double deltaPhi = m_tileDeltaPhiCenter[iR] + m_aerogelAriGapDeltaPhiCenter[iR];
57 for (Int_t j = 0; j < m_nTiles[iR]; j++) {
58 phi = phi0 + deltaPhi * j;
59 n = m.second.size();
60 x = new double [n];
61 y = new double [n];
62 //Loop over polygonal points which defines bins.
63 for (unsigned int i = 0; i < n; i++) {
64 xold = m.second[i].X();
65 yold = m.second[i].Y();
66 makeRotation(xold, yold, xnew, ynew, phi);
67 x[i] = xnew;
68 y[i] = ynew;
69 }
70 AddBin(n, x, y);
71 delete []x;
72 delete []y;
73 }
74 iR++;
75 }// for (auto& m: m_verticesMap) {
76
77}
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 101 of file ARICHAerogelHist.cc.

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

◆ dumpVerticesMap()

void dumpVerticesMap ( )
protected

Function to print vertices for one aerogel tile.

Definition at line 136 of file ARICHAerogelHist.cc.

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

◆ 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 79 of file ARICHAerogelHist.cc.

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

◆ 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.
phirotation angle.

Definition at line 95 of file ARICHAerogelHist.cc.

96{
97 xnew = xold * std::cos(phi) - yold * std::sin(phi);
98 ynew = xold * std::sin(phi) + yold * std::cos(phi);
99}

◆ SetInitialParametersByDefault()

void SetInitialParametersByDefault ( )
protected

Function which set initial values of input parameters.

Definition at line 203 of file ARICHAerogelHist.cc.

204{
205
206 //Verbose level
207 m_verboseLevel = 0;
208 //Distance between aerogel tiles
210 //Number of circular points
212 //Number of tiles per ring.
213 m_nTiles.push_back(22);
214 m_nTiles.push_back(28);
215 m_nTiles.push_back(34);
216 m_nTiles.push_back(40);
217 //Minimum radius of aerogel ring
218 m_tileRmin.push_back(441.0);
219 m_tileRmin.push_back(614.0);
220 m_tileRmin.push_back(787.0);
221 m_tileRmin.push_back(960.0);
222 //Maximum radius of aerogel ring
223 m_tileRmax.push_back(613.0);
224 m_tileRmax.push_back(786.0);
225 m_tileRmax.push_back(959.0);
226 m_tileRmax.push_back(1117.0);
227 for (unsigned int i = 0; i < m_nTiles.size(); i++) {
228 m_tileRcenter.push_back((m_tileRmax[i] + m_tileRmin[i]) / 2.0);
229 double r = m_tileRcenter[i];
230 double l = 2 * TMath::Pi() * r / m_nTiles[i] - m_aerogelTileGap;
231 double phi = l / r;
232 m_tileDeltaPhiCenter.push_back(phi);
234 }
235
236}
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 149 of file ARICHAerogelHist.cc.

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