Geant4 11.4.0
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// Author: Paul Kent (CERN), 23.01.1994 - First version
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
74/**
75 * @brief G4Tubs is a tube or tube segment with curved sides parallel to
76 * the Z-axis. The tube has a specified half-length along the Z-axis, about
77 * which it is centered, and a given minimum and maximum radius. A minimum
78 * radius of 0 corresponds to filled tube/cylinder. The tube segment is
79 * specified by starting and delta angles for phi, with 0 being the +x axis,
80 * PI/2 the +y axis. A delta angle of 2PI signifies a complete, unsegmented
81 * tube/cylinder.
82 */
83
84class G4Tubs : public G4CSGSolid
85{
86 public:
87
88 /**
89 * Constructs a tubs with the given name and dimensions.
90 * It checks the input parameters, converting angles so 0<sphi+dpshi<=2_PI
91 * if pdphi>2PI then reset it to 2PI.
92 * @param[in] pName The name of the solid.
93 * @param[in] pRMin Inner radius.
94 * @param[in] pRMax Outer radius.
95 * @param[in] pDz Half length in Z.
96 * @param[in] pSPhi Starting phi angle in radians.
97 * @param[in] pDPhi Angle of the segment in radians.
98 */
99 G4Tubs( const G4String& pName,
100 G4double pRMin,
101 G4double pRMax,
102 G4double pDz,
103 G4double pSPhi,
104 G4double pDPhi );
105
106 /**
107 * Default destructor.
108 */
109 ~G4Tubs() override = default;
110
111 /**
112 * Accessors.
113 */
114 inline G4double GetInnerRadius () const;
115 inline G4double GetOuterRadius () const;
116 inline G4double GetZHalfLength () const;
117 inline G4double GetStartPhiAngle () const;
118 inline G4double GetDeltaPhiAngle () const;
119 inline G4double GetSinStartPhi () const;
120 inline G4double GetCosStartPhi () const;
121 inline G4double GetSinEndPhi () const;
122 inline G4double GetCosEndPhi () const;
123
124 /**
125 * Modifiers.
126 */
127 inline void SetInnerRadius (G4double newRMin);
128 inline void SetOuterRadius (G4double newRMax);
129 inline void SetZHalfLength (G4double newDz);
130 inline void SetStartPhiAngle (G4double newSPhi, G4bool trig=true);
131 inline void SetDeltaPhiAngle (G4double newDPhi);
132
133 /**
134 * Returning an estimation of the solid volume (capacity) and
135 * surface area, in internal units.
136 */
137 G4double GetCubicVolume() override;
138 G4double GetSurfaceArea() override;
139
140 /**
141 * Dispatch method for parameterisation replication mechanism and
142 * dimension computation.
143 */
145 const G4int n,
146 const G4VPhysicalVolume* pRep ) override;
147
148 /**
149 * Computes the bounding limits of the solid.
150 * @param[out] pMin The minimum bounding limit point.
151 * @param[out] pMax The maximum bounding limit point.
152 */
153 void BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const override;
154
155 /**
156 * Calculates the minimum and maximum extent of the solid, when under the
157 * specified transform, and within the specified limits.
158 * @param[in] pAxis The axis along which compute the extent.
159 * @param[in] pVoxelLimit The limiting space dictated by voxels.
160 * @param[in] pTransform The internal transformation applied to the solid.
161 * @param[out] pMin The minimum extent value.
162 * @param[out] pMax The maximum extent value.
163 * @returns True if the solid is intersected by the extent region.
164 */
165 G4bool CalculateExtent(const EAxis pAxis,
166 const G4VoxelLimits& pVoxelLimit,
167 const G4AffineTransform& pTransform,
168 G4double& pmin, G4double& pmax) const override;
169
170 /**
171 * Concrete implementations of the expected query interfaces for
172 * solids, as defined in the base class G4VSolid.
173 */
174 EInside Inside( const G4ThreeVector& p ) const override;
175 G4ThreeVector SurfaceNormal( const G4ThreeVector& p ) const override;
177 const G4ThreeVector& v) const override;
178 G4double DistanceToIn(const G4ThreeVector& p) const override;
180 const G4bool calcNorm = false,
181 G4bool* validNorm = nullptr,
182 G4ThreeVector* n = nullptr) const override;
183 G4double DistanceToOut(const G4ThreeVector& p) const override;
184
185 /**
186 * Returns the type ID, "G4Tubs" of the solid.
187 */
188 G4GeometryType GetEntityType() const override;
189
190 /**
191 * Returns a random point located and uniformly distributed on the
192 * surface of the solid.
193 */
194 G4ThreeVector GetPointOnSurface() const override;
195
196 /**
197 * Makes a clone of the object for use in multi-treading.
198 * @returns A pointer to the new cloned allocated solid.
199 */
200 G4VSolid* Clone() const override;
201
202 /**
203 * Streams the object contents to an output stream.
204 */
205 std::ostream& StreamInfo( std::ostream& os ) const override;
206
207 /**
208 * Methods for creating graphical representations (i.e. for visualisation).
209 */
210 void DescribeYourselfTo (G4VGraphicsScene& scene) const override;
211 G4Polyhedron* CreatePolyhedron () const override;
212
213 /**
214 * Fake default constructor for usage restricted to direct object
215 * persistency for clients requiring preallocation of memory for
216 * persistifiable objects.
217 */
218 G4Tubs(__void__&);
219
220 /**
221 * Copy constructor and assignment operator.
222 */
223 G4Tubs(const G4Tubs& rhs) = default;
224 G4Tubs& operator=(const G4Tubs& rhs);
225
226 protected:
227
228 /**
229 * Resets the relevant values to zero.
230 */
231 inline void Initialize();
232
233 /**
234 * Methods resetting relevant flags and angle values.
235 */
236 inline void CheckSPhiAngle(G4double sPhi);
237 inline void CheckDPhiAngle(G4double dPhi);
238 inline void CheckPhiAngles(G4double sPhi, G4double dPhi);
239
240 /**
241 * Recomputes relevant trigonometric values and caches them.
242 */
244
245 /**
246 * Computes fast inverse cylindrical (Rxy) radius for points expected to
247 * be on a cylindrical surface. Ensures that surface normal vector
248 * produced has magnitude with 'normalTolerance' of unit.
249 */
250 inline G4double FastInverseRxy( const G4ThreeVector& pos, G4double invRad,
251 G4double normalTolerance ) const;
252
253 /**
254 * Algorithm for SurfaceNormal() following the original specification
255 * for points not on the surface.
256 */
258
259 protected:
260
261 /** Radial and angular tolerances. */
263
264 /** Tolerance of unity for surface normal. */
265 static constexpr G4double kNormTolerance = 1.0e-6;
266
267 /** Radial and angular dimensions. */
269
270 /** Cached trigonometric values. */
273
274 /** Flag for identification of section or full tube. */
276
277 /** More cached values - inverse of Rmax, Rmin. */
279
280 /** Cached half tolerance values. */
282};
283
284#include "G4Tubs.icc"
285
286#endif
287
288#endif
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4String G4GeometryType
Definition G4VSolid.hh:70
G4AffineTransform is a class for geometric affine transformations. It supports efficient arbitrary ro...
G4CSGSolid(const G4String &pName)
Definition G4CSGSolid.cc:49
void CheckPhiAngles(G4double sPhi, G4double dPhi)
G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector &p) const
Definition G4Tubs.cc:585
G4Tubs(const G4String &pName, G4double pRMin, G4double pRMax, G4double pDz, G4double pSPhi, G4double pDPhi)
Definition G4Tubs.cc:73
G4double GetZHalfLength() const
void CheckSPhiAngle(G4double sPhi)
void SetDeltaPhiAngle(G4double newDPhi)
~G4Tubs() override=default
G4double cosHDPhiIT
Definition G4Tubs.hh:271
void InitializeTrigonometry()
G4Tubs(const G4Tubs &rhs)=default
G4double sinCPhi
Definition G4Tubs.hh:271
G4double cosEPhi
Definition G4Tubs.hh:272
G4ThreeVector GetPointOnSurface() const override
Definition G4Tubs.cc:1650
void SetStartPhiAngle(G4double newSPhi, G4bool trig=true)
G4double GetSurfaceArea() override
Definition G4Tubs.cc:1746
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override
Definition G4Tubs.cc:168
G4double fRMin
Definition G4Tubs.hh:268
G4double kAngTolerance
Definition G4Tubs.hh:262
G4double fDPhi
Definition G4Tubs.hh:268
G4double GetCosStartPhi() const
void DescribeYourselfTo(G4VGraphicsScene &scene) const override
Definition G4Tubs.cc:1765
G4double GetCosEndPhi() const
G4double fRMax
Definition G4Tubs.hh:268
G4double GetInnerRadius() const
G4double fInvRmin
Definition G4Tubs.hh:278
G4double GetCubicVolume() override
Definition G4Tubs.cc:1731
G4Tubs & operator=(const G4Tubs &rhs)
Definition G4Tubs.cc:123
G4double GetOuterRadius() const
G4double sinSPhi
Definition G4Tubs.hh:272
G4double fInvRmax
Definition G4Tubs.hh:278
static constexpr G4double kNormTolerance
Definition G4Tubs.hh:265
void CheckDPhiAngle(G4double dPhi)
G4double cosCPhi
Definition G4Tubs.hh:271
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const override
Definition G4Tubs.cc:494
G4double GetStartPhiAngle() const
EInside Inside(const G4ThreeVector &p) const override
Definition G4Tubs.cc:327
G4double cosHDPhi
Definition G4Tubs.hh:271
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const override
Definition G4Tubs.cc:211
G4Polyhedron * CreatePolyhedron() const override
Definition G4Tubs.cc:1770
G4double cosSPhi
Definition G4Tubs.hh:272
G4double GetSinEndPhi() const
G4double sinEPhi
Definition G4Tubs.hh:272
std::ostream & StreamInfo(std::ostream &os) const override
Definition G4Tubs.cc:1627
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const override
Definition G4Tubs.cc:1127
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep) override
Definition G4Tubs.cc:157
void SetInnerRadius(G4double newRMin)
void SetOuterRadius(G4double newRMax)
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override
Definition G4Tubs.cc:717
void Initialize()
G4double kRadTolerance
Definition G4Tubs.hh:262
G4double FastInverseRxy(const G4ThreeVector &pos, G4double invRad, G4double normalTolerance) const
G4double fDz
Definition G4Tubs.hh:268
G4bool fPhiFullTube
Definition G4Tubs.hh:275
G4GeometryType GetEntityType() const override
Definition G4Tubs.cc:1609
G4double halfRadTolerance
Definition G4Tubs.hh:281
G4double halfAngTolerance
Definition G4Tubs.hh:281
G4double halfCarTolerance
Definition G4Tubs.hh:281
G4double GetDeltaPhiAngle() const
G4double GetSinStartPhi() const
G4VSolid * Clone() const override
Definition G4Tubs.cc:1618
G4double cosHDPhiOT
Definition G4Tubs.hh:271
void SetZHalfLength(G4double newDz)
G4double fSPhi
Definition G4Tubs.hh:268
G4VPVParameterisation ia an abstract base class for Parameterisation, able to compute the transformat...
G4VPhysicalVolume is an abstract base class for the representation of a positioned volume....
G4VSolid(const G4String &name)
Definition G4VSolid.cc:59
G4VoxelLimits represents limitation/restrictions of space, where restrictions are only made perpendic...
EAxis
Definition geomdefs.hh:54
EInside
Definition geomdefs.hh:67