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