Geant4 9.6.0
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// $Id$
27//
28// class G4GDMLParameterisation Implementation
29//
30// History:
31// - Created. Zoltan Torzsok, November 2007
32// -------------------------------------------------------------------------
33
35
37{
38 return (G4int)parameterList.size();
39}
40
42{
43 parameterList.push_back(newParameter);
44}
45
46void G4GDMLParameterisation::
47ComputeTransformation(const G4int index,G4VPhysicalVolume* physvol) const
48{
49 physvol->SetTranslation(parameterList[index].position);
50 physvol->SetRotation(parameterList[index].pRot);
51}
52
53void G4GDMLParameterisation::
54ComputeDimensions(G4Box& box,const G4int index,const G4VPhysicalVolume*) const
55{
56 box.SetXHalfLength(parameterList[index].dimension[0]);
57 box.SetYHalfLength(parameterList[index].dimension[1]);
58 box.SetZHalfLength(parameterList[index].dimension[2]);
59}
60
61void G4GDMLParameterisation::
62ComputeDimensions(G4Trd& trd,const G4int index,const G4VPhysicalVolume*) const
63{
64 trd.SetXHalfLength1(parameterList[index].dimension[0]);
65 trd.SetXHalfLength2(parameterList[index].dimension[1]);
66 trd.SetYHalfLength1(parameterList[index].dimension[2]);
67 trd.SetYHalfLength2(parameterList[index].dimension[3]);
68 trd.SetZHalfLength(parameterList[index].dimension[4]);
69}
70
71void G4GDMLParameterisation::
72ComputeDimensions(G4Trap& trap,const G4int index,const G4VPhysicalVolume*) const
73{
74 trap.SetAllParameters(parameterList[index].dimension[0], // Dz
75 parameterList[index].dimension[1], // Theta
76 parameterList[index].dimension[2], // Phi
77 parameterList[index].dimension[3], // Dy1
78 parameterList[index].dimension[4], // Dx1
79 parameterList[index].dimension[5], // Dx2
80 parameterList[index].dimension[6], // pAlp1,
81 parameterList[index].dimension[7], // pDy2,
82 parameterList[index].dimension[8], // pDx3,
83 parameterList[index].dimension[9], // pDx4,
84 parameterList[index].dimension[10]); // pAlp2
85}
86
87void G4GDMLParameterisation::
88ComputeDimensions(G4Tubs& tubs,const G4int index,const G4VPhysicalVolume*) const
89{
90 tubs.SetInnerRadius(parameterList[index].dimension[0]);
91 tubs.SetOuterRadius(parameterList[index].dimension[1]);
92 tubs.SetZHalfLength(parameterList[index].dimension[2]);
93 tubs.SetStartPhiAngle(parameterList[index].dimension[3]);
94 tubs.SetDeltaPhiAngle(parameterList[index].dimension[4]);
95}
96
97void G4GDMLParameterisation::
98ComputeDimensions(G4Cons& cons,const G4int index,const G4VPhysicalVolume*) const
99{
100 cons.SetInnerRadiusMinusZ(parameterList[index].dimension[0]);
101 cons.SetOuterRadiusMinusZ(parameterList[index].dimension[1]);
102 cons.SetInnerRadiusPlusZ(parameterList[index].dimension[2]);
103 cons.SetOuterRadiusPlusZ(parameterList[index].dimension[3]);
104 cons.SetZHalfLength(parameterList[index].dimension[4]);
105 cons.SetStartPhiAngle(parameterList[index].dimension[5]);
106 cons.SetDeltaPhiAngle(parameterList[index].dimension[6]);
107}
108
109void G4GDMLParameterisation::
110ComputeDimensions(G4Sphere& sphere,const G4int index,const G4VPhysicalVolume*) const
111{
112 sphere.SetInsideRadius(parameterList[index].dimension[0]);
113 sphere.SetOuterRadius(parameterList[index].dimension[1]);
114 sphere.SetStartPhiAngle(parameterList[index].dimension[2]);
115 sphere.SetDeltaPhiAngle(parameterList[index].dimension[3]);
116 sphere.SetStartThetaAngle(parameterList[index].dimension[4]);
117 sphere.SetDeltaThetaAngle(parameterList[index].dimension[5]);
118}
119
120void G4GDMLParameterisation::
121ComputeDimensions(G4Orb& orb,const G4int index,const G4VPhysicalVolume*) const
122{
123 orb.SetRadius(parameterList[index].dimension[0]);
124}
125
126void G4GDMLParameterisation::
127ComputeDimensions(G4Torus& torus,const G4int index,const G4VPhysicalVolume*) const
128{
129 torus.SetAllParameters(parameterList[index].dimension[0], // pRmin
130 parameterList[index].dimension[1], // pRmax
131 parameterList[index].dimension[2], // pRtor
132 parameterList[index].dimension[3], // pSPhi
133 parameterList[index].dimension[4]); // pDPhi
134}
135
136void G4GDMLParameterisation::
137ComputeDimensions(G4Para& para,const G4int index,const G4VPhysicalVolume*) const
138{
139 para.SetXHalfLength(parameterList[index].dimension[0]);
140 para.SetYHalfLength(parameterList[index].dimension[1]);
141 para.SetZHalfLength(parameterList[index].dimension[2]);
142 para.SetAlpha(parameterList[index].dimension[3]);
143 para.SetTanAlpha(std::tan(parameterList[index].dimension[3]));
144 para.SetThetaAndPhi(parameterList[index].dimension[4],parameterList[index].dimension[5]);
145}
146
147void G4GDMLParameterisation::
148ComputeDimensions(G4Hype& hype,const G4int index,const G4VPhysicalVolume*) const
149{
150 hype.SetInnerRadius(parameterList[index].dimension[0]);
151 hype.SetOuterRadius(parameterList[index].dimension[1]);
152 hype.SetZHalfLength(parameterList[index].dimension[4]);
153 hype.SetInnerStereo(parameterList[index].dimension[2]);
154 hype.SetOuterStereo(parameterList[index].dimension[3]);
155}
156
157void G4GDMLParameterisation::
158ComputeDimensions(G4Polycone&,const G4int,const G4VPhysicalVolume*) const
159{
160 G4Exception("G4GDMLParameterisation::ComputeDimensions()",
161 "InvalidSetup", FatalException,
162 "Parameterisation of G4Polycone not implemented yet. Sorry!");
163}
164
165void G4GDMLParameterisation::
166ComputeDimensions(G4Polyhedra&,const G4int,const G4VPhysicalVolume*) const
167{
168 G4Exception("G4GDMLParameterisation::ComputeDimensions()",
169 "InvalidSetup", FatalException,
170 "Parameterisation of G4Polyhedra not implemented yet. Sorry!");
171}
@ FatalException
int G4int
Definition: G4Types.hh:66
Definition: G4Box.hh:55
void SetZHalfLength(G4double dz)
Definition: G4Box.cc:170
void SetYHalfLength(G4double dy)
Definition: G4Box.cc:150
void SetXHalfLength(G4double dx)
Definition: G4Box.cc:130
Definition: G4Cons.hh:75
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 AddParameter(const PARAMETER &)
Definition: G4Hype.hh:67
void SetOuterStereo(G4double newOSte)
void SetOuterRadius(G4double newORad)
void SetZHalfLength(G4double newHLZ)
void SetInnerStereo(G4double newISte)
void SetInnerRadius(G4double newIRad)
Definition: G4Orb.hh:52
void SetRadius(G4double newRmax)
Definition: G4Para.hh:77
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 SetDeltaPhiAngle(G4double newDphi)
void SetStartThetaAngle(G4double newSTheta)
void SetOuterRadius(G4double newRmax)
void SetDeltaThetaAngle(G4double newDTheta)
void SetInsideRadius(G4double newRmin)
void SetStartPhiAngle(G4double newSphi, G4bool trig=true)
void SetAllParameters(G4double pRmin, G4double pRmax, G4double pRtor, G4double pSPhi, G4double pDPhi)
Definition: G4Torus.cc:100
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:603
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:77
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 *)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41