Geant4  10.01
G4Hype.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: G4Hype.hh 83572 2014-09-01 15:23:27Z 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 #include "G4VSolid.hh"
61 #include "G4ThreeVector.hh"
62 #include "G4Polyhedron.hh"
63 
64 class G4SolidExtentList;
65 class G4ClippablePolygon;
66 
67 class G4Hype : public G4VSolid
68 {
69  public: // with description
70 
71  G4Hype(const G4String& pName,
72  G4double newInnerRadius,
73  G4double newOuterRadius,
74  G4double newInnerStereo,
75  G4double newOuterStereo,
76  G4double newHalfLenZ);
77 
78  virtual ~G4Hype();
79 
81  const G4int n,
82  const G4VPhysicalVolume* pRep);
83 
84  G4bool CalculateExtent(const EAxis pAxis,
85  const G4VoxelLimits& pVoxelLimit,
86  const G4AffineTransform& pTransform,
87  G4double& pmin, G4double& pmax) const;
88 
89  inline G4double GetInnerRadius () const;
90  inline G4double GetOuterRadius () const;
91  inline G4double GetZHalfLength () const;
92  inline G4double GetInnerStereo () const;
93  inline G4double GetOuterStereo () const;
94 
95  inline void SetInnerRadius (G4double newIRad);
96  inline void SetOuterRadius (G4double newORad);
97  inline void SetZHalfLength (G4double newHLZ);
98  inline void SetInnerStereo (G4double newISte);
99  inline void SetOuterStereo (G4double newOSte);
100 
101  EInside Inside(const G4ThreeVector& p) const;
102 
103  G4ThreeVector SurfaceNormal(const G4ThreeVector& p) const;
104 
105  G4double DistanceToIn(const G4ThreeVector& p, const G4ThreeVector& v) const;
106  G4double DistanceToIn(const G4ThreeVector& p) const;
108  const G4bool calcNorm=G4bool(false),
109  G4bool *validNorm=0, G4ThreeVector *n=0) const;
110  G4double DistanceToOut(const G4ThreeVector& p) const;
111 
113 
114  G4VSolid* Clone() const;
115 
116  std::ostream& StreamInfo(std::ostream& os) const;
117 
120 
122 
123  void DescribeYourselfTo (G4VGraphicsScene& scene) const;
124  G4VisExtent GetExtent () const;
125  G4Polyhedron* CreatePolyhedron () const;
126  G4Polyhedron* GetPolyhedron () const;
127 
128  public: // without description
129 
130  G4Hype(__void__&);
131  // Fake default constructor for usage restricted to direct object
132  // persistency for clients requiring preallocation of memory for
133  // persistifiable objects.
134 
135  G4Hype(const G4Hype& rhs);
136  G4Hype& operator=(const G4Hype& rhs);
137  // Copy constructor and assignment operator.
138 
139  protected: // without description
140 
141  inline G4bool InnerSurfaceExists() const;
142  // whether we have an inner surface or not
143 
145  G4double r0, G4double tanPhi );
147  G4double r0, G4double tan2Phi );
148  // approximate isotropic distance to hyperbolic surface
149 
150  inline G4double HypeInnerRadius2(G4double zVal) const;
151  inline G4double HypeOuterRadius2(G4double zVal) const;
152  // values of hype radius at a given Z
153 
154  static G4int IntersectHype( const G4ThreeVector &p, const G4ThreeVector &v,
155  G4double r2, G4double tan2Phi, G4double s[2] );
156  // intersection with hyperbolic surface
157 
158  static void AddPolyToExtent( const G4ThreeVector &v0,
159  const G4ThreeVector &v1,
160  const G4ThreeVector &w1,
161  const G4ThreeVector &w0,
162  const G4VoxelLimits &voxelLimit,
163  const EAxis axis,
164  G4SolidExtentList &extentList );
165 
166  protected:
167 
173 
174  // precalculated parameters, squared quantities
175 
178  G4double tanInnerStereo2; // squared tan of Inner Stereo angle
179  G4double tanOuterStereo2; // squared tan of Outer Stereo angle
180  G4double innerRadius2; // squared Inner Radius
181  G4double outerRadius2; // squared Outer Radius
182  G4double endInnerRadius2; // squared endcap Inner Radius
183  G4double endOuterRadius2; // squared endcap Outer Radius
184  G4double endInnerRadius; // endcap Inner Radius
185  G4double endOuterRadius; // endcap Outer Radius
186 
187  // Used by distanceToOut
188 
190 
191  private:
192 
193  G4double asinh(G4double arg);
194 
195  private:
196 
199 
201 
204 };
205 
206 #include "G4Hype.icc"
207 
208 #endif
G4double outerStereo
Definition: G4Hype.hh:172
G4bool fRebuildPolyhedron
Definition: G4Hype.hh:202
G4double GetCubicVolume()
Definition: G4Hype.cc:1380
G4double outerRadius2
Definition: G4Hype.hh:181
G4double fHalfTol
Definition: G4Hype.hh:200
G4double GetOuterStereo() const
G4Hype(const G4String &pName, G4double newInnerRadius, G4double newOuterRadius, G4double newInnerStereo, G4double newOuterStereo, G4double newHalfLenZ)
Definition: G4Hype.cc:74
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Hype.cc:603
G4double GetInnerStereo() const
CLHEP::Hep3Vector G4ThreeVector
G4double GetInnerRadius() const
G4double GetSurfaceArea()
Definition: G4Hype.cc:1391
G4double innerRadius2
Definition: G4Hype.hh:180
G4GeometryType GetEntityType() const
Definition: G4Hype.cc:1362
G4double innerStereo
Definition: G4Hype.hh:171
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Hype.cc:212
static G4double ApproxDistInside(G4double pr, G4double pz, G4double r0, G4double tan2Phi)
Definition: G4Hype.cc:1338
G4double fSurfaceArea
Definition: G4Hype.hh:198
G4Polyhedron * CreatePolyhedron() const
Definition: G4Hype.cc:1551
G4double endOuterRadius2
Definition: G4Hype.hh:183
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Hype.cc:1529
G4double tanOuterStereo
Definition: G4Hype.hh:177
int G4int
Definition: G4Types.hh:78
G4bool InnerSurfaceExists() const
virtual ~G4Hype()
Definition: G4Hype.cc:152
G4double HypeOuterRadius2(G4double zVal) const
static const double s
Definition: G4SIunits.hh:150
G4double tanOuterStereo2
Definition: G4Hype.hh:179
G4VSolid * Clone() const
Definition: G4Hype.cc:1371
static G4double ApproxDistOutside(G4double pr, G4double pz, G4double r0, G4double tanPhi)
Definition: G4Hype.cc:1280
G4double tanInnerStereo2
Definition: G4Hype.hh:178
void SetOuterRadius(G4double newORad)
G4double outerRadius
Definition: G4Hype.hh:169
G4double HypeInnerRadius2(G4double zVal) const
Definition: G4Hype.hh:67
G4double endOuterRadius
Definition: G4Hype.hh:185
bool G4bool
Definition: G4Types.hh:79
void SetInnerRadius(G4double newIRad)
G4double GetOuterRadius() const
const G4int n
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=G4bool(false), G4bool *validNorm=0, G4ThreeVector *n=0) const
Definition: G4Hype.cc:979
G4Hype & operator=(const G4Hype &rhs)
Definition: G4Hype.cc:179
G4double halfLenZ
Definition: G4Hype.hh:170
G4ThreeVector GetPointOnSurface() const
Definition: G4Hype.cc:1426
static void AddPolyToExtent(const G4ThreeVector &v0, const G4ThreeVector &v1, const G4ThreeVector &w1, const G4ThreeVector &w0, const G4VoxelLimits &voxelLimit, const EAxis axis, G4SolidExtentList &extentList)
Definition: G4Hype.cc:487
G4VisExtent GetExtent() const
Definition: G4Hype.cc:1538
static G4int IntersectHype(const G4ThreeVector &p, const G4ThreeVector &v, G4double r2, G4double tan2Phi, G4double s[2])
Definition: G4Hype.cc:1213
EInside
Definition: geomdefs.hh:58
G4double fCubicVolume
Definition: G4Hype.hh:197
EAxis
Definition: geomdefs.hh:54
G4double tanInnerStereo
Definition: G4Hype.hh:176
EInside Inside(const G4ThreeVector &p) const
Definition: G4Hype.cc:513
void SetZHalfLength(G4double newHLZ)
G4double GetZHalfLength() const
G4Polyhedron * fpPolyhedron
Definition: G4Hype.hh:203
G4double asinh(G4double arg)
Definition: G4Hype.cc:1581
G4double endInnerRadius2
Definition: G4Hype.hh:182
void SetOuterStereo(G4double newOSte)
double G4double
Definition: G4Types.hh:76
G4double innerRadius
Definition: G4Hype.hh:168
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4Hype.cc:223
G4double endInnerRadius
Definition: G4Hype.hh:184
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Hype.cc:1402
G4Polyhedron * GetPolyhedron() const
Definition: G4Hype.cc:1561
const G4double r0
void SetInnerStereo(G4double newISte)
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Hype.cc:557