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