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