BOSS 7.1.2
BESIII Offline Software System
Loading...
Searching...
No Matches
BesEmcEndGeometry.cc
Go to the documentation of this file.
1#include "G4ThreeVector.hh"
2
4#include "BesEmcParameter.hh"
5#include "G4ios.hh"
6
7#include <iomanip>
8#include <fstream>
9#include <sstream>
10#include <iostream>
11#include <assert.h>
12
13#include "G4SystemOfUnits.hh"
14
15using namespace std;
16using namespace CLHEP;
17
20
23
25{
27
28 WorldRmin1 = emcPara.GetWorldRmin1();
29 WorldRmax1 = emcPara.GetWorldRmax1();
30 WorldRmin2 = emcPara.GetWorldRmin2();//+5.*mm; //add 5mm to avoid warning
31 WorldRmax2 = emcPara.GetWorldRmax2();
32 WorldDz = emcPara.GetWorldDz();
33 WorldZPosition = emcPara.GetWorldZPosition();
34 CrystalLength = emcPara.GetCrystalLength();
35
36 SectorRmin1 = 489.27*mm;
37 SectorRmax1 = 880.*mm;
38 SectorRmin2 = 621.57*mm;
39 SectorRmax2 = 1113.53*mm;
40 SectorDz = WorldDz-44.*mm;
41 SectorZPosition = -18.*mm;
42
43 for(G4int i=0;i<6;i++)
44 cryNumInOneLayer[i] = emcPara.GetCryInOneLayer(i);
45 for(G4int i=0;i<5;i++)
46 pentaInOneSector[i] = emcPara.GetPentaInOneSector(i);
47
48 fTyvekThickness = emcPara.GetTyvekThickness();
49 fAlThickness = emcPara.GetAlThickness();
50 fMylarThickness = emcPara.GetMylarThickness();
51 totalThickness=fTyvekThickness+fAlThickness+fMylarThickness;
52
53 G4String ParaPath = getenv("EMCSIMROOT");
54 if(!ParaPath){
55 G4cout<<"BOOST environment not set!"<<G4endl;
56 exit(-1);
57 }
58 ParaPath += "/dat/EmcEndGeometry.dat";
59
60 ifstream fin;
61 fin.open(ParaPath);
62 assert(fin);
63 for(G4int i=0;i<30;i++)
64 for (G4int j=0;j<24;j++)
65 fin>>param[i][j];
66 for(G4int i=0;i<5;i++)
67 for (G4int j=0;j<6;j++)
68 fin>>penta[i][j];
69 fin.close();
70}
71
72
74{
75
76 ////////////////////////////////////////////////////////////////////////
77 // emc endcap crystal //
78 //////////////////////////////////////////////////////////////////////////
79 // 1 //
80 // 0 __ -- //
81 // -- \ //
82 // | \ \ //
83 // | \ __ - 2 //
84 // | \ -- | //
85 // | | 3 | //
86 // | | | //
87 // | | | //
88 // | | | //
89 // | | 5 | //
90 // 4 \ | | //
91 // \ | __ - 6 //
92 // \| -- //
93 // 7 //
94 /////////////////////////////////////////////////////////////////////////
95
97
98 G4int pentaNb=0;
99 for(G4int i=0;i<30;i++)
100 {
101 for (G4int j=0;j<8;j++)
102 {
103 //use 24 parameters to construct 8 point of one crystal
104 G4ThreeVector* pPnt = new G4ThreeVector(param[i][j],param[i][j+8],param[i][j+16]-WorldZPosition-SectorZPosition);
105 fPnt[i][j] = *pPnt;
106 }
107
108 if(i==5||i==6||i==14||i==15||i==16)
109 //if(0) //no pentagonal crystal now
110 {
111 for(G4int j=0;j<8;j++)
112 fPnt[30+pentaNb][j] = fPnt[i][j];
113
114 //compute the 5th point of the pentagonal crystal
115 G4double y0,y1,y4,y5;
116 G4ThreeVector v0,v1,v4,v5;
117 y0 = penta[pentaNb][0];
118 y1 = penta[pentaNb][1];
119 y4 = penta[pentaNb][2];
120 y5 = penta[pentaNb][3];
121 v0 = (y0-fPnt[i][3].getY())*(fPnt[i][0]-fPnt[i][3])/(fPnt[i][0].getY()-fPnt[i][3].getY())
122 +fPnt[i][3];
123 v1 = (y1-fPnt[i][2].getY())*(fPnt[i][1]-fPnt[i][2])/(fPnt[i][1].getY()-fPnt[i][2].getY())
124 +fPnt[i][2];
125 v4 = (y4-fPnt[i][7].getY())*(fPnt[i][4]-fPnt[i][7])/(fPnt[i][4].getY()-fPnt[i][7].getY())
126 +fPnt[i][7];
127 v5 = (y5-fPnt[i][6].getY())*(fPnt[i][5]-fPnt[i][6])/(fPnt[i][5].getY()-fPnt[i][6].getY())
128 +fPnt[i][6];
129
130 G4double x1,x5;
131 x1 = penta[pentaNb][4];
132 x5 = penta[pentaNb][5];
133 v1 = (x1-v0.getX())*(v1-v0)/(v1.getX()-v0.getX())+v0; //v1', the fifth point
134 v5 = (x5-v4.getX())*(v5-v4)/(v5.getX()-v4.getX())+v4; //v5'
135
136 fPnt[i][0] = v0;
137 fPnt[i][1] = v1;
138 fPnt[i][4] = v4;
139 fPnt[i][5] = v5;
140
141 fPnt[30+pentaNb][2] = v1;
142 fPnt[30+pentaNb][3] = v0;
143 fPnt[30+pentaNb][6] = v5;
144 fPnt[30+pentaNb][7] = v4;
145
146 pentaNb++;
147 }
148 }
149
150 //reflect point in one sector according to new design
151 G4ThreeVector temp[35][8];
152 for(G4int i=0;i<35;i++)
153 {
154 for (G4int j=0;j<8;j++)
155 {
156 temp[i][j]=fPnt[i][j];
157 fPnt[i][j].rotateZ(157.5*deg);
158 fPnt[i][j].setX(-fPnt[i][j].getX());
159 }
160
161 // point 0<-->3, 1<-->2, 4<-->7, 5<-->6
162 for (G4int j=0;j<8;j++)
163 {
164 if(j<2)
165 {
166 G4ThreeVector v=fPnt[i][j];
167 fPnt[i][j]=fPnt[i][3-j];
168 fPnt[i][3-j]=v;
169 }
170 else if(j>=4&&j<6)
171 {
172 G4ThreeVector v=fPnt[i][j];
173 fPnt[i][j]=fPnt[i][11-j];
174 fPnt[i][11-j]=v;
175 }
176 }
177 }
178
179 //exchange sequence in the same layer
180 //Exchange(0,3); Exchange(1,2); Exchange(4,7); Exchange(5,6);
181 //Exchange(8,12); Exchange(9,11); Exchange(13,17); Exchange(14,16);
182 //Exchange(18,23); Exchange(19,22); Exchange(20,21); Exchange(24,29); Exchange(25,28); Exchange(26,27);
183
184 Exchange(0,3); Exchange(1,2); Exchange(4,7); Exchange(5,31); Exchange(6,30);
185 Exchange(8,12); Exchange(9,11); Exchange(13,17); Exchange(14,34); Exchange(15,33); Exchange(16,32);
186 Exchange(18,23); Exchange(19,22); Exchange(20,21); Exchange(24,29); Exchange(25,28); Exchange(26,27);
187
188 /*for(G4int i=0;i<35;i++)
189 {
190 G4cout<<"crystal number: "<<i<<G4endl;
191 for(G4int j=0;j<8;j++)
192 G4cout<<fPnt[i][j]<<G4endl;
193 }*/
194
195
196 //compute the 6 crystal beside the 20mm gap
197 for(G4int i=0;i<35;i++)
198 {
199 for (G4int j=0;j<8;j++)
200 {
201 G4ThreeVector pPnt1 = temp[i][j];
202 fPnt1[i][j] = pPnt1.rotateZ(67.5*deg);
203 }
204 if((i==3)||(i==7)||(i==12)||(i==17)||(i==23)||(i==29))
205 {
206 fPnt1[i][0].setX(10);
207 fPnt1[i][1].setX(10);
208 fPnt1[i][4].setX(10);
209 fPnt1[i][5].setX(10);
210
211 G4double y0 = fPnt1[i][0].getY()+10*(fPnt1[i][3].getY()-fPnt1[i][0].getY())/fPnt1[i][3].getX();
212 G4double y1 = fPnt1[i][1].getY()+10*(fPnt1[i][2].getY()-fPnt1[i][1].getY())/fPnt1[i][2].getX();
213 G4double y4 = fPnt1[i][4].getY()+10*(fPnt1[i][7].getY()-fPnt1[i][4].getY())/fPnt1[i][7].getX();
214 G4double y5 = fPnt1[i][5].getY()+10*(fPnt1[i][6].getY()-fPnt1[i][5].getY())/fPnt1[i][6].getX();
215
216 G4double z0 = fPnt1[i][0].getZ()+10*(fPnt1[i][3].getZ()-fPnt1[i][0].getZ())/fPnt1[i][3].getX();
217 G4double z1 = fPnt1[i][1].getZ()+10*(fPnt1[i][2].getZ()-fPnt1[i][1].getZ())/fPnt1[i][2].getX();
218 G4double z4 = fPnt1[i][4].getZ()+10*(fPnt1[i][7].getZ()-fPnt1[i][4].getZ())/fPnt1[i][7].getX();
219 G4double z5 = fPnt1[i][5].getZ()+10*(fPnt1[i][6].getZ()-fPnt1[i][5].getZ())/fPnt1[i][6].getX();
220
221 fPnt1[i][0].setY(y0);
222 fPnt1[i][1].setY(y1);
223 fPnt1[i][4].setY(y4);
224 fPnt1[i][5].setY(y5);
225
226 fPnt1[i][0].setZ(z0);
227 fPnt1[i][1].setZ(z1);
228 fPnt1[i][4].setZ(z4);
229 fPnt1[i][5].setZ(z5);
230 }
231 }
232}
233
234void BesEmcEndGeometry::Exchange(G4int cry1, G4int cry2)
235{
236 G4ThreeVector v;
237 for(G4int i=0;i<8;i++)
238 {
239 v = fPnt[cry1][i];
240 fPnt[cry1][i] = fPnt[cry2][i];
241 fPnt[cry2][i] = v;
242 }
243}
244
245void BesEmcEndGeometry::ExchangeSector7(G4int cry1, G4int cry2)
246{
247 G4ThreeVector v;
248 for(G4int i=0;i<8;i++)
249 {
250 v = fPnt1[cry1][i];
251 fPnt1[cry1][i] = fPnt1[cry2][i];
252 fPnt1[cry2][i] = v;
253 }
254}
255
256//reflect x axis to construct sectors on the left
258{
259 for(G4int i=0;i<35;i++)
260 {
261 for(G4int j=0;j<8;j++)
262 {
263 fPnt1[i][j].setX(-fPnt1[i][j].getX());
264 }
265
266 // point 0<-->3, 1<-->2, 4<-->7, 5<-->6
267 for (G4int j=0;j<8;j++)
268 {
269 if(j<2)
270 {
271 G4ThreeVector v=fPnt1[i][j];
272 fPnt1[i][j]=fPnt1[i][3-j];
273 fPnt1[i][3-j]=v;
274 }
275 else if(j>=4&&j<6)
276 {
277 G4ThreeVector v=fPnt1[i][j];
278 fPnt1[i][j]=fPnt1[i][11-j];
279 fPnt1[i][11-j]=v;
280 }
281 }
282
283 //change the sequence of eight points according to the requirment of IrregBox
284 // point 0<-->1, 2<-->3, 4<-->5, 6<-->7
285 /*for(G4int j=0;j<8;j++)
286 {
287 if((j%2)==0)
288 {
289 G4ThreeVector v2=fPnt1[i][j];
290 fPnt1[i][j]=fPnt1[i][j+1];
291 fPnt1[i][j+1]=v2;
292 }
293 }*/
294 }
297 ExchangeSector7(9,11); ExchangeSector7(13,17); ExchangeSector7(14,34);
298 ExchangeSector7(15,33); ExchangeSector7(16,32); ExchangeSector7(18,23);
299 ExchangeSector7(19,22); ExchangeSector7(20,21); ExchangeSector7(24,29);
300 ExchangeSector7(25,28); ExchangeSector7(26,27);
301}
302
303void BesEmcEndGeometry::Zoom(const G4ThreeVector pos[8], const G4double factor)
304{
305 G4ThreeVector center1(0,0,0);
306 G4ThreeVector center2(0,0,0);
307 for(G4int i=0;i<8;i++)
308 zoomPoint[i]=G4ThreeVector(0,0,0);
309
310 for(G4int i=0;i<8;i++)
311 {
312 if(i<4) center1+=pos[i];
313 else center2+=pos[i];
314 }
315 center1/=4;
316 center2/=4;
317
318 for(G4int i=0;i<8;i++)
319 {
320 if(i<4) zoomPoint[i]=factor*pos[i]+(1-factor)*center1;
321 else zoomPoint[i]=factor*pos[i]+(1-factor)*center2;
322 }
323}
324
325//subtract the casing
326void BesEmcEndGeometry::ModifyForCasing(G4ThreeVector pos[8], G4int CryNb)
327{
328 G4ThreeVector center1(0,0,0); //the center of large surface
329 G4ThreeVector center2(0,0,0); //the center of small surface
330
331 const G4double dt=1e-5; //avoid overlap between crystal and casing
332
333 if(CryNb==5||CryNb==6||CryNb==14||CryNb==15||CryNb==16)
334 {
335 center1 = (pos[0]+pos[1])*(1-dt)/2+(pos[2]+pos[3])*dt/2;
336 center2 = (pos[4]+pos[5])*(1-dt)/2+(pos[6]+pos[7])*dt/2;
337 }
338 else if(CryNb>=30&&CryNb<35)
339 {
340 center1 = (pos[0]+pos[1])*dt/2+(pos[2]+pos[3])*(1-dt)/2;
341 center2 = (pos[4]+pos[5])*dt/2+(pos[6]+pos[7])*(1-dt)/2;
342 }
343 else
344 {
345 center1 = (pos[0]+pos[1]+pos[2]+pos[3])/4;
346 center2 = (pos[4]+pos[5]+pos[6]+pos[7])/4;
347 }
348
349 G4double r1=(pos[1]-center1).r();
350 G4double r2=(pos[2]-center1).r();
351 G4double r12=(pos[1]-pos[2]).r();
352 G4double theta=acos((r2*r2+r12*r12-r1*r1)/(2*r2*r12));
353 G4double h=r2*sin(theta); //distance from the center to the vertical side
354 G4double t1=totalThickness/h;
355
356 r1=(pos[5]-center2).r();
357 r2=(pos[6]-center2).r();
358 r12=(pos[5]-pos[6]).r();
359 theta=acos((r2*r2+r12*r12-r1*r1)/(2*r2*r12));
360 h=r2*sin(theta);
361 G4double t2=totalThickness/h;
362
363 for(G4int i=0;i<8;i++)
364 {
365 if(i<4)
366 {
367 cryPoint[i] = (1-t1)*pos[i]+t1*center1;
368 }
369 else
370 {
371 G4ThreeVector temp = (1-t2)*pos[i]+t2*center2;
372 cryPoint[i] = (1-totalThickness/CrystalLength)*temp+(totalThickness/CrystalLength)*pos[i-4];
373 }
374 //G4cout<<"casing="<<pos[i]<<"\t"<<"crystal="<<cryPoint[i]<<G4endl;
375 }
376}
377
379(const G4int partId, const G4int numTheta, const G4int numPhi)
380{
381 //ComputeParameters();
382 G4int sector=-1, cryNb=-1;
383 G4int leftFlag=0;
384 G4int downFlag=0;
385 G4int pentaFlag=0;
386 G4int flag=0;
387 G4double A1=0,a1=0,B1=0,b1=0,C1=0,c1=0,D1=0,d1=0,E1=0,e1=0; //dimension of pentagonal crystal
388 G4int m_numPhi=0;
389 G4ThreeVector position(0,0,0);
390 G4int cryInOneSector = cryNumInOneLayer[numTheta]/16; //number of crystal in one layer in one sector
391 G4ThreeVector pos[8];
392
393 if(partId==2) //west end
394 {
395 if(numPhi>=0&&numPhi<8*cryInOneSector)
396 m_numPhi=8*cryInOneSector-1-numPhi;
397 else if(numPhi>=8*cryInOneSector&&numPhi<16*cryInOneSector)
398 m_numPhi=16*cryInOneSector-1-numPhi;
399 }
400 else //east end
401 m_numPhi=numPhi;
402
403 if(numPhi>=4*cryInOneSector&&numPhi<5*cryInOneSector) //crystal in 4th sector
404 {
405 leftFlag = 1;
406 m_numPhi=8*cryInOneSector-1-numPhi;
407 }
408 else if(numPhi>=11*cryInOneSector&&numPhi<12*cryInOneSector) //crystal in 11th sector
409 {
410 leftFlag = 1;
411 downFlag = 1;
412 m_numPhi=numPhi-8*cryInOneSector;
413 }
414 if(numPhi>=12*cryInOneSector&&numPhi<13*cryInOneSector) //crystal in 12th sector
415 {
416 downFlag = 1;
417 m_numPhi=16*cryInOneSector-1-numPhi;
418 }
419
420 //translate numTheta and numPhi to sector and cryNb
421 G4int cryNbOffset = 0;
422 for(G4int i=0;i<numTheta;i++)
423 cryNbOffset += cryNumInOneLayer[i]/16;
424 if(cryInOneSector)
425 sector = m_numPhi/cryInOneSector;
426 cryNb = m_numPhi-cryInOneSector*sector+cryNbOffset;
427 sector += 3;
428 if(sector>15&&sector<=18)
429 sector -= 16;
430
431 //sector beside the gap
432 if(sector==6)
433 for(G4int i=0;i<8;i++)
434 pos[i]=fPnt1[cryNb][i];
435 else
436 for(G4int i=0;i<8;i++)
437 {
438 pos[i]=fPnt[cryNb][i];
439 pos[i].rotateZ(-67.5*deg+sector*22.5*deg);
440 }
441
442 //crystal dimension
443 Zoom(pos,0.999);
444 ModifyForCasing(zoomPoint,cryNb);
445 G4double A = (cryPoint[0]-cryPoint[3]).r();
446 G4double a = (cryPoint[4]-cryPoint[7]).r();
447 G4double B = (cryPoint[1]-cryPoint[2]).r();
448 G4double b = (cryPoint[5]-cryPoint[6]).r();
449 G4double C = (cryPoint[0]-cryPoint[1]).r();
450 G4double c = (cryPoint[4]-cryPoint[5]).r();
451 G4double D = (cryPoint[2]-cryPoint[3]).r();
452 G4double d = (cryPoint[6]-cryPoint[7]).r();
453
454 //reflect the axis according to the flag
455 for(G4int i=0;i<8;i++)
456 {
457 pos[i].setZ(pos[i].getZ()+WorldZPosition+SectorZPosition); //give the absolute z coordinate
458 if(leftFlag==1)
459 pos[i].setX(-pos[i].getX());
460 if(downFlag==1)
461 pos[i].setY(-pos[i].getY());
462 if(partId==2)
463 {
464 pos[i].setX(-pos[i].getX());
465 pos[i].setZ(-pos[i].getZ());
466 }
467 }
468
469 //compute the position
470 for(G4int j=4;j<8;j++)
471 position += pos[j];
472 position /= 4;
473
474 //compute pentagonal crystal
475 for(G4int i=0;i<5;i++)
476 {
477 if(cryNb==pentaInOneSector[i])
478 {
479 pentaFlag = 1;
480 G4ThreeVector penta[8];
481
482 //sector beside the gap
483 if(sector==6)
484 for(G4int j=0;j<8;j++)
485 penta[j]=fPnt1[30+i][j];
486 else
487 for(G4int j=0;j<8;j++)
488 {
489 penta[j]=fPnt[30+i][j];
490 penta[j].rotateZ(-67.5*deg+sector*22.5*deg);
491 }
492
493 //crystal dimension
494 ModifyForCasing(penta,30+i);
495 A1 = (cryPoint[1]-cryPoint[2]).r();
496 a1 = (cryPoint[5]-cryPoint[6]).r();
497 B1 = (cryPoint[1]-cryPoint[0]).r();
498 b1 = (cryPoint[5]-cryPoint[4]).r();
499 C1 = (cryPoint[0]-cryPoint[3]).r()+A;
500 c1 = (cryPoint[4]-cryPoint[7]).r()+a;
501 D1 = D;
502 d1 = d;
503 E1 = B;
504 e1 = b;
505
506 //reflect the axis according to the flag
507 for(G4int j=0;j<8;j++)
508 {
509 penta[j].setZ(penta[j].getZ()+WorldZPosition+SectorZPosition); //give the absolute z coordinate
510 if(leftFlag==1)
511 penta[j].setX(-penta[j].getX());
512 if(downFlag==1)
513 penta[j].setY(-penta[j].getY());
514 if(partId==2)
515 {
516 penta[j].setX(-penta[j].getX());
517 penta[j].setZ(-penta[j].getZ());
518 }
519 }
520
521 //compute the position
522 for(G4int j=4;j<8;j++)
523 {
524 if(j!=0&&j!=4)
525 position += pos[j];
526 if(j==0||j==1||j==4||j==5)
527 position += penta[j];
528 }
529 position /= 10;
530
531 flag = leftFlag+downFlag;
532 if(flag==1)
533 {
534 G4double temp1 = B1; B1 = D1; D1 = temp1;
535 temp1 = A1; A1 = E1; E1 = temp1;
536 temp1 = b1; b1 = d1; d1 = temp1;
537 temp1 = a1; a1 = e1; e1 = temp1;
538 }
539
540 break;
541 }
542 }
543
544 flag = leftFlag+downFlag+partId/2;
545 if(flag==1||flag==3)
546 {
547 G4double temp = C; C = D; D = temp;
548 temp = c; c = d; d = temp;
549 }
550
551 /*G4cout<<"##########################################################################"<<G4endl;
552 G4cout<<"###sector="<<sector<<",cryNb="<<cryNb<<",left flag="<<leftFlag<<",down flag="<<downFlag<<G4endl;
553 G4cout<<"partId="<<partId<<"\t"<<"theta number="<<numTheta<<"\t"<<"phi number="<<numPhi<<G4endl;
554 G4cout<<"crystal point:"<<G4endl;
555 for(G4int i=0;i<8;i++)
556 G4cout<<pos[i]<<G4endl;
557 G4cout<<"crystal position:"<<"\t"<<position<<"\t"<<"phi="<<position.phi()/deg<<G4endl;
558 G4cout<<"crystal dimension:"<<G4endl;
559 if(pentaFlag==1)
560 G4cout<<"A="<<A1<<",a="<<a1<<",B="<<B1<<",b="<<b1<<",C="<<C1<<",c="<<c1<<",D="<<D1<<",d="<<d1<<",E="<<E1<<",e="<<e1<<G4endl;
561 else
562 G4cout<<"A="<<A<<",a="<<a<<",B="<<B<<",b="<<b<<",C="<<C<<",c="<<c<<",D="<<D<<",d="<<d<<G4endl;
563 G4cout<<"##########################################################################"<<G4endl;*/
564
565 return position;
566}
double sin(const BesAngle a)
Definition BesAngle.h:210
Double_t e1
**********Class see also m_nmax DOUBLE PRECISION m_amel DOUBLE PRECISION m_x2 DOUBLE PRECISION m_alfinv DOUBLE PRECISION m_Xenph INTEGER m_KeyWtm INTEGER m_idyfs DOUBLE PRECISION m_zini DOUBLE PRECISION m_q2 DOUBLE PRECISION m_Wt_KF DOUBLE PRECISION m_WtCut INTEGER m_KFfin *COMMON c_KarLud $ !Input CMS energy[GeV] $ !CMS energy after beam spread beam strahlung[GeV] $ !Beam energy spread[GeV] $ !z boost due to beam spread $ !electron beam mass *ff pair spectrum $ !minimum v
Definition KarLud.h:35
TGraph2DErrors * dt
Definition McCor.cxx:45
***************************************************************************************Pseudo Class RRes *****************************************************************************************Parameters and physical constants **Maarten sept ************************************************************************DOUBLE PRECISION xsmu **************************************************************************PARTICLE DATA all others are from PDG *Only resonances with known widths into electron pairs are sept ************************************************************************C Declarations C
Definition RRes.h:29
void Exchange(G4int cry1, G4int cry2)
void ExchangeSector7(G4int cry1, G4int cry2)
void ModifyForCasing(G4ThreeVector pos[8], G4int CryNb)
void Zoom(const G4ThreeVector pos[8], const G4double factor)
G4ThreeVector ComputeDimAndPos(const G4int, const G4int, const G4int)
G4double GetWorldZPosition()
G4double GetCrystalLength()
G4int GetCryInOneLayer(G4int nb)
G4double GetWorldRmax1()
static BesEmcParameter & GetInstance()
G4double GetTyvekThickness()
G4double GetMylarThickness()
G4double GetWorldRmin1()
G4double GetWorldRmin2()
G4double GetWorldRmax2()
G4double GetWorldDz()
G4int GetPentaInOneSector(G4int nb)
G4double GetAlThickness()
const double b
Definition slope.cxx:9