FDsimG4VCorrectorRing.cc
Go to the documentation of this file.
2 
3 #include "G4Polycone.hh"
4 
5 #include "G4LogicalVolume.hh"
6 #include "G4PVPlacement.hh"
7 #include "G4ThreeVector.hh"
8 
9 using namespace TelescopeSimulatorLX ;
10 
11 FDsimG4VCorrectorRing::FDsimG4VCorrectorRing(G4double rmin, G4double rmax, G4double minPhi, G4double deltaPhi, G4double thickness, G4Material* material, G4int npoints)
12 {
13 
14  fMinRadius = rmin ;
15  fMaxRadius = rmax ;
16 
17  fStartPhi = minPhi ;
18  fDeltaPhi = deltaPhi ;
19  fThickness = thickness ;
20 
21  fNpoints = npoints ;
22 
23  fMaterial = material ;
24  fLogicalVol = NULL ;
25 
26 }
27 
29 
30  delete [] fRminVec ;
31  delete [] fRmaxVec ;
32  delete [] fZVec ;
33 
34  fRminVec = NULL ;
35  fRmaxVec = NULL ;
36  fZVec = NULL ;
37 
38 }
39 
40 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
42 
43 
44  static const G4int NUMZPLANES = fNpoints ;
45  static const G4double z0 = GetSagita(fMinRadius) ;
46 
47  fRminVec = new G4double[NUMZPLANES] ;
48  fRmaxVec = new G4double[NUMZPLANES] ;
49  fZVec = new G4double[NUMZPLANES] ;
50 
51  G4double SagMin = GetSagita(fMinRadius) ;
52  G4double SagMax = GetSagita(fMaxRadius) ;
53 
54  G4double Sagita ;
55 
56  for (G4int iz=0 ; iz < NUMZPLANES-1 ; iz++){
57  fRmaxVec[iz] = fMaxRadius ;
58  Sagita = SagMax - (iz)*((SagMax-SagMin)/(NUMZPLANES-2)) ;
59  fRminVec[iz] = GetRadius(Sagita) ;
60  fZVec[iz] = Sagita - z0 ;
61 
62  }
63 
64 
65  fRminVec[NUMZPLANES-1] = fMinRadius ;
66  fRmaxVec[NUMZPLANES-1] = fMaxRadius ;
67  fZVec[NUMZPLANES-1] = fZVec[NUMZPLANES-2] - fThickness ;
68 
69 
70  G4Polycone *ring_solid
71  = new G4Polycone("CorrectorRing",fStartPhi,fDeltaPhi,NUMZPLANES,fZVec,fRminVec,fRmaxVec);
72 
74  = new G4LogicalVolume(ring_solid,fMaterial,"CorrectorRing",0,0,0);
75 
76 
77 }
78 
79 
80 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
82 
83  static const G4int NUMZPLANES = fNpoints ;
84 
85  static const G4double z0 = GetSagita(fMinRadius) ;
86 
87  fRminVec = new G4double[NUMZPLANES] ;
88  fRmaxVec = new G4double[NUMZPLANES] ;
89  fZVec = new G4double[NUMZPLANES] ;
90 
91  for (G4int iRad=0 ; iRad < NUMZPLANES-1 ; iRad++){
92  fRmaxVec[iRad] = fMaxRadius ;
93  fRminVec[iRad] = fMaxRadius - (iRad)*((fMaxRadius-fMinRadius)/(NUMZPLANES-2)) ;
94  fZVec[iRad] = GetSagita(fRminVec[iRad]) - z0;
95  }
96 
97 
98  fRminVec[NUMZPLANES-1] = fMinRadius ;
99  fRmaxVec[NUMZPLANES-1] = fMaxRadius ;
100  fZVec[NUMZPLANES-1] = fZVec[NUMZPLANES-2] - fThickness ;
101 
102 
103  G4Polycone *ring_solid
104  = new G4Polycone("CorrectorRing",fStartPhi,fDeltaPhi,NUMZPLANES,fZVec,fRminVec,fRmaxVec);
105 
106  fLogicalVol
107  = new G4LogicalVolume(ring_solid,fMaterial,"CorrectorRing",0,0,0);
108 
109 
110 }
111 
112 
113 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
114 
116 
117  G4cerr << " ================================================================ " << G4endl ;
118  G4cerr << " = Corrector Lens parameters " << G4endl ;
119  G4cerr << " = " << G4endl ;
120  G4cerr << " = Inner radius : " << fMinRadius/mm << " (mm)" << G4endl;
121  G4cerr << " = Outer radius : " << fMaxRadius/mm << " (mm)" << G4endl;
122  G4cerr << " = " << G4endl ;
123  G4cerr << " = Azimuthal span : " << fStartPhi*180./pi << " < Phi < " << (fStartPhi+fDeltaPhi)*180.0/pi << " (deg.) " << G4endl ;
124  G4cerr << " = " << G4endl ;
125  G4cerr << " = Thickness : " << fThickness/mm << " (mm) " << G4endl ;
126  G4cerr << " = " << G4endl ;
127  G4cerr << " = Material : " << fMaterial->GetName() << G4endl ;
128  G4cerr << " = " << G4endl ;
129  G4cerr << " ================================================================ " << G4endl ;
130 
131 #ifdef VERBOSE_USE
132  for (G4int iz=0 ; iz < 10 ; iz++){
133  G4cerr << "Rmin,Rmax,z " << fRminVec[iz]/mm << " " << fRmaxVec[iz]/mm << " " <<
134  fZVec[iz]/mm << G4endl ;
135  }
136 
137  G4cerr << "." << G4endl << "." << G4endl << "." << G4endl << "." << G4endl ;
138 
139  for (G4int iz= (fNpoints<10) ? fNpoints : fNpoints-10 ; iz < fNpoints ; iz++){
140  G4cerr << "Rmin,Rmax,z " << fRminVec[iz]/mm << " " << fRmaxVec[iz]/mm << " " <<
141  fZVec[iz]/mm << G4endl ;
142  }
143 
144 #endif
145 
146 }
constexpr double mm
Definition: AugerUnits.h:113
virtual G4double GetSagita(G4double)=0
virtual G4double GetRadius(G4double)=0
FDsimG4VCorrectorRing(G4double rMin, G4double rMax, G4double minPhi, G4double deltaPhi, G4double Thickness, G4Material *Material, G4int Npoints)

, generated on Tue Sep 26 2023.