Geant4 9.6.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//
27// $Id$
28//
29//
30// --------------------------------------------------------------------
31// GEANT 4 class header file
32//
33//
34// G4Tubs
35//
36// Class description:
37//
38// A tube or tube segment with curved sides parallel to
39// the z-axis. The tube has a specified half-length along
40// the z-axis, about which it is centered, and a given
41// minimum and maximum radius. A minimum radius of 0
42// corresponds to filled tube /cylinder. The tube segment is
43// specified by starting and delta angles for phi, with 0
44// being the +x axis, PI/2 the +y axis.
45// A delta angle of 2PI signifies a complete, unsegmented
46// tube/cylinder.
47//
48// Member Data:
49//
50// fRMin Inner radius
51// fRMax Outer radius
52// fDz half length in z
53//
54// fSPhi The starting phi angle in radians,
55// adjusted such that fSPhi+fDPhi<=2PI, fSPhi>-2PI
56//
57// fDPhi Delta angle of the segment.
58//
59// fPhiFullTube Boolean variable used for indicate the Phi Section
60
61// History:
62// 10.08.95 P.Kent: General cleanup, use G4VSolid extent helper functions
63// to CalculateExtent()
64// 23.01.94 P.Kent: Converted to `tolerant' geometry
65// 19.07.96 J.Allison: G4GraphicsScene - see G4Box
66// 22.07.96 J.Allison: Changed SendPolyhedronTo to CreatePolyhedron
67// --------------------------------------------------------------------
68
69#ifndef G4TUBS_HH
70#define G4TUBS_HH
71
73
74#include "G4CSGSolid.hh"
75
76class G4Tubs : public G4CSGSolid
77{
78 public: // with description
79
80 G4Tubs( const G4String& pName,
81 G4double pRMin,
82 G4double pRMax,
83 G4double pDz,
84 G4double pSPhi,
85 G4double pDPhi );
86 //
87 // Constructs a tubs with the given name and dimensions
88
89 virtual ~G4Tubs();
90 //
91 // Destructor
92
93 // Accessors
94
95 inline G4double GetInnerRadius () const;
96 inline G4double GetOuterRadius () const;
97 inline G4double GetZHalfLength () const;
98 inline G4double GetStartPhiAngle () const;
99 inline G4double GetDeltaPhiAngle () const;
100
101 // Modifiers
102
103 inline void SetInnerRadius (G4double newRMin);
104 inline void SetOuterRadius (G4double newRMax);
105 inline void SetZHalfLength (G4double newDz);
106 inline void SetStartPhiAngle (G4double newSPhi, G4bool trig=true);
107 inline void SetDeltaPhiAngle (G4double newDPhi);
108
109 // Methods for solid
110
113
115 const G4int n,
116 const G4VPhysicalVolume* pRep );
117
118 G4bool CalculateExtent( const EAxis pAxis,
119 const G4VoxelLimits& pVoxelLimit,
120 const G4AffineTransform& pTransform,
121 G4double& pmin, G4double& pmax ) const;
122
123 EInside Inside( const G4ThreeVector& p ) const;
124
125 G4ThreeVector SurfaceNormal( const G4ThreeVector& p ) const;
126
127 G4double DistanceToIn(const G4ThreeVector& p, const G4ThreeVector& v) const;
128 G4double DistanceToIn(const G4ThreeVector& p) const;
130 const G4bool calcNorm=G4bool(false),
131 G4bool *validNorm=0, G4ThreeVector *n=0) const;
132 G4double DistanceToOut(const G4ThreeVector& p) const;
133
135
137
138 G4VSolid* Clone() const;
139
140 std::ostream& StreamInfo( std::ostream& os ) const;
141
142 // Visualisation functions
143
144 void DescribeYourselfTo ( G4VGraphicsScene& scene ) const;
146 G4NURBS* CreateNURBS () const;
147
148 public: // without description
149
150 G4Tubs(__void__&);
151 //
152 // Fake default constructor for usage restricted to direct object
153 // persistency for clients requiring preallocation of memory for
154 // persistifiable objects.
155
156 G4Tubs(const G4Tubs& rhs);
157 G4Tubs& operator=(const G4Tubs& rhs);
158 // Copy constructor and assignment operator.
159
160 // Older names for access functions
161
162 inline G4double GetRMin() const;
163 inline G4double GetRMax() const;
164 inline G4double GetDz () const;
165 inline G4double GetSPhi() const;
166 inline G4double GetDPhi() const;
167
168 protected:
169
171 CreateRotatedVertices( const G4AffineTransform& pTransform ) const;
172 //
173 // Creates the List of transformed vertices in the format required
174 // for G4VSolid:: ClipCrossSection and ClipBetweenSections
175
176 inline void Initialize();
177 //
178 // Reset relevant values to zero
179
180 inline void CheckSPhiAngle(G4double sPhi);
181 inline void CheckDPhiAngle(G4double dPhi);
182 inline void CheckPhiAngles(G4double sPhi, G4double dPhi);
183 //
184 // Reset relevant flags and angle values
185
187 //
188 // Recompute relevant trigonometric values and cache them
189
190 virtual G4ThreeVector ApproxSurfaceNormal( const G4ThreeVector& p ) const;
191 //
192 // Algorithm for SurfaceNormal() following the original
193 // specification for points not on the surface
194
195 protected:
196
197 // Used by distanceToOut
198 //
200
201 // Used by normal
202 //
204
206 //
207 // Radial and angular tolerances
208
210 //
211 // Radial and angular dimensions
212
215 //
216 // Cached trigonometric values
217
219 //
220 // Flag for identification of section or full tube
221};
222
223#include "G4Tubs.icc"
224
225#endif
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4VSolid.hh:85
Definition: G4Tubs.hh:77
G4GeometryType GetEntityType() const
Definition: G4Tubs.cc:1807
void CheckPhiAngles(G4double sPhi, G4double dPhi)
virtual G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector &p) const
Definition: G4Tubs.cc:678
G4double GetZHalfLength() const
G4Polyhedron * CreatePolyhedron() const
Definition: G4Tubs.cc:1923
G4double GetSurfaceArea()
void CheckSPhiAngle(G4double sPhi)
void SetDeltaPhiAngle(G4double newDPhi)
EInside Inside(const G4ThreeVector &p) const
Definition: G4Tubs.cc:414
G4ThreeVector GetPointOnSurface() const
Definition: G4Tubs.cc:1848
G4double cosHDPhiIT
Definition: G4Tubs.hh:213
G4double GetRMin() const
void InitializeTrigonometry()
G4double sinCPhi
Definition: G4Tubs.hh:213
G4double cosEPhi
Definition: G4Tubs.hh:214
void SetStartPhiAngle(G4double newSPhi, G4bool trig=true)
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Tubs.cc:189
G4VSolid * Clone() const
Definition: G4Tubs.cc:1816
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Tubs.cc:1825
G4double GetSPhi() const
G4double fRMin
Definition: G4Tubs.hh:209
G4ThreeVectorList * CreateRotatedVertices(const G4AffineTransform &pTransform) const
Definition: G4Tubs.cc:1721
G4double kAngTolerance
Definition: G4Tubs.hh:205
@ kEPhi
Definition: G4Tubs.hh:199
@ kRMax
Definition: G4Tubs.hh:199
@ kPZ
Definition: G4Tubs.hh:199
@ kMZ
Definition: G4Tubs.hh:199
@ kRMin
Definition: G4Tubs.hh:199
@ kSPhi
Definition: G4Tubs.hh:199
@ kNull
Definition: G4Tubs.hh:199
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4Tubs.cc:200
G4double fDPhi
Definition: G4Tubs.hh:209
G4double fRMax
Definition: G4Tubs.hh:209
G4double GetInnerRadius() const
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=G4bool(false), G4bool *validNorm=0, G4ThreeVector *n=0) const
Definition: G4Tubs.cc:1233
G4Tubs & operator=(const G4Tubs &rhs)
Definition: G4Tubs.cc:160
G4double GetOuterRadius() const
G4NURBS * CreateNURBS() const
Definition: G4Tubs.cc:1928
G4double sinSPhi
Definition: G4Tubs.hh:214
void CheckDPhiAngle(G4double dPhi)
G4double cosCPhi
Definition: G4Tubs.hh:213
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Tubs.cc:1918
@ kNRMax
Definition: G4Tubs.hh:203
@ kNZ
Definition: G4Tubs.hh:203
@ kNRMin
Definition: G4Tubs.hh:203
@ kNEPhi
Definition: G4Tubs.hh:203
@ kNSPhi
Definition: G4Tubs.hh:203
G4double GetStartPhiAngle() const
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Tubs.cc:810
G4double GetRMax() const
G4double cosSPhi
Definition: G4Tubs.hh:214
G4double sinEPhi
Definition: G4Tubs.hh:214
virtual ~G4Tubs()
Definition: G4Tubs.cc:136
void SetInnerRadius(G4double newRMin)
void SetOuterRadius(G4double newRMax)
void Initialize()
G4double GetDz() const
G4double GetDPhi() const
G4double kRadTolerance
Definition: G4Tubs.hh:205
G4double fDz
Definition: G4Tubs.hh:209
G4bool fPhiFullTube
Definition: G4Tubs.hh:218
G4double GetCubicVolume()
G4double GetDeltaPhiAngle() const
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Tubs.cc:584
G4double cosHDPhiOT
Definition: G4Tubs.hh:213
void SetZHalfLength(G4double newDz)
G4double fSPhi
Definition: G4Tubs.hh:209
EAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:58