Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TwistTrapParallelSide.hh
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// G4TwistTrapParallelSide
27//
28// Class description:
29//
30// Class describing a twisted boundary surface for a trapezoid.
31
32// Author: 27-Oct-2004 - O.Link ([email protected])
33// --------------------------------------------------------------------
34#ifndef G4TWISTTRAPPARALLELSIDE_HH
35#define G4TWISTTRAPPARALLELSIDE_HH
36
37#include "G4VTwistSurface.hh"
38
39#include <vector>
40
42{
43 public: // with description
44
46 G4double PhiTwist, // twist angle
47 G4double pDz, // half z lenght
48 G4double pTheta, // direction between end planes
49 G4double pPhi, // by polar and azimutal angles
50 G4double pDy1, // half y length at -pDz
51 G4double pDx1, // half x length at -pDz,-pDy
52 G4double pDx2, // half x length at -pDz,+pDy
53 G4double pDy2, // half y length at +pDz
54 G4double pDx3, // half x length at +pDz,-pDy
55 G4double pDx4, // half x length at +pDz,+pDy
56 G4double pAlph, // tilt angle at +pDz
57 G4double AngleSide // parity
58 );
59
61
62 virtual G4ThreeVector GetNormal(const G4ThreeVector& xx,
63 G4bool isGlobal = false) ;
64
65 virtual G4int DistanceToSurface(const G4ThreeVector& gp,
66 const G4ThreeVector& gv,
67 G4ThreeVector gxx[],
68 G4double distance[],
69 G4int areacode[],
70 G4bool isvalid[],
71 EValidate validate = kValidateWithTol);
72
73 virtual G4int DistanceToSurface(const G4ThreeVector& gp,
74 G4ThreeVector gxx[],
75 G4double distance[],
76 G4int areacode[]);
77
78 public: // without description
79
80 G4TwistTrapParallelSide(__void__&);
81 // Fake default constructor for usage restricted to direct object
82 // persistency for clients requiring preallocation of memory for
83 // persistifiable objects.
84
85 private:
86
87 virtual G4int GetAreaCode(const G4ThreeVector& xx,
88 G4bool withTol = true);
89 virtual void SetCorners();
90 virtual void SetBoundaries();
91
92 void GetPhiUAtX(G4ThreeVector p, G4double& phi, G4double& u);
93 G4ThreeVector ProjectPoint(const G4ThreeVector& p,
94 G4bool isglobal = false);
95
96 virtual G4ThreeVector SurfacePoint(G4double phi, G4double u,
97 G4bool isGlobal = false);
98 virtual G4double GetBoundaryMin(G4double phi);
99 virtual G4double GetBoundaryMax(G4double phi);
100 virtual G4double GetSurfaceArea();
101 virtual void GetFacets( G4int m, G4int n, G4double xyz[][3],
102 G4int faces[][4], G4int iside );
103
104 inline G4ThreeVector NormAng(G4double phi, G4double u);
105 inline G4double GetValueB(G4double phi) ;
106 inline G4double Xcoef(G4double u);
107 // To calculate the w(u) function
108
109 private:
110
111 G4double fTheta;
112 G4double fPhi ;
113
114 G4double fDy1;
115 G4double fDx1;
116 G4double fDx2;
117
118 G4double fDy2;
119 G4double fDx3;
120 G4double fDx4;
121
122 G4double fDz; // Half-length along the z axis
123
124 G4double fAlph;
125 G4double fTAlph; // std::tan(fAlph)
126
127 G4double fPhiTwist; // twist angle ( dphi in surface equation)
128
129 G4double fAngleSide;
130
131 G4double fdeltaX;
132 G4double fdeltaY;
133
134 G4double fDx4plus2; // fDx4 + fDx2 == a2/2 + a1/2
135 G4double fDx4minus2; // fDx4 - fDx2 -
136 G4double fDx3plus1; // fDx3 + fDx1 == d2/2 + d1/2
137 G4double fDx3minus1; // fDx3 - fDx1 -
138 G4double fDy2plus1; // fDy2 + fDy1 == b2/2 + b1/2
139 G4double fDy2minus1; // fDy2 - fDy1 -
140 G4double fa1md1; // 2 fDx2 - 2 fDx1 == a1 - d1
141 G4double fa2md2; // 2 fDx4 - 2 fDx3
142};
143
144//========================================================
145// inline functions
146//========================================================
147
148
149inline
150G4double G4TwistTrapParallelSide::GetValueB(G4double phi)
151{
152 return ( fDy2plus1 + fDy2minus1 * ( 2 * phi ) / fPhiTwist ) ;
153}
154
155inline
156G4double G4TwistTrapParallelSide::Xcoef(G4double phi)
157{
158 return GetValueB(phi)/2. ;
159}
160
161inline G4ThreeVector
162G4TwistTrapParallelSide::
163SurfacePoint( G4double phi, G4double u, G4bool isGlobal )
164{
165 // function to calculate a point on the surface, given by parameters phi,u
166
167 G4ThreeVector SurfPoint ( u*std::cos(phi) - Xcoef(phi)*std::sin(phi)
168 + fdeltaX*phi/fPhiTwist,
169 u*std::sin(phi) + Xcoef(phi)*std::cos(phi)
170 + fdeltaY*phi/fPhiTwist,
171 2*fDz*phi/fPhiTwist );
172 if (isGlobal) { return (fRot * SurfPoint + fTrans); }
173 return SurfPoint;
174}
175
176
177inline
178G4double G4TwistTrapParallelSide::GetBoundaryMin(G4double phi)
179{
180 return -(fPhiTwist*(fDx2 + fDx4 - fDy2plus1*fTAlph)
181 + 2*fDx4minus2*phi - 2*fDy2minus1*fTAlph*phi)/(2.*fPhiTwist) ;
182}
183
184inline
185G4double G4TwistTrapParallelSide::GetBoundaryMax(G4double phi)
186{
187 return (fDx2 + fDx4 + fDy2plus1*fTAlph)/ 2.
188 + ((fDx4minus2 + fDy2minus1*fTAlph)*phi)/fPhiTwist ;
189}
190
191inline
192G4double G4TwistTrapParallelSide::GetSurfaceArea()
193{
194 return 2*fDx4plus2*fDz ;
195}
196
197inline
198G4ThreeVector G4TwistTrapParallelSide::NormAng( G4double phi, G4double u )
199{
200 // function to calculate the norm at a given point on the surface
201 // replace a1-d1
202
203 G4ThreeVector nvec(-2*fDz*std::sin(phi) ,
204 2*fDz*std::cos(phi) ,
205 -(fDy2minus1 + fPhiTwist*u + fdeltaY*std::cos(phi)
206 -fdeltaX*std::sin(phi))) ;
207 return nvec.unit();
208}
209
210#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
virtual G4ThreeVector GetNormal(const G4ThreeVector &xx, G4bool isGlobal=false)
virtual G4int DistanceToSurface(const G4ThreeVector &gp, const G4ThreeVector &gv, G4ThreeVector gxx[], G4double distance[], G4int areacode[], G4bool isvalid[], EValidate validate=kValidateWithTol)
G4RotationMatrix fRot
G4ThreeVector fTrans