Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4GDMLParameterisation.cc
이 파일의 문서화 페이지로 가기
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: G4GDMLParameterisation.cc 77556 2013-11-26 08:56:14Z gcosmo $
27 //
28 // class G4GDMLParameterisation Implementation
29 //
30 // History:
31 // - Created. Zoltan Torzsok, November 2007
32 // -------------------------------------------------------------------------
33 
35 #include "G4PolyconeHistorical.hh"
36 #include "G4PolyhedraHistorical.hh"
37 
39 {
40  return (G4int)parameterList.size();
41 }
42 
44 {
45  parameterList.push_back(newParameter);
46 }
47 
49 ComputeTransformation(const G4int index,G4VPhysicalVolume* physvol) const
50 {
51  physvol->SetTranslation(parameterList[index].position);
52  physvol->SetRotation(parameterList[index].pRot);
53 }
54 
56 ComputeDimensions(G4Box& box,const G4int index,const G4VPhysicalVolume*) const
57 {
58  box.SetXHalfLength(parameterList[index].dimension[0]);
59  box.SetYHalfLength(parameterList[index].dimension[1]);
60  box.SetZHalfLength(parameterList[index].dimension[2]);
61 }
62 
64 ComputeDimensions(G4Trd& trd,const G4int index,const G4VPhysicalVolume*) const
65 {
66  trd.SetXHalfLength1(parameterList[index].dimension[0]);
67  trd.SetXHalfLength2(parameterList[index].dimension[1]);
68  trd.SetYHalfLength1(parameterList[index].dimension[2]);
69  trd.SetYHalfLength2(parameterList[index].dimension[3]);
70  trd.SetZHalfLength(parameterList[index].dimension[4]);
71 }
72 
74 ComputeDimensions(G4Trap& trap,const G4int index,const G4VPhysicalVolume*) const
75 {
76  trap.SetAllParameters(parameterList[index].dimension[0], // Dz
77  parameterList[index].dimension[1], // Theta
78  parameterList[index].dimension[2], // Phi
79  parameterList[index].dimension[3], // Dy1
80  parameterList[index].dimension[4], // Dx1
81  parameterList[index].dimension[5], // Dx2
82  parameterList[index].dimension[6], // pAlp1,
83  parameterList[index].dimension[7], // pDy2,
84  parameterList[index].dimension[8], // pDx3,
85  parameterList[index].dimension[9], // pDx4,
86  parameterList[index].dimension[10]); // pAlp2
87 }
88 
90 ComputeDimensions(G4Tubs& tubs,const G4int index,const G4VPhysicalVolume*) const
91 {
92  tubs.SetInnerRadius(parameterList[index].dimension[0]);
93  tubs.SetOuterRadius(parameterList[index].dimension[1]);
94  tubs.SetZHalfLength(parameterList[index].dimension[2]);
95  tubs.SetStartPhiAngle(parameterList[index].dimension[3]);
96  tubs.SetDeltaPhiAngle(parameterList[index].dimension[4]);
97 }
98 
100 ComputeDimensions(G4Cons& cons,const G4int index,const G4VPhysicalVolume*) const
101 {
102  cons.SetInnerRadiusMinusZ(parameterList[index].dimension[0]);
103  cons.SetOuterRadiusMinusZ(parameterList[index].dimension[1]);
104  cons.SetInnerRadiusPlusZ(parameterList[index].dimension[2]);
105  cons.SetOuterRadiusPlusZ(parameterList[index].dimension[3]);
106  cons.SetZHalfLength(parameterList[index].dimension[4]);
107  cons.SetStartPhiAngle(parameterList[index].dimension[5]);
108  cons.SetDeltaPhiAngle(parameterList[index].dimension[6]);
109 }
110 
112 ComputeDimensions(G4Sphere& sphere,const G4int index,const G4VPhysicalVolume*) const
113 {
114  sphere.SetInnerRadius(parameterList[index].dimension[0]);
115  sphere.SetOuterRadius(parameterList[index].dimension[1]);
116  sphere.SetStartPhiAngle(parameterList[index].dimension[2]);
117  sphere.SetDeltaPhiAngle(parameterList[index].dimension[3]);
118  sphere.SetStartThetaAngle(parameterList[index].dimension[4]);
119  sphere.SetDeltaThetaAngle(parameterList[index].dimension[5]);
120 }
121 
123 ComputeDimensions(G4Orb& orb,const G4int index,const G4VPhysicalVolume*) const
124 {
125  orb.SetRadius(parameterList[index].dimension[0]);
126 }
127 
129 ComputeDimensions(G4Ellipsoid& ellipsoid,const G4int index,const G4VPhysicalVolume*) const
130 {
131  ellipsoid.SetSemiAxis(parameterList[index].dimension[0],
132  parameterList[index].dimension[1],
133  parameterList[index].dimension[2]);
134  ellipsoid.SetZCuts(parameterList[index].dimension[3],
135  parameterList[index].dimension[4]);
136 }
137 
139 ComputeDimensions(G4Torus& torus,const G4int index,const G4VPhysicalVolume*) const
140 {
141  torus.SetAllParameters(parameterList[index].dimension[0], // pRmin
142  parameterList[index].dimension[1], // pRmax
143  parameterList[index].dimension[2], // pRtor
144  parameterList[index].dimension[3], // pSPhi
145  parameterList[index].dimension[4]); // pDPhi
146 }
147 
149 ComputeDimensions(G4Para& para,const G4int index,const G4VPhysicalVolume*) const
150 {
151  para.SetXHalfLength(parameterList[index].dimension[0]);
152  para.SetYHalfLength(parameterList[index].dimension[1]);
153  para.SetZHalfLength(parameterList[index].dimension[2]);
154  para.SetAlpha(parameterList[index].dimension[3]);
155  para.SetTanAlpha(std::tan(parameterList[index].dimension[3]));
156  para.SetThetaAndPhi(parameterList[index].dimension[4],parameterList[index].dimension[5]);
157 }
158 
160 ComputeDimensions(G4Hype& hype,const G4int index,const G4VPhysicalVolume*) const
161 {
162  hype.SetInnerRadius(parameterList[index].dimension[0]);
163  hype.SetOuterRadius(parameterList[index].dimension[1]);
164  hype.SetZHalfLength(parameterList[index].dimension[4]);
165  hype.SetInnerStereo(parameterList[index].dimension[2]);
166  hype.SetOuterStereo(parameterList[index].dimension[3]);
167 }
168 
170 ComputeDimensions(G4Polycone& pcone,const G4int index,const G4VPhysicalVolume*) const
171 {
172 
173  G4PolyconeHistorical origparam( *(pcone.GetOriginalParameters()) );
174  origparam.Start_angle = parameterList[index].dimension[0];
175  origparam.Opening_angle = parameterList[index].dimension[1];
176  origparam.Num_z_planes = (G4int) parameterList[index].dimension[2];
177  G4int nZplanes = origparam.Num_z_planes;
178 
179  for( G4int ii = 0; ii < nZplanes; ii++ )
180  {
181  origparam.Rmin[ii] = parameterList[index].dimension[3+ii*3] ;
182  origparam.Rmax[ii] = parameterList[index].dimension[4+ii*3] ;
183  origparam.Z_values[ii] = parameterList[index].dimension[5+ii*3] ;
184  }
185 
186  pcone.SetOriginalParameters(&origparam); // copy values & transfer pointers
187  pcone.Reset(); // reset to new solid parameters
188 
189 }
190 
192 ComputeDimensions(G4Polyhedra& polyhedra,const G4int index,const G4VPhysicalVolume*) const
193 {
194  G4PolyhedraHistorical origparam( *(polyhedra.GetOriginalParameters()) );
195  origparam.Start_angle = parameterList[index].dimension[0];
196  origparam.Opening_angle = parameterList[index].dimension[1];
197  origparam.Num_z_planes = (G4int) parameterList[index].dimension[2];
198  origparam.numSide = (G4int) parameterList[index].dimension[3];
199 
200  G4int nZplanes = origparam.Num_z_planes;
201 
202  for( G4int ii = 0; ii < nZplanes; ii++ )
203  {
204 
205  origparam.Rmin[ii] = parameterList[index].dimension[4+ii*3] ;
206  origparam.Rmax[ii] = parameterList[index].dimension[5+ii*3] ;
207  origparam.Z_values[ii] = parameterList[index].dimension[6+ii*3] ;
208 
209  }
210  polyhedra.SetOriginalParameters(&origparam); // copy values & transfer pointers
211  polyhedra.Reset(); // reset to new solid parameters
212 }
void SetXHalfLength2(G4double val)
void SetAllParameters(G4double pRmin, G4double pRmax, G4double pRtor, G4double pSPhi, G4double pDPhi)
Definition: G4Torus.cc:103
void SetXHalfLength(G4double dx)
Definition: G4Box.cc:136
G4bool Reset()
Definition: G4Polyhedra.cc:486
void SetYHalfLength1(G4double val)
void SetInnerRadius(G4double newRMin)
void SetAlpha(G4double alpha)
void SetOuterRadius(G4double newRMax)
void SetOriginalParameters(G4PolyconeHistorical *pars)
void SetDeltaThetaAngle(G4double newDTheta)
void SetYHalfLength(G4double val)
void SetYHalfLength2(G4double val)
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:298
Definition: G4Tubs.hh:85
void SetDeltaPhiAngle(G4double newDphi)
void SetZHalfLength(G4double dz)
Definition: G4Box.cc:183
void SetOuterStereo(G4double newOSte)
void AddParameter(const PARAMETER &)
void SetRadius(G4double newRmax)
void SetXHalfLength1(G4double val)
void SetTranslation(const G4ThreeVector &v)
void SetOuterRadiusPlusZ(G4double Rmax2)
void SetStartPhiAngle(G4double newSphi, G4bool trig=true)
void SetXHalfLength(G4double val)
G4PolyhedraHistorical * GetOriginalParameters() const
void SetTanAlpha(G4double val)
void SetOuterRadiusMinusZ(G4double Rmax1)
void SetStartThetaAngle(G4double newSTheta)
void SetOriginalParameters(G4PolyhedraHistorical *pars)
void SetStartPhiAngle(G4double newSPhi, G4bool trig=true)
void SetSemiAxis(G4double x, G4double y, G4double z)
void SetZHalfLength(G4double newDz)
void SetZHalfLength(G4double val)
void SetStartPhiAngle(G4double newSPhi, G4bool trig=true)
Definition: G4Cons.hh:83
void SetInnerRadius(G4double newIRad)
void SetZHalfLength(G4double newDz)
Definition: G4Box.hh:64
void ComputeTransformation(const G4int, G4VPhysicalVolume *) const
void SetZHalfLength(G4double val)
Definition: G4Orb.hh:62
std::vector< PARAMETER > parameterList
void SetZCuts(G4double newzBottomCut, G4double newzTopCut)
void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
G4bool Reset()
Definition: G4Polycone.cc:447
void SetInnerRadius(G4double newRMin)
void SetYHalfLength(G4double dy)
Definition: G4Box.cc:160
void SetOuterRadius(G4double newORad)
void SetDeltaPhiAngle(G4double newDPhi)
void SetDeltaPhiAngle(G4double newDPhi)
int G4int
Definition: G4Types.hh:78
void SetInnerRadiusPlusZ(G4double Rmin2)
Definition: G4Trd.hh:72
void SetInnerRadiusMinusZ(G4double Rmin1)
void SetRotation(G4RotationMatrix *)
G4PolyconeHistorical * GetOriginalParameters() const
void SetOuterRadius(G4double newRmax)
void SetInnerStereo(G4double newISte)
void SetZHalfLength(G4double newHLZ)
Definition: G4Para.hh:86
void SetThetaAndPhi(double pTheta, double pPhi)
Definition: G4Hype.hh:76