Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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:
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 ~G4Tubs() override;
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
113 inline G4double GetCubicVolume() override;
114 inline G4double GetSurfaceArea() override;
115
117 const G4int n,
118 const G4VPhysicalVolume* pRep ) override;
119
120 void BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const override;
121
122 G4bool CalculateExtent(const EAxis pAxis,
123 const G4VoxelLimits& pVoxelLimit,
124 const G4AffineTransform& pTransform,
125 G4double& pmin, G4double& pmax) const override;
126
127 EInside Inside( const G4ThreeVector& p ) const override;
128
129 G4ThreeVector SurfaceNormal( const G4ThreeVector& p ) const override;
130
132 const G4ThreeVector& v) const override;
133 G4double DistanceToIn(const G4ThreeVector& p) const override;
135 const G4bool calcNorm = false,
136 G4bool* validNorm = nullptr,
137 G4ThreeVector* n = nullptr) const override;
138 G4double DistanceToOut(const G4ThreeVector& p) const override;
139
140 G4GeometryType GetEntityType() const override;
141
142 G4ThreeVector GetPointOnSurface() const override;
143
144 G4VSolid* Clone() const override;
145
146 std::ostream& StreamInfo( std::ostream& os ) const override;
147
148 // Visualisation functions
149
150 void DescribeYourselfTo (G4VGraphicsScene& scene) const override;
151 G4Polyhedron* CreatePolyhedron () const override;
152
153 G4Tubs(__void__&);
154 //
155 // Fake default constructor for usage restricted to direct object
156 // persistency for clients requiring preallocation of memory for
157 // persistifiable objects.
158
159 G4Tubs(const G4Tubs& rhs);
160 G4Tubs& operator=(const G4Tubs& rhs);
161 // Copy constructor and assignment operator.
162
163 protected:
164
165 inline void Initialize();
166 //
167 // Reset relevant values to zero
168
169 inline void CheckSPhiAngle(G4double sPhi);
170 inline void CheckDPhiAngle(G4double dPhi);
171 inline void CheckPhiAngles(G4double sPhi, G4double dPhi);
172 //
173 // Reset relevant flags and angle values
174
176 //
177 // Recompute relevant trigonometric values and cache them
178
179 inline G4double FastInverseRxy( const G4ThreeVector& pos, G4double invRad,
180 G4double normalTolerance ) const;
181 //
182 // Compute fast inverse cylindrical (Rxy) radius for points expected to
183 // be on a cylindrical surface. Ensures that surface normal vector
184 // produced has magnitude with 'normalTolerance' of unit
185
186 virtual G4ThreeVector ApproxSurfaceNormal( const G4ThreeVector& p ) const;
187 //
188 // Algorithm for SurfaceNormal() following the original
189 // specification for points not on the surface
190
191 protected:
192
193 // Used by distanceToOut
194 //
196
197 // Used by normal
198 //
200
202 //
203 // Radial and angular tolerances
204
205 static constexpr G4double kNormTolerance = 1.0e-6;
206 //
207 // Tolerance of unity for surface normal
208 // (for speedup - use fInvRmax if possible )
209
211 //
212 // Radial and angular dimensions
213
216 //
217 // Cached trigonometric values
218
220 //
221 // Flag for identification of section or full tube
222
224 //
225 // More cached values - inverse of Rmax, Rmin.
226
228 //
229 // Cached half tolerance values
230};
231
232#include "G4Tubs.icc"
233
234#endif
235
236#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void CheckPhiAngles(G4double sPhi, G4double dPhi)
virtual G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector &p) const
Definition G4Tubs.cc:576
G4Tubs(const G4String &pName, G4double pRMin, G4double pRMax, G4double pDz, G4double pSPhi, G4double pDPhi)
Definition G4Tubs.cc:58
G4double GetZHalfLength() const
void CheckSPhiAngle(G4double sPhi)
void SetDeltaPhiAngle(G4double newDPhi)
G4double cosHDPhiIT
Definition G4Tubs.hh:214
void InitializeTrigonometry()
G4double sinCPhi
Definition G4Tubs.hh:214
G4double cosEPhi
Definition G4Tubs.hh:215
G4ThreeVector GetPointOnSurface() const override
Definition G4Tubs.cc:1651
void SetStartPhiAngle(G4double newSPhi, G4bool trig=true)
G4double GetSurfaceArea() override
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override
Definition G4Tubs.cc:165
G4double fRMin
Definition G4Tubs.hh:210
G4double kAngTolerance
Definition G4Tubs.hh:201
@ kEPhi
Definition G4Tubs.hh:195
@ kRMax
Definition G4Tubs.hh:195
@ kPZ
Definition G4Tubs.hh:195
@ kMZ
Definition G4Tubs.hh:195
@ kRMin
Definition G4Tubs.hh:195
@ kSPhi
Definition G4Tubs.hh:195
@ kNull
Definition G4Tubs.hh:195
G4double fDPhi
Definition G4Tubs.hh:210
G4double GetCosStartPhi() const
~G4Tubs() override
void DescribeYourselfTo(G4VGraphicsScene &scene) const override
Definition G4Tubs.cc:1732
G4double GetCosEndPhi() const
G4double fRMax
Definition G4Tubs.hh:210
G4double GetInnerRadius() const
G4double fInvRmin
Definition G4Tubs.hh:223
G4double GetCubicVolume() override
G4Tubs & operator=(const G4Tubs &rhs)
Definition G4Tubs.cc:120
G4double GetOuterRadius() const
G4double sinSPhi
Definition G4Tubs.hh:215
G4double fInvRmax
Definition G4Tubs.hh:223
static constexpr G4double kNormTolerance
Definition G4Tubs.hh:205
void CheckDPhiAngle(G4double dPhi)
G4double cosCPhi
Definition G4Tubs.hh:214
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const override
Definition G4Tubs.cc:485
@ kNRMax
Definition G4Tubs.hh:199
@ kNZ
Definition G4Tubs.hh:199
@ kNRMin
Definition G4Tubs.hh:199
@ kNEPhi
Definition G4Tubs.hh:199
@ kNSPhi
Definition G4Tubs.hh:199
G4double GetStartPhiAngle() const
EInside Inside(const G4ThreeVector &p) const override
Definition G4Tubs.cc:318
G4double cosHDPhi
Definition G4Tubs.hh:214
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const override
Definition G4Tubs.cc:208
G4Polyhedron * CreatePolyhedron() const override
Definition G4Tubs.cc:1737
G4double cosSPhi
Definition G4Tubs.hh:215
G4double GetSinEndPhi() const
G4double sinEPhi
Definition G4Tubs.hh:215
std::ostream & StreamInfo(std::ostream &os) const override
Definition G4Tubs.cc:1628
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const override
Definition G4Tubs.cc:1128
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep) override
Definition G4Tubs.cc:154
void SetInnerRadius(G4double newRMin)
void SetOuterRadius(G4double newRMax)
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override
Definition G4Tubs.cc:708
void Initialize()
G4double kRadTolerance
Definition G4Tubs.hh:201
G4double FastInverseRxy(const G4ThreeVector &pos, G4double invRad, G4double normalTolerance) const
G4double fDz
Definition G4Tubs.hh:210
G4Tubs(const G4Tubs &rhs)
G4bool fPhiFullTube
Definition G4Tubs.hh:219
G4GeometryType GetEntityType() const override
Definition G4Tubs.cc:1610
G4double halfRadTolerance
Definition G4Tubs.hh:227
G4double halfAngTolerance
Definition G4Tubs.hh:227
G4double halfCarTolerance
Definition G4Tubs.hh:227
G4double GetDeltaPhiAngle() const
G4double GetSinStartPhi() const
G4VSolid * Clone() const override
Definition G4Tubs.cc:1619
G4double cosHDPhiOT
Definition G4Tubs.hh:214
void SetZHalfLength(G4double newDz)
G4double fSPhi
Definition G4Tubs.hh:210
EAxis
Definition geomdefs.hh:54
EInside
Definition geomdefs.hh:67