Geant4-11
G4Tubs.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// G4Tubs
27//
28// Class description:
29//
30// A tube or tube segment with curved sides parallel to
31// the z-axis. The tube has a specified half-length along
32// the z-axis, about which it is centered, and a given
33// minimum and maximum radius. A minimum radius of 0
34// corresponds to filled tube /cylinder. The tube segment is
35// specified by starting and delta angles for phi, with 0
36// being the +x axis, PI/2 the +y axis.
37// A delta angle of 2PI signifies a complete, unsegmented
38// tube/cylinder.
39//
40// Member Data:
41//
42// fRMin Inner radius
43// fRMax Outer radius
44// fDz half length in z
45//
46// fSPhi The starting phi angle in radians,
47// adjusted such that fSPhi+fDPhi<=2PI, fSPhi>-2PI
48//
49// fDPhi Delta angle of the segment.
50//
51// fPhiFullTube Boolean variable used for indicate the Phi Section
52
53// 23.01.94 P.Kent: First version. Converted to `tolerant' geometry
54// --------------------------------------------------------------------
55#ifndef G4TUBS_HH
56#define G4TUBS_HH
57
58#include "G4GeomTypes.hh"
59
60#if defined(G4GEOM_USE_USOLIDS)
61#define G4GEOM_USE_UTUBS 1
62#endif
63
64#if defined(G4GEOM_USE_UTUBS)
65 #define G4UTubs G4Tubs
66 #include "G4UTubs.hh"
67#else
68
70
71#include "G4CSGSolid.hh"
72#include "G4Polyhedron.hh"
73
74class G4Tubs : public G4CSGSolid
75{
76 public: // with description
77
78 G4Tubs( const G4String& pName,
79 G4double pRMin,
80 G4double pRMax,
81 G4double pDz,
82 G4double pSPhi,
83 G4double pDPhi );
84 //
85 // Constructs a tubs with the given name and dimensions
86
87 virtual ~G4Tubs();
88 //
89 // Destructor
90
91 // Accessors
92
93 inline G4double GetInnerRadius () const;
94 inline G4double GetOuterRadius () const;
95 inline G4double GetZHalfLength () const;
96 inline G4double GetStartPhiAngle () const;
97 inline G4double GetDeltaPhiAngle () const;
98 inline G4double GetSinStartPhi () const;
99 inline G4double GetCosStartPhi () const;
100 inline G4double GetSinEndPhi () const;
101 inline G4double GetCosEndPhi () const;
102
103 // Modifiers
104
105 inline void SetInnerRadius (G4double newRMin);
106 inline void SetOuterRadius (G4double newRMax);
107 inline void SetZHalfLength (G4double newDz);
108 inline void SetStartPhiAngle (G4double newSPhi, G4bool trig=true);
109 inline void SetDeltaPhiAngle (G4double newDPhi);
110
111 // Methods for solid
112
115
117 const G4int n,
118 const G4VPhysicalVolume* pRep );
119
121
122 G4bool CalculateExtent( const EAxis pAxis,
123 const G4VoxelLimits& pVoxelLimit,
124 const G4AffineTransform& pTransform,
125 G4double& pmin, G4double& pmax ) const;
126
127 EInside Inside( const G4ThreeVector& p ) const;
128
129 G4ThreeVector SurfaceNormal( const G4ThreeVector& p ) const;
130
131 G4double DistanceToIn(const G4ThreeVector& p, const G4ThreeVector& v) const;
132 G4double DistanceToIn(const G4ThreeVector& p) const;
134 const G4bool calcNorm = false,
135 G4bool* validNorm = nullptr,
136 G4ThreeVector* n = nullptr) const;
137 G4double DistanceToOut(const G4ThreeVector& p) const;
138
140
142
143 G4VSolid* Clone() const;
144
145 std::ostream& StreamInfo( std::ostream& os ) const;
146
147 // Visualisation functions
148
149 void DescribeYourselfTo ( G4VGraphicsScene& scene ) const;
151
152 public: // without description
153
154 G4Tubs(__void__&);
155 //
156 // Fake default constructor for usage restricted to direct object
157 // persistency for clients requiring preallocation of memory for
158 // persistifiable objects.
159
160 G4Tubs(const G4Tubs& rhs);
161 G4Tubs& operator=(const G4Tubs& rhs);
162 // Copy constructor and assignment operator.
163
164 protected:
165
166 inline void Initialize();
167 //
168 // Reset relevant values to zero
169
170 inline void CheckSPhiAngle(G4double sPhi);
171 inline void CheckDPhiAngle(G4double dPhi);
172 inline void CheckPhiAngles(G4double sPhi, G4double dPhi);
173 //
174 // Reset relevant flags and angle values
175
177 //
178 // Recompute relevant trigonometric values and cache them
179
181 G4double normalTolerance ) const;
182 //
183 // Compute fast inverse cylindrical (Rxy) radius for points expected to
184 // be on a cylindrical surface. Ensures that surface normal vector
185 // produced has magnitude with 'normalTolerance' of unit
186
187 virtual G4ThreeVector ApproxSurfaceNormal( const G4ThreeVector& p ) const;
188 //
189 // Algorithm for SurfaceNormal() following the original
190 // specification for points not on the surface
191
192 protected:
193
194 // Used by distanceToOut
195 //
197
198 // Used by normal
199 //
201
203 //
204 // Radial and angular tolerances
205
206 static constexpr G4double kNormTolerance = 1.0e-6;
207 //
208 // Tolerance of unity for surface normal
209 // (for speedup - use fInvRmax if possible )
210
212 //
213 // Radial and angular dimensions
214
217 //
218 // Cached trigonometric values
219
221 //
222 // Flag for identification of section or full tube
223
225 //
226 // More cached values - inverse of Rmax, Rmin.
227
229 //
230 // Cached half tolerance values
231};
232
233#include "G4Tubs.icc"
234
235#endif
236
237#endif
static const G4double pos
static const G4double pMax
static const G4double pMin
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
Definition: G4Tubs.hh:75
G4GeometryType GetEntityType() const
Definition: G4Tubs.cc:1632
void CheckPhiAngles(G4double sPhi, G4double dPhi)
virtual G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector &p) const
Definition: G4Tubs.cc:598
G4Tubs(const G4String &pName, G4double pRMin, G4double pRMax, G4double pDz, G4double pSPhi, G4double pDPhi)
Definition: G4Tubs.cc:58
G4double GetZHalfLength() const
G4Polyhedron * CreatePolyhedron() const
Definition: G4Tubs.cc:1759
G4double GetSurfaceArea()
void CheckSPhiAngle(G4double sPhi)
void SetDeltaPhiAngle(G4double newDPhi)
EInside Inside(const G4ThreeVector &p) const
Definition: G4Tubs.cc:340
G4ThreeVector GetPointOnSurface() const
Definition: G4Tubs.cc:1673
G4double cosHDPhiIT
Definition: G4Tubs.hh:215
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4Tubs.cc:187
void InitializeTrigonometry()
G4double sinCPhi
Definition: G4Tubs.hh:215
G4double cosEPhi
Definition: G4Tubs.hh:216
void SetStartPhiAngle(G4double newSPhi, G4bool trig=true)
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Tubs.cc:176
G4VSolid * Clone() const
Definition: G4Tubs.cc:1641
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Tubs.cc:1650
G4double fRMin
Definition: G4Tubs.hh:211
G4double kAngTolerance
Definition: G4Tubs.hh:202
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const
Definition: G4Tubs.cc:1150
@ kEPhi
Definition: G4Tubs.hh:196
@ kRMax
Definition: G4Tubs.hh:196
@ kPZ
Definition: G4Tubs.hh:196
@ kMZ
Definition: G4Tubs.hh:196
@ kRMin
Definition: G4Tubs.hh:196
@ kSPhi
Definition: G4Tubs.hh:196
@ kNull
Definition: G4Tubs.hh:196
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4Tubs.cc:230
G4double fDPhi
Definition: G4Tubs.hh:211
G4double GetCosStartPhi() const
G4double GetCosEndPhi() const
G4double fRMax
Definition: G4Tubs.hh:211
G4double GetInnerRadius() const
G4double fInvRmin
Definition: G4Tubs.hh:224
G4Tubs & operator=(const G4Tubs &rhs)
Definition: G4Tubs.cc:142
G4double GetOuterRadius() const
G4double sinSPhi
Definition: G4Tubs.hh:216
G4double fInvRmax
Definition: G4Tubs.hh:224
static constexpr G4double kNormTolerance
Definition: G4Tubs.hh:206
void CheckDPhiAngle(G4double dPhi)
G4double cosCPhi
Definition: G4Tubs.hh:215
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Tubs.cc:1754
@ kNRMax
Definition: G4Tubs.hh:200
@ kNZ
Definition: G4Tubs.hh:200
@ kNRMin
Definition: G4Tubs.hh:200
@ kNEPhi
Definition: G4Tubs.hh:200
@ kNSPhi
Definition: G4Tubs.hh:200
G4double GetStartPhiAngle() const
G4double cosHDPhi
Definition: G4Tubs.hh:215
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Tubs.cc:730
G4double cosSPhi
Definition: G4Tubs.hh:216
G4double GetSinEndPhi() const
G4double sinEPhi
Definition: G4Tubs.hh:216
virtual ~G4Tubs()
Definition: G4Tubs.cc:114
void SetInnerRadius(G4double newRMin)
void SetOuterRadius(G4double newRMax)
void Initialize()
G4double kRadTolerance
Definition: G4Tubs.hh:202
G4double FastInverseRxy(const G4ThreeVector &pos, G4double invRad, G4double normalTolerance) const
G4double fDz
Definition: G4Tubs.hh:211
G4bool fPhiFullTube
Definition: G4Tubs.hh:220
G4double halfRadTolerance
Definition: G4Tubs.hh:228
G4double halfAngTolerance
Definition: G4Tubs.hh:228
G4double halfCarTolerance
Definition: G4Tubs.hh:228
G4double GetCubicVolume()
G4double GetDeltaPhiAngle() const
G4double GetSinStartPhi() const
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Tubs.cc:507
G4double cosHDPhiOT
Definition: G4Tubs.hh:215
void SetZHalfLength(G4double newDz)
G4double fSPhi
Definition: G4Tubs.hh:211
EAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:67