Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GDMLParameterisation.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// G4GDMLParameterisation implementation
27//
28// Author: Zoltan Torzsok, November 2007
29// --------------------------------------------------------------------
30
34
35// --------------------------------------------------------------------
37{
38 return (G4int) parameterList.size();
39}
40
41// --------------------------------------------------------------------
43{
44 parameterList.push_back(newParameter);
45}
46
47// --------------------------------------------------------------------
48void G4GDMLParameterisation::ComputeTransformation(
49 const G4int index, G4VPhysicalVolume* physvol) const
50{
51 physvol->SetTranslation(parameterList[index].position);
52 physvol->SetRotation(parameterList[index].pRot);
53}
54
55// --------------------------------------------------------------------
56void G4GDMLParameterisation::ComputeDimensions(G4Box& box, const G4int index,
57 const G4VPhysicalVolume*) const
58{
59 box.SetXHalfLength(parameterList[index].dimension[0]);
60 box.SetYHalfLength(parameterList[index].dimension[1]);
61 box.SetZHalfLength(parameterList[index].dimension[2]);
62}
63
64// --------------------------------------------------------------------
65void G4GDMLParameterisation::ComputeDimensions(G4Trd& trd, const G4int index,
66 const G4VPhysicalVolume*) const
67{
68 trd.SetXHalfLength1(parameterList[index].dimension[0]);
69 trd.SetXHalfLength2(parameterList[index].dimension[1]);
70 trd.SetYHalfLength1(parameterList[index].dimension[2]);
71 trd.SetYHalfLength2(parameterList[index].dimension[3]);
72 trd.SetZHalfLength(parameterList[index].dimension[4]);
73}
74
75// --------------------------------------------------------------------
76void G4GDMLParameterisation::ComputeDimensions(G4Trap& trap, const G4int index,
77 const G4VPhysicalVolume*) const
78{
79 trap.SetAllParameters(parameterList[index].dimension[0], // Dz
80 parameterList[index].dimension[1], // Theta
81 parameterList[index].dimension[2], // Phi
82 parameterList[index].dimension[3], // Dy1
83 parameterList[index].dimension[4], // Dx1
84 parameterList[index].dimension[5], // Dx2
85 parameterList[index].dimension[6], // pAlp1,
86 parameterList[index].dimension[7], // pDy2,
87 parameterList[index].dimension[8], // pDx3,
88 parameterList[index].dimension[9], // pDx4,
89 parameterList[index].dimension[10]); // pAlp2
90}
91
92// --------------------------------------------------------------------
93void G4GDMLParameterisation::ComputeDimensions(G4Tubs& tubs, const G4int index,
94 const G4VPhysicalVolume*) const
95{
96 tubs.SetInnerRadius(parameterList[index].dimension[0]);
97 tubs.SetOuterRadius(parameterList[index].dimension[1]);
98 tubs.SetZHalfLength(parameterList[index].dimension[2]);
99 tubs.SetStartPhiAngle(parameterList[index].dimension[3]);
100 tubs.SetDeltaPhiAngle(parameterList[index].dimension[4]);
101}
102
103// --------------------------------------------------------------------
104void G4GDMLParameterisation::ComputeDimensions(G4Cons& cons, const G4int index,
105 const G4VPhysicalVolume*) const
106{
107 cons.SetInnerRadiusMinusZ(parameterList[index].dimension[0]);
108 cons.SetOuterRadiusMinusZ(parameterList[index].dimension[1]);
109 cons.SetInnerRadiusPlusZ(parameterList[index].dimension[2]);
110 cons.SetOuterRadiusPlusZ(parameterList[index].dimension[3]);
111 cons.SetZHalfLength(parameterList[index].dimension[4]);
112 cons.SetStartPhiAngle(parameterList[index].dimension[5]);
113 cons.SetDeltaPhiAngle(parameterList[index].dimension[6]);
114}
115
116// --------------------------------------------------------------------
117void G4GDMLParameterisation::ComputeDimensions(G4Sphere& sphere,
118 const G4int index,
119 const G4VPhysicalVolume*) const
120{
121 sphere.SetInnerRadius(parameterList[index].dimension[0]);
122 sphere.SetOuterRadius(parameterList[index].dimension[1]);
123 sphere.SetStartPhiAngle(parameterList[index].dimension[2]);
124 sphere.SetDeltaPhiAngle(parameterList[index].dimension[3]);
125 sphere.SetStartThetaAngle(parameterList[index].dimension[4]);
126 sphere.SetDeltaThetaAngle(parameterList[index].dimension[5]);
127}
128
129// --------------------------------------------------------------------
130void G4GDMLParameterisation::ComputeDimensions(G4Orb& orb, const G4int index,
131 const G4VPhysicalVolume*) const
132{
133 orb.SetRadius(parameterList[index].dimension[0]);
134}
135
136// --------------------------------------------------------------------
137void G4GDMLParameterisation::ComputeDimensions(G4Ellipsoid& ellipsoid,
138 const G4int index,
139 const G4VPhysicalVolume*) const
140{
141 ellipsoid.SetSemiAxis(parameterList[index].dimension[0],
142 parameterList[index].dimension[1],
143 parameterList[index].dimension[2]);
144 ellipsoid.SetZCuts(parameterList[index].dimension[3],
145 parameterList[index].dimension[4]);
146}
147
148// --------------------------------------------------------------------
149void G4GDMLParameterisation::ComputeDimensions(G4Torus& torus,
150 const G4int index,
151 const G4VPhysicalVolume*) const
152{
153 torus.SetAllParameters(parameterList[index].dimension[0], // pRmin
154 parameterList[index].dimension[1], // pRmax
155 parameterList[index].dimension[2], // pRtor
156 parameterList[index].dimension[3], // pSPhi
157 parameterList[index].dimension[4]); // pDPhi
158}
159
160// --------------------------------------------------------------------
161void G4GDMLParameterisation::ComputeDimensions(G4Para& para, const G4int index,
162 const G4VPhysicalVolume*) const
163{
164 para.SetXHalfLength(parameterList[index].dimension[0]);
165 para.SetYHalfLength(parameterList[index].dimension[1]);
166 para.SetZHalfLength(parameterList[index].dimension[2]);
167 para.SetAlpha(parameterList[index].dimension[3]);
168 para.SetTanAlpha(std::tan(parameterList[index].dimension[3]));
169 para.SetThetaAndPhi(parameterList[index].dimension[4],
170 parameterList[index].dimension[5]);
171}
172
173// --------------------------------------------------------------------
174void G4GDMLParameterisation::ComputeDimensions(G4Hype& hype, const G4int index,
175 const G4VPhysicalVolume*) const
176{
177 hype.SetInnerRadius(parameterList[index].dimension[0]);
178 hype.SetOuterRadius(parameterList[index].dimension[1]);
179 hype.SetZHalfLength(parameterList[index].dimension[4]);
180 hype.SetInnerStereo(parameterList[index].dimension[2]);
181 hype.SetOuterStereo(parameterList[index].dimension[3]);
182}
183
184// --------------------------------------------------------------------
185void G4GDMLParameterisation::ComputeDimensions(G4Polycone& pcone,
186 const G4int index,
187 const G4VPhysicalVolume*) const
188{
189 G4PolyconeHistorical origparam(*(pcone.GetOriginalParameters()));
190 origparam.Start_angle = parameterList[index].dimension[0];
191 origparam.Opening_angle = parameterList[index].dimension[1];
192 origparam.Num_z_planes = (G4int) parameterList[index].dimension[2];
193 G4int nZplanes = origparam.Num_z_planes;
194
195 for(G4int ii = 0; ii < nZplanes; ++ii)
196 {
197 origparam.Rmin[ii] = parameterList[index].dimension[3 + ii * 3];
198 origparam.Rmax[ii] = parameterList[index].dimension[4 + ii * 3];
199 origparam.Z_values[ii] = parameterList[index].dimension[5 + ii * 3];
200 }
201
202 pcone.SetOriginalParameters(&origparam); // copy values & transfer pointers
203 pcone.Reset(); // reset to new solid parameters
204}
205
206// --------------------------------------------------------------------
207void G4GDMLParameterisation::ComputeDimensions(G4Polyhedra& polyhedra,
208 const G4int index,
209 const G4VPhysicalVolume*) const
210{
211 G4PolyhedraHistorical origparam(*(polyhedra.GetOriginalParameters()));
212 origparam.Start_angle = parameterList[index].dimension[0];
213 origparam.Opening_angle = parameterList[index].dimension[1];
214 origparam.Num_z_planes = (G4int) parameterList[index].dimension[2];
215 origparam.numSide = (G4int) parameterList[index].dimension[3];
216
217 G4int nZplanes = origparam.Num_z_planes;
218
219 for(G4int ii = 0; ii < nZplanes; ++ii)
220 {
221 origparam.Rmin[ii] = parameterList[index].dimension[4 + ii * 3];
222 origparam.Rmax[ii] = parameterList[index].dimension[5 + ii * 3];
223 origparam.Z_values[ii] = parameterList[index].dimension[6 + ii * 3];
224 }
225 polyhedra.SetOriginalParameters(
226 &origparam); // copy values & transfer pointers
227 polyhedra.Reset(); // reset to new solid parameters
228}
int G4int
Definition: G4Types.hh:85
Definition: G4Box.hh:56
void SetZHalfLength(G4double dz)
Definition: G4Box.cc:172
void SetYHalfLength(G4double dy)
Definition: G4Box.cc:149
void SetXHalfLength(G4double dx)
Definition: G4Box.cc:125
Definition: G4Cons.hh:78
void SetInnerRadiusPlusZ(G4double Rmin2)
void SetZHalfLength(G4double newDz)
void SetStartPhiAngle(G4double newSPhi, G4bool trig=true)
void SetOuterRadiusMinusZ(G4double Rmax1)
void SetOuterRadiusPlusZ(G4double Rmax2)
void SetDeltaPhiAngle(G4double newDPhi)
void SetInnerRadiusMinusZ(G4double Rmin1)
void SetZCuts(G4double newzBottomCut, G4double newzTopCut)
void SetSemiAxis(G4double x, G4double y, G4double z)
void AddParameter(const PARAMETER &)
Definition: G4Hype.hh:69
void SetOuterStereo(G4double newOSte)
void SetOuterRadius(G4double newORad)
void SetZHalfLength(G4double newHLZ)
void SetInnerStereo(G4double newISte)
void SetInnerRadius(G4double newIRad)
Definition: G4Orb.hh:56
void SetRadius(G4double newRmax)
Definition: G4Para.hh:79
void SetYHalfLength(G4double val)
void SetTanAlpha(G4double val)
void SetThetaAndPhi(double pTheta, double pPhi)
void SetZHalfLength(G4double val)
void SetXHalfLength(G4double val)
void SetAlpha(G4double alpha)
void SetOriginalParameters(G4PolyconeHistorical *pars)
G4PolyconeHistorical * GetOriginalParameters() const
G4bool Reset()
Definition: G4Polycone.cc:436
void SetOriginalParameters(G4PolyhedraHistorical *pars)
G4PolyhedraHistorical * GetOriginalParameters() const
G4bool Reset()
Definition: G4Polyhedra.cc:463
void SetDeltaPhiAngle(G4double newDphi)
void SetStartThetaAngle(G4double newSTheta)
void SetOuterRadius(G4double newRmax)
void SetDeltaThetaAngle(G4double newDTheta)
void SetInnerRadius(G4double newRMin)
void SetStartPhiAngle(G4double newSphi, G4bool trig=true)
void SetAllParameters(G4double pRmin, G4double pRmax, G4double pRtor, G4double pSPhi, G4double pDPhi)
Definition: G4Torus.cc:81
void SetAllParameters(G4double pDz, G4double pTheta, G4double pPhi, G4double pDy1, G4double pDx1, G4double pDx2, G4double pAlp1, G4double pDy2, G4double pDx3, G4double pDx4, G4double pAlp2)
Definition: G4Trap.cc:282
Definition: G4Trd.hh:63
void SetYHalfLength2(G4double val)
void SetXHalfLength1(G4double val)
void SetYHalfLength1(G4double val)
void SetXHalfLength2(G4double val)
void SetZHalfLength(G4double val)
Definition: G4Tubs.hh:75
void SetDeltaPhiAngle(G4double newDPhi)
void SetStartPhiAngle(G4double newSPhi, G4bool trig=true)
void SetInnerRadius(G4double newRMin)
void SetOuterRadius(G4double newRMax)
void SetZHalfLength(G4double newDz)
void SetTranslation(const G4ThreeVector &v)
void SetRotation(G4RotationMatrix *)