Geant4  10.01.p02
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: G4Tubs.hh 79491 2014-03-05 15:24:29Z gcosmo $
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 #ifndef G4TUBS_HH
69 #define G4TUBS_HH
70 
71 #if defined(G4GEOM_USE_USOLIDS)
72 #define G4GEOM_USE_UTUBS 1
73 #endif
74 
75 #if defined(G4GEOM_USE_UTUBS)
76  #define G4UTubs G4Tubs
77  #include "G4UTubs.hh"
78 #else
79 
80 #include <CLHEP/Units/PhysicalConstants.h>
81 
82 #include "G4CSGSolid.hh"
83 #include "G4Polyhedron.hh"
84 
85 class G4Tubs : public G4CSGSolid
86 {
87  public: // with description
88 
89  G4Tubs( const G4String& pName,
90  G4double pRMin,
91  G4double pRMax,
92  G4double pDz,
93  G4double pSPhi,
94  G4double pDPhi );
95  //
96  // Constructs a tubs with the given name and dimensions
97 
98  virtual ~G4Tubs();
99  //
100  // Destructor
101 
102  // Accessors
103 
104  inline G4double GetInnerRadius () const;
105  inline G4double GetOuterRadius () const;
106  inline G4double GetZHalfLength () const;
107  inline G4double GetStartPhiAngle () const;
108  inline G4double GetDeltaPhiAngle () const;
109 
110  // Modifiers
111 
112  inline void SetInnerRadius (G4double newRMin);
113  inline void SetOuterRadius (G4double newRMax);
114  inline void SetZHalfLength (G4double newDz);
115  inline void SetStartPhiAngle (G4double newSPhi, G4bool trig=true);
116  inline void SetDeltaPhiAngle (G4double newDPhi);
117 
118  // Methods for solid
119 
120  inline G4double GetCubicVolume();
121  inline G4double GetSurfaceArea();
122 
124  const G4int n,
125  const G4VPhysicalVolume* pRep );
126 
127  G4bool CalculateExtent( const EAxis pAxis,
128  const G4VoxelLimits& pVoxelLimit,
129  const G4AffineTransform& pTransform,
130  G4double& pmin, G4double& pmax ) const;
131 
132  EInside Inside( const G4ThreeVector& p ) const;
133 
134  G4ThreeVector SurfaceNormal( const G4ThreeVector& p ) const;
135 
136  G4double DistanceToIn(const G4ThreeVector& p, const G4ThreeVector& v) const;
137  G4double DistanceToIn(const G4ThreeVector& p) const;
139  const G4bool calcNorm=G4bool(false),
140  G4bool *validNorm=0, G4ThreeVector *n=0) const;
141  G4double DistanceToOut(const G4ThreeVector& p) const;
142 
144 
146 
147  G4VSolid* Clone() const;
148 
149  std::ostream& StreamInfo( std::ostream& os ) const;
150 
151  // Visualisation functions
152 
153  void DescribeYourselfTo ( G4VGraphicsScene& scene ) const;
154  G4Polyhedron* CreatePolyhedron () const;
155 
156  public: // without description
157 
158  G4Tubs(__void__&);
159  //
160  // Fake default constructor for usage restricted to direct object
161  // persistency for clients requiring preallocation of memory for
162  // persistifiable objects.
163 
164  G4Tubs(const G4Tubs& rhs);
165  G4Tubs& operator=(const G4Tubs& rhs);
166  // Copy constructor and assignment operator.
167 
168  // Older names for access functions
169 
170  inline G4double GetRMin() const;
171  inline G4double GetRMax() const;
172  inline G4double GetDz () const;
173  inline G4double GetSPhi() const;
174  inline G4double GetDPhi() const;
175 
176  protected:
177 
179  CreateRotatedVertices( const G4AffineTransform& pTransform ) const;
180  //
181  // Creates the List of transformed vertices in the format required
182  // for G4VSolid:: ClipCrossSection and ClipBetweenSections
183 
184  inline void Initialize();
185  //
186  // Reset relevant values to zero
187 
188  inline void CheckSPhiAngle(G4double sPhi);
189  inline void CheckDPhiAngle(G4double dPhi);
190  inline void CheckPhiAngles(G4double sPhi, G4double dPhi);
191  //
192  // Reset relevant flags and angle values
193 
194  inline void InitializeTrigonometry();
195  //
196  // Recompute relevant trigonometric values and cache them
197 
198  virtual G4ThreeVector ApproxSurfaceNormal( const G4ThreeVector& p ) const;
199  //
200  // Algorithm for SurfaceNormal() following the original
201  // specification for points not on the surface
202 
203  protected:
204 
205  // Used by distanceToOut
206  //
208 
209  // Used by normal
210  //
212 
214  //
215  // Radial and angular tolerances
216 
218  //
219  // Radial and angular dimensions
220 
223  //
224  // Cached trigonometric values
225 
227  //
228  // Flag for identification of section or full tube
229 
231  //
232  // Cached half tolerance values
233 };
234 
235 #include "G4Tubs.icc"
236 
237 #endif
238 
239 #endif
G4bool fPhiFullTube
Definition: G4Tubs.hh:226
G4double GetRMax() const
CLHEP::Hep3Vector G4ThreeVector
void CheckSPhiAngle(G4double sPhi)
G4double fDPhi
Definition: G4Tubs.hh:217
G4double GetSurfaceArea()
G4double halfCarTolerance
Definition: G4Tubs.hh:230
G4Tubs(const G4String &pName, G4double pRMin, G4double pRMax, G4double pDz, G4double pSPhi, G4double pDPhi)
Definition: G4Tubs.cc:85
Definition: G4Tubs.hh:85
G4double cosEPhi
Definition: G4Tubs.hh:221
void SetStartPhiAngle(G4double newSPhi, G4bool trig=true)
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Tubs.cc:812
G4double GetCubicVolume()
G4double GetDz() const
void Initialize()
G4ThreeVector GetPointOnSurface() const
Definition: G4Tubs.cc:1844
virtual G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector &p) const
Definition: G4Tubs.cc:680
int G4int
Definition: G4Types.hh:78
G4double cosHDPhiOT
Definition: G4Tubs.hh:221
G4double GetSPhi() const
void CheckDPhiAngle(G4double dPhi)
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Tubs.cc:1821
void SetDeltaPhiAngle(G4double newDPhi)
G4double GetDPhi() const
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Tubs.cc:1914
G4double sinSPhi
Definition: G4Tubs.hh:221
G4GeometryType GetEntityType() const
Definition: G4Tubs.cc:1803
G4double GetDeltaPhiAngle() const
G4double fRMin
Definition: G4Tubs.hh:217
G4Polyhedron * CreatePolyhedron() const
Definition: G4Tubs.cc:1919
bool G4bool
Definition: G4Types.hh:79
G4VSolid * Clone() const
Definition: G4Tubs.cc:1812
G4double cosCPhi
Definition: G4Tubs.hh:221
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4VSolid.hh:79
G4double fSPhi
Definition: G4Tubs.hh:217
G4double GetStartPhiAngle() const
const G4int n
G4double GetRMin() const
G4double GetInnerRadius() const
virtual ~G4Tubs()
Definition: G4Tubs.cc:138
void SetInnerRadius(G4double newRMin)
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4Tubs.cc:208
G4double sinEPhi
Definition: G4Tubs.hh:221
G4double halfAngTolerance
Definition: G4Tubs.hh:230
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=G4bool(false), G4bool *validNorm=0, G4ThreeVector *n=0) const
Definition: G4Tubs.cc:1232
void SetOuterRadius(G4double newRMax)
EInside Inside(const G4ThreeVector &p) const
Definition: G4Tubs.cc:422
G4double halfRadTolerance
Definition: G4Tubs.hh:230
EInside
Definition: geomdefs.hh:58
EAxis
Definition: geomdefs.hh:54
G4double GetZHalfLength() const
G4double fRMax
Definition: G4Tubs.hh:217
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Tubs.cc:589
G4double kRadTolerance
Definition: G4Tubs.hh:213
G4double fDz
Definition: G4Tubs.hh:217
G4double cosSPhi
Definition: G4Tubs.hh:221
double G4double
Definition: G4Types.hh:76
void InitializeTrigonometry()
void CheckPhiAngles(G4double sPhi, G4double dPhi)
G4double cosHDPhiIT
Definition: G4Tubs.hh:221
G4double sinCPhi
Definition: G4Tubs.hh:221
G4double GetOuterRadius() const
G4double kAngTolerance
Definition: G4Tubs.hh:213
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Tubs.cc:197
G4Tubs & operator=(const G4Tubs &rhs)
Definition: G4Tubs.cc:165
void SetZHalfLength(G4double newDz)
G4ThreeVectorList * CreateRotatedVertices(const G4AffineTransform &pTransform) const
Definition: G4Tubs.cc:1717