Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4Hype.hh
이 파일의 문서화 페이지로 가기
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: G4Hype.hh 108078 2017-12-20 08:15:44Z gcosmo $
28 // $Original: G4Hype.hh,v 1.0 1998/06/09 16:57:50 safai Exp $
29 //
30 //
31 // --------------------------------------------------------------------
32 // GEANT 4 class header file
33 //
34 //
35 // G4Hype
36 //
37 // Class description:
38 //
39 // This class implements a tube with hyperbolic profile.
40 //
41 // It describes an hyperbolic volume with curved sides parallel to
42 // the z-axis. The solid has a specified half-length along the z axis,
43 // about which it is centered, and a given minimum and maximum radius.
44 // A minimum radius of 0 signifies a filled Hype (with hyperbolical
45 // inner surface). To have a filled Hype the user must specify
46 // inner radius = 0 AND inner stereo angle = 0.
47 //
48 // The inner and outer hyperbolical surfaces can have different
49 // stereo angles. A stereo angle of 0 gives a cylindrical surface.
50 
51 // Authors:
52 // Ernesto Lamanna (Ernesto.Lamanna@roma1.infn.it) &
53 // Francesco Safai Tehrani (Francesco.SafaiTehrani@roma1.infn.it)
54 // Rome, INFN & University of Rome "La Sapienza", 9 June 1998.
55 //
56 // --------------------------------------------------------------------
57 #ifndef G4HYPE_HH
58 #define G4HYPE_HH
59 
60 #if defined(G4GEOM_USE_USOLIDS)
61 #define G4GEOM_USE_UHYPE 1
62 #endif
63 
64 #if (defined(G4GEOM_USE_UHYPE) && defined(G4GEOM_USE_SYS_USOLIDS))
65  #define G4UHype G4Hype
66  #include "G4UHype.hh"
67 #else
68 
69 #include "G4VSolid.hh"
70 #include "G4ThreeVector.hh"
71 #include "G4Polyhedron.hh"
72 
73 class G4SolidExtentList;
74 class G4ClippablePolygon;
75 
76 class G4Hype : public G4VSolid
77 {
78  public: // with description
79 
80  G4Hype(const G4String& pName,
81  G4double newInnerRadius,
82  G4double newOuterRadius,
83  G4double newInnerStereo,
84  G4double newOuterStereo,
85  G4double newHalfLenZ);
86 
87  virtual ~G4Hype();
88 
90  const G4int n,
91  const G4VPhysicalVolume* pRep);
92 
93  void BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const;
94 
95  G4bool CalculateExtent(const EAxis pAxis,
96  const G4VoxelLimits& pVoxelLimit,
97  const G4AffineTransform& pTransform,
98  G4double& pMin, G4double& pMax) const;
99 
100  inline G4double GetInnerRadius () const;
101  inline G4double GetOuterRadius () const;
102  inline G4double GetZHalfLength () const;
103  inline G4double GetInnerStereo () const;
104  inline G4double GetOuterStereo () const;
105 
106  inline void SetInnerRadius (G4double newIRad);
107  inline void SetOuterRadius (G4double newORad);
108  inline void SetZHalfLength (G4double newHLZ);
109  inline void SetInnerStereo (G4double newISte);
110  inline void SetOuterStereo (G4double newOSte);
111 
112  EInside Inside(const G4ThreeVector& p) const;
113 
114  G4ThreeVector SurfaceNormal(const G4ThreeVector& p) const;
115 
116  G4double DistanceToIn(const G4ThreeVector& p, const G4ThreeVector& v) const;
117  G4double DistanceToIn(const G4ThreeVector& p) const;
119  const G4bool calcNorm=G4bool(false),
120  G4bool *validNorm=0, G4ThreeVector *n=0) const;
121  G4double DistanceToOut(const G4ThreeVector& p) const;
122 
124 
125  G4VSolid* Clone() const;
126 
127  std::ostream& StreamInfo(std::ostream& os) const;
128 
131 
133 
134  void DescribeYourselfTo (G4VGraphicsScene& scene) const;
135  G4VisExtent GetExtent () const;
136  G4Polyhedron* CreatePolyhedron () const;
137  G4Polyhedron* GetPolyhedron () const;
138 
139  public: // without description
140 
141  G4Hype(__void__&);
142  // Fake default constructor for usage restricted to direct object
143  // persistency for clients requiring preallocation of memory for
144  // persistifiable objects.
145 
146  G4Hype(const G4Hype& rhs);
147  G4Hype& operator=(const G4Hype& rhs);
148  // Copy constructor and assignment operator.
149 
150  protected: // without description
151 
152  inline G4bool InnerSurfaceExists() const;
153  // whether we have an inner surface or not
154 
156  G4double r0, G4double tanPhi );
158  G4double r0, G4double tan2Phi );
159  // approximate isotropic distance to hyperbolic surface
160 
161  inline G4double HypeInnerRadius2(G4double zVal) const;
162  inline G4double HypeOuterRadius2(G4double zVal) const;
163  // values of hype radius at a given Z
164 
165  static G4int IntersectHype( const G4ThreeVector &p, const G4ThreeVector &v,
166  G4double r2, G4double tan2Phi, G4double s[2] );
167  // intersection with hyperbolic surface
168 
169  protected:
170 
176 
177  // precalculated parameters, squared quantities
178 
181  G4double tanInnerStereo2; // squared tan of Inner Stereo angle
182  G4double tanOuterStereo2; // squared tan of Outer Stereo angle
183  G4double innerRadius2; // squared Inner Radius
184  G4double outerRadius2; // squared Outer Radius
185  G4double endInnerRadius2; // squared endcap Inner Radius
186  G4double endOuterRadius2; // squared endcap Outer Radius
187  G4double endInnerRadius; // endcap Inner Radius
188  G4double endOuterRadius; // endcap Outer Radius
189 
190  // Used by distanceToOut
191 
193 
194  private:
195 
196  G4double asinh(G4double arg);
197 
198  private:
199 
202 
204 
207 };
208 
209 #include "G4Hype.icc"
210 
211 #endif // defined(G4GEOM_USE_UHYPE) && defined(G4GEOM_USE_SYS_USOLIDS)
212 
213 #endif // G4HYPE_HH
G4double HypeInnerRadius2(G4double zVal) const
G4double GetOuterStereo() const
G4Polyhedron * CreatePolyhedron() const
Definition: G4Hype.cc:1309
G4double tanInnerStereo2
Definition: G4Hype.hh:181
G4double fHalfTol
Definition: G4Hype.hh:203
G4double GetOuterRadius() const
G4double tanInnerStereo
Definition: G4Hype.hh:179
G4double outerRadius
Definition: G4Hype.hh:172
static G4double ApproxDistOutside(G4double pr, G4double pz, G4double r0, G4double tanPhi)
Definition: G4Hype.cc:1036
G4double innerRadius2
Definition: G4Hype.hh:183
const char * p
Definition: xmltok.h:285
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Hype.cc:312
void SetOuterStereo(G4double newOSte)
G4VisExtent GetExtent() const
Definition: G4Hype.cc:1296
G4double fCubicVolume
Definition: G4Hype.hh:200
G4double HypeOuterRadius2(G4double zVal) const
G4double endOuterRadius
Definition: G4Hype.hh:188
G4double GetSurfaceArea()
Definition: G4Hype.cc:1147
G4double asinh(G4double arg)
Definition: G4Hype.cc:1339
G4VSolid * Clone() const
Definition: G4Hype.cc:1127
G4Hype(const G4String &pName, G4double newInnerRadius, G4double newOuterRadius, G4double newInnerStereo, G4double newOuterStereo, G4double newHalfLenZ)
Definition: G4Hype.cc:76
G4double GetCubicVolume()
Definition: G4Hype.cc:1136
static G4double ApproxDistInside(G4double pr, G4double pz, G4double r0, G4double tan2Phi)
Definition: G4Hype.cc:1094
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Hype.cc:1287
G4Hype & operator=(const G4Hype &rhs)
Definition: G4Hype.cc:181
const XML_Char * s
Definition: expat.h:262
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
Definition: G4Hype.cc:249
G4double halfLenZ
Definition: G4Hype.hh:173
void SetInnerRadius(G4double newIRad)
G4GeometryType GetEntityType() const
Definition: G4Hype.cc:1118
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4Hype.cc:225
G4double innerRadius
Definition: G4Hype.hh:171
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Hype.cc:358
G4double endOuterRadius2
Definition: G4Hype.hh:186
EInside Inside(const G4ThreeVector &p) const
Definition: G4Hype.cc:268
G4double fSurfaceArea
Definition: G4Hype.hh:201
G4double outerRadius2
Definition: G4Hype.hh:184
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=G4bool(false), G4bool *validNorm=0, G4ThreeVector *n=0) const
Definition: G4Hype.cc:734
static G4int IntersectHype(const G4ThreeVector &p, const G4ThreeVector &v, G4double r2, G4double tan2Phi, G4double s[2])
Definition: G4Hype.cc:968
void SetOuterRadius(G4double newORad)
G4Polyhedron * GetPolyhedron() const
Definition: G4Hype.cc:1319
virtual ~G4Hype()
Definition: G4Hype.cc:154
int G4int
Definition: G4Types.hh:78
EInside
Definition: geomdefs.hh:58
G4bool fRebuildPolyhedron
Definition: G4Hype.hh:205
EAxis
Definition: geomdefs.hh:54
G4bool InnerSurfaceExists() const
G4double endInnerRadius
Definition: G4Hype.hh:187
G4double tanOuterStereo
Definition: G4Hype.hh:180
G4double GetInnerRadius() const
G4double GetInnerStereo() const
G4double endInnerRadius2
Definition: G4Hype.hh:185
G4double outerStereo
Definition: G4Hype.hh:175
Char_t n[5]
G4double tanOuterStereo2
Definition: G4Hype.hh:182
G4ThreeVector GetPointOnSurface() const
Definition: G4Hype.cc:1182
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Hype.cc:214
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Hype.cc:1158
void SetInnerStereo(G4double newISte)
void SetZHalfLength(G4double newHLZ)
G4double GetZHalfLength() const
G4Polyhedron * fpPolyhedron
Definition: G4Hype.hh:206
G4double innerStereo
Definition: G4Hype.hh:174
Definition: G4Hype.hh:76