Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4TwistTubsHypeSide.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// G4TwistTubsHypeSide
27//
28// Class description:
29//
30// Class describing a hyperbolic boundary surface for a cylinder.
31
32// 01-Aug-2002 - Kotoyo Hoshina (hoshina@hepburn.s.chiba-u.ac.jp), created.
33// 13-Nov-2003 - O.Link (Oliver.Link@cern.ch), Integration in Geant4
34// from original version in Jupiter-2.5.02 application.
35// --------------------------------------------------------------------
36#ifndef G4TWISTTUBSHYPESIDE_HH
37#define G4TWISTTUBSHYPESIDE_HH
38
39#include "G4VTwistSurface.hh"
40#include "G4Integrator.hh"
42
44{
45 public: // with description
46
47 G4TwistTubsHypeSide(const G4String& name,
48 const G4RotationMatrix& rot, // 0.5*(phi-width segment)
49 const G4ThreeVector& tlate,
50 const G4int handedness,// R-hand = 1, L-hand = -1
51 const G4double kappa, // tan(TwistAngle/2)/fZHalfLen
52 const G4double tanstereo, // tan(stereo angle)
53 const G4double r0, // radius at z = 0
54 const EAxis axis0 = kPhi,
55 const EAxis axis1 = kZAxis,
56 G4double axis0min = -kInfinity,
57 G4double axis1min = -kInfinity,
58 G4double axis0max = kInfinity,
59 G4double axis1max = kInfinity);
60
61 G4TwistTubsHypeSide(const G4String& name,
62 G4double EndInnerRadius[2],
63 G4double EndOuterRadius[2],
64 G4double DPhi,
65 G4double EndPhi[2],
66 G4double EndZ[2],
67 G4double InnerRadius,
68 G4double OuterRadius,
69 G4double Kappa,
70 G4double TanInnerStereo,
71 G4double TanOuterStereo,
72 G4int handedness) ;
73
74 virtual ~G4TwistTubsHypeSide();
75
76 virtual G4int DistanceToSurface(const G4ThreeVector& gp,
77 const G4ThreeVector& gv,
78 G4ThreeVector gxx[],
79 G4double distance[],
80 G4int areacode[],
81 G4bool isvalid[],
82 EValidate validate = kValidateWithTol);
83
84 virtual G4int DistanceToSurface(const G4ThreeVector& gp,
85 G4ThreeVector gxx[],
86 G4double distance[],
87 G4int areacode[]);
88
89 virtual G4ThreeVector GetNormal(const G4ThreeVector& xx,
90 G4bool isGlobal = false) ;
91 virtual EInside Inside(const G4ThreeVector& gp) ;
92
93 virtual G4double GetRhoAtPZ(const G4ThreeVector& p,
94 G4bool isglobal = false) const ;
95
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 public: // without description
105
106 G4TwistTubsHypeSide(__void__&);
107 // Fake default constructor for usage restricted to direct object
108 // persistency for clients requiring preallocation of memory for
109 // persistifiable objects.
110
111 private:
112
113 virtual G4int GetAreaCode(const G4ThreeVector& xx,
114 G4bool withTol = true);
115 virtual G4int GetAreaCodeInPhi(const G4ThreeVector& xx,
116 G4bool withTol = true);
117 virtual void SetCorners();
118
119 virtual void SetCorners(G4double EndInnerRadius[2],
120 G4double EndOuterRadius[2],
121 G4double DPhi,
122 G4double EndPhi[2],
123 G4double EndZ[2]);
124 virtual void SetBoundaries();
125
126 private:
127
128 G4double fKappa; // std::tan(TwistedAngle/2)/HalfLenZ;
129 G4double fTanStereo; // std::tan(StereoAngle)
130 G4double fTan2Stereo; // std::tan(StereoAngle)**2
131 G4double fR0; // radius at z = 0
132 G4double fR02; // radius**2 at z = 0
133 G4double fDPhi ; // segment
134
135 class Insidetype
136 {
137 public:
138 G4ThreeVector gp;
139 EInside inside;
140 };
141 Insidetype fInside;
142};
143
144//========================================================
145// inline functions
146//========================================================
147
148inline
150 G4bool isglobal) const
151{
152 // Get Rho at p.z() on Hyperbolic Surface.
153 G4ThreeVector tmpp;
154 if (isglobal) { tmpp = fRot.inverse()*p - fTrans; }
155 else { tmpp = p; }
156
157 return std::sqrt(fR02 + tmpp.z() * tmpp.z() * fTan2Stereo);
158}
159
160inline
162SurfacePoint(G4double phi , G4double z , G4bool isGlobal)
163{
164 G4double rho = std::sqrt(fR02 + z * z * fTan2Stereo) ;
165
166 G4ThreeVector SurfPoint (rho*std::cos(phi), rho*std::sin(phi), z) ;
167
168 if (isGlobal) { return (fRot * SurfPoint + fTrans); }
169 return SurfPoint;
170}
171
172inline
174{
175 G4ThreeVector ptmp(0,0,z) ; // temporary point with z Komponent only
176 G4ThreeVector lowerlimit; // lower phi-boundary limit at z = ptmp.z()
177 lowerlimit = GetBoundaryAtPZ(sAxis0 & sAxisMin, ptmp);
178 return std::atan2( lowerlimit.y(), lowerlimit.x() ) ;
179}
180
181inline
183{
184 G4ThreeVector ptmp(0,0,z) ; // temporary point with z Komponent only
185 G4ThreeVector upperlimit; // upper phi-boundary limit at z = ptmp.z()
186 upperlimit = GetBoundaryAtPZ(sAxis0 & sAxisMax, ptmp);
187 return std::atan2( upperlimit.y(), upperlimit.x() ) ;
188}
189
190inline
192{
193 // approximation with tube surface
194
195 return ( fAxisMax[1] - fAxisMin[1] ) * fR0 * fDPhi ;
196}
197
198#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
double z() const
double x() const
double y() const
HepRotation inverse() const
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)
virtual G4double GetBoundaryMin(G4double phi)
virtual EInside Inside(const G4ThreeVector &gp)
virtual G4ThreeVector SurfacePoint(G4double, G4double, G4bool isGlobal=false)
virtual void GetFacets(G4int m, G4int n, G4double xyz[][3], G4int faces[][4], G4int iside)
virtual G4double GetRhoAtPZ(const G4ThreeVector &p, G4bool isglobal=false) const
virtual G4double GetSurfaceArea()
virtual G4double GetBoundaryMax(G4double phi)
static const G4int sAxisMax
static const G4int sAxis0
G4double fAxisMax[2]
G4RotationMatrix fRot
static const G4int sAxisMin
virtual G4ThreeVector GetBoundaryAtPZ(G4int areacode, const G4ThreeVector &p) const
G4ThreeVector fTrans
G4double fAxisMin[2]
EAxis
Definition: geomdefs.hh:54
@ kPhi
Definition: geomdefs.hh:60
@ kZAxis
Definition: geomdefs.hh:57
EInside
Definition: geomdefs.hh:67