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
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// G4Hype
27//
28// Class description:
29//
30// This class implements a tube with hyperbolic profile.
31//
32// It describes an hyperbolic volume with curved sides parallel to
33// the z-axis. The solid has a specified half-length along the z axis,
34// about which it is centered, and a given minimum and maximum radius.
35// A minimum radius of 0 signifies a filled Hype (with hyperbolical
36// inner surface). To have a filled Hype the user must specify
37// inner radius = 0 AND inner stereo angle = 0.
38//
39// The inner and outer hyperbolical surfaces can have different
40// stereo angles. A stereo angle of 0 gives a cylindrical surface.
41
42// Authors:
43// Ernesto Lamanna (Ernesto.Lamanna@roma1.infn.it) &
44// Francesco Safai Tehrani (Francesco.SafaiTehrani@roma1.infn.it)
45// Rome, INFN & University of Rome "La Sapienza", 9 June 1998.
46// --------------------------------------------------------------------
47#ifndef G4HYPE_HH
48#define G4HYPE_HH
49
50#include "G4GeomTypes.hh"
51
52#if defined(G4GEOM_USE_USOLIDS)
53#define G4GEOM_USE_UHYPE 1
54#endif
55
56#if (defined(G4GEOM_USE_UHYPE) && defined(G4GEOM_USE_SYS_USOLIDS))
57 #define G4UHype G4Hype
58 #include "G4UHype.hh"
59#else
60
61#include "G4VSolid.hh"
62#include "G4ThreeVector.hh"
63#include "G4Polyhedron.hh"
64
67
68class G4Hype : public G4VSolid
69{
70 public: // with description
71
72 G4Hype(const G4String& pName,
73 G4double newInnerRadius,
74 G4double newOuterRadius,
75 G4double newInnerStereo,
76 G4double newOuterStereo,
77 G4double newHalfLenZ);
78
79 virtual ~G4Hype();
80
82 const G4int n,
83 const G4VPhysicalVolume* pRep);
84
85 void BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const;
86
87 G4bool CalculateExtent(const EAxis pAxis,
88 const G4VoxelLimits& pVoxelLimit,
89 const G4AffineTransform& pTransform,
90 G4double& pMin, G4double& pMax) const;
91
92 inline G4double GetInnerRadius () const;
93 inline G4double GetOuterRadius () const;
94 inline G4double GetZHalfLength () const;
95 inline G4double GetInnerStereo () const;
96 inline G4double GetOuterStereo () const;
97
98 inline void SetInnerRadius (G4double newIRad);
99 inline void SetOuterRadius (G4double newORad);
100 inline void SetZHalfLength (G4double newHLZ);
101 inline void SetInnerStereo (G4double newISte);
102 inline void SetOuterStereo (G4double newOSte);
103
104 EInside Inside(const G4ThreeVector& p) const;
105
107
108 G4double DistanceToIn(const G4ThreeVector& p, const G4ThreeVector& v) const;
109 G4double DistanceToIn(const G4ThreeVector& p) const;
111 const G4bool calcNorm = false,
112 G4bool* validNorm = nullptr,
113 G4ThreeVector* n = nullptr) const;
114 G4double DistanceToOut(const G4ThreeVector& p) const;
115
117
118 G4VSolid* Clone() const;
119
120 std::ostream& StreamInfo(std::ostream& os) const;
121
124
126
127 void DescribeYourselfTo (G4VGraphicsScene& scene) const;
128 G4VisExtent GetExtent () const;
130 G4Polyhedron* GetPolyhedron () const;
131
132 public: // without description
133
134 G4Hype(__void__&);
135 // Fake default constructor for usage restricted to direct object
136 // persistency for clients requiring preallocation of memory for
137 // persistifiable objects.
138
139 G4Hype(const G4Hype& rhs);
140 G4Hype& operator=(const G4Hype& rhs);
141 // Copy constructor and assignment operator.
142
143 protected: // without description
144
146 // whether we have an inner surface or not
147
149 G4double r0, G4double tanPhi );
151 G4double r0, G4double tan2Phi );
152 // approximate isotropic distance to hyperbolic surface
153
156 // values of hype radius at a given Z
157
158 static G4int IntersectHype( const G4ThreeVector &p, const G4ThreeVector &v,
159 G4double r2, G4double tan2Phi, G4double s[2] );
160 // intersection with hyperbolic surface
161
162 protected:
163
169
170 // precalculated parameters, squared quantities
171
174 G4double tanInnerStereo2; // squared tan of Inner Stereo angle
175 G4double tanOuterStereo2; // squared tan of Outer Stereo angle
176 G4double innerRadius2; // squared Inner Radius
177 G4double outerRadius2; // squared Outer Radius
178 G4double endInnerRadius2; // squared endcap Inner Radius
179 G4double endOuterRadius2; // squared endcap Outer Radius
180 G4double endInnerRadius; // endcap Inner Radius
181 G4double endOuterRadius; // endcap Outer Radius
182
183 // Used by distanceToOut
184
186
187 private:
188
189 G4double asinh(G4double arg);
190
191 private:
192
193 G4double fCubicVolume = 0.0;
194 G4double fSurfaceArea = 0.0;
195
196 G4double fHalfTol;
197
198 mutable G4bool fRebuildPolyhedron = false;
199 mutable G4Polyhedron* fpPolyhedron = nullptr;
200};
201
202#include "G4Hype.icc"
203
204#endif // defined(G4GEOM_USE_UHYPE) && defined(G4GEOM_USE_SYS_USOLIDS)
205
206#endif // G4HYPE_HH
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
Definition: G4Hype.hh:69
G4ThreeVector GetPointOnSurface() const
Definition: G4Hype.cc:1148
G4double endInnerRadius
Definition: G4Hype.hh:180
G4double tanOuterStereo
Definition: G4Hype.hh:173
G4Polyhedron * CreatePolyhedron() const
Definition: G4Hype.cc:1269
G4VisExtent GetExtent() const
Definition: G4Hype.cc:1258
G4Polyhedron * GetPolyhedron() const
Definition: G4Hype.cc:1277
G4double tanOuterStereo2
Definition: G4Hype.hh:175
static G4int IntersectHype(const G4ThreeVector &p, const G4ThreeVector &v, G4double r2, G4double tan2Phi, G4double s[2])
Definition: G4Hype.cc:921
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Hype.cc:279
G4double tanInnerStereo2
Definition: G4Hype.hh:174
G4double halfLenZ
Definition: G4Hype.hh:166
G4Hype & operator=(const G4Hype &rhs)
Definition: G4Hype.cc:159
G4double endInnerRadius2
Definition: G4Hype.hh:178
EInside Inside(const G4ThreeVector &p) const
Definition: G4Hype.cc:238
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
Definition: G4Hype.cc:221
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Hype.cc:190
G4VSolid * Clone() const
Definition: G4Hype.cc:1072
G4double innerRadius2
Definition: G4Hype.hh:176
static G4double ApproxDistOutside(G4double pr, G4double pz, G4double r0, G4double tanPhi)
Definition: G4Hype.cc:986
void SetOuterStereo(G4double newOSte)
G4double endOuterRadius
Definition: G4Hype.hh:181
G4double outerStereo
Definition: G4Hype.hh:168
G4double outerRadius
Definition: G4Hype.hh:165
G4double GetCubicVolume()
Definition: G4Hype.cc:1081
G4double tanInnerStereo
Definition: G4Hype.hh:172
G4double endOuterRadius2
Definition: G4Hype.hh:179
@ innerFace
Definition: G4Hype.hh:185
@ leftCap
Definition: G4Hype.hh:185
@ rightCap
Definition: G4Hype.hh:185
@ outerFace
Definition: G4Hype.hh:185
G4bool InnerSurfaceExists() const
void SetOuterRadius(G4double newORad)
G4double GetInnerStereo() const
G4double outerRadius2
Definition: G4Hype.hh:177
G4double HypeOuterRadius2(G4double zVal) const
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Hype.cc:323
G4double innerRadius
Definition: G4Hype.hh:164
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const
Definition: G4Hype.cc:693
virtual ~G4Hype()
Definition: G4Hype.cc:136
void SetZHalfLength(G4double newHLZ)
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4Hype.cc:199
G4double GetZHalfLength() const
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Hype.cc:1127
G4double innerStereo
Definition: G4Hype.hh:167
static G4double ApproxDistInside(G4double pr, G4double pz, G4double r0, G4double tan2Phi)
Definition: G4Hype.cc:1043
G4GeometryType GetEntityType() const
Definition: G4Hype.cc:1065
G4double HypeInnerRadius2(G4double zVal) const
G4double GetSurfaceArea()
Definition: G4Hype.cc:1093
G4double GetOuterStereo() const
G4double GetOuterRadius() const
void SetInnerStereo(G4double newISte)
void SetInnerRadius(G4double newIRad)
G4double GetInnerRadius() const
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Hype.cc:1251
EAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:67