Geant4  10.00.p01
G4Sphere.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: G4Sphere.hh 76263 2013-11-08 11:41:52Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------------
31 // GEANT 4 class header file
32 //
33 // G4Sphere
34 //
35 // Class description:
36 //
37 // A G4Sphere is, in the general case, a section of a spherical shell,
38 // between specified phi and theta angles
39 //
40 // The phi and theta segments are described by a starting angle,
41 // and the +ve delta angle for the shape.
42 // If the delta angle is >=2*pi, or >=pi the shape is treated as
43 // continuous in phi or theta respectively.
44 //
45 // Theta must lie between 0-pi (incl).
46 //
47 // Member Data:
48 //
49 // fRmin inner radius
50 // fRmax outer radius
51 //
52 // fSPhi starting angle of the segment in radians
53 // fDPhi delta angle of the segment in radians
54 //
55 // fSTheta starting angle of the segment in radians
56 // fDTheta delta angle of the segment in radians
57 //
58 //
59 // Note:
60 // Internally fSPhi & fDPhi are adjusted so that fDPhi<=2PI,
61 // and fDPhi+fSPhi<=2PI. This enables simpler comparisons to be
62 // made with (say) Phi of a point.
63 
64 // History:
65 // 28.3.94 P.Kent: old C++ code converted to tolerant geometry
66 // 17.9.96 V.Grichine: final modifications to commit
67 // --------------------------------------------------------------------
68 
69 #ifndef G4Sphere_HH
70 #define G4Sphere_HH
71 
72 #if defined(G4GEOM_USE_USOLIDS)
73 #define G4GEOM_USE_USPHERE 1
74 #endif
75 
76 #if defined(G4GEOM_USE_USPHERE)
77  #define G4USphere G4Sphere
78  #include "G4USphere.hh"
79 #else
80 
81 #include <CLHEP/Units/PhysicalConstants.h>
82 #include "G4CSGSolid.hh"
83 
84 class G4VisExtent;
85 
86 class G4Sphere : public G4CSGSolid
87 {
88  public: // with description
89 
90  G4Sphere(const G4String& pName,
91  G4double pRmin, G4double pRmax,
92  G4double pSPhi, G4double pDPhi,
93  G4double pSTheta, G4double pDTheta);
94  //
95  // Constructs a sphere or sphere shell section
96  // with the given name and dimensions
97 
98  ~G4Sphere();
99  //
100  // Destructor
101 
102  // Accessors
103 
104  inline G4double GetInnerRadius () const;
105  inline G4double GetOuterRadius () const;
106  inline G4double GetStartPhiAngle () const;
107  inline G4double GetDeltaPhiAngle () const;
108  inline G4double GetStartThetaAngle() const;
109  inline G4double GetDeltaThetaAngle() const;
110 
111  // Modifiers
112 
113  inline void SetInnerRadius (G4double newRMin);
114  inline void SetOuterRadius (G4double newRmax);
115  inline void SetStartPhiAngle (G4double newSphi, G4bool trig=true);
116  inline void SetDeltaPhiAngle (G4double newDphi);
117  inline void SetStartThetaAngle(G4double newSTheta);
118  inline void SetDeltaThetaAngle(G4double newDTheta);
119 
120  // Methods for solid
121 
122  inline G4double GetCubicVolume();
124 
126  const G4int n,
127  const G4VPhysicalVolume* pRep);
128 
129  G4bool CalculateExtent(const EAxis pAxis,
130  const G4VoxelLimits& pVoxelLimit,
131  const G4AffineTransform& pTransform,
132  G4double& pmin, G4double& pmax) const;
133 
134  EInside Inside(const G4ThreeVector& p) const;
135 
136  G4ThreeVector SurfaceNormal( const G4ThreeVector& p) const;
137 
139  const G4ThreeVector& v) const;
140 
141  G4double DistanceToIn(const G4ThreeVector& p) const;
142 
144  const G4ThreeVector& v,
145  const G4bool calcNorm=G4bool(false),
146  G4bool *validNorm=0,
147  G4ThreeVector *n=0) const;
148 
149  G4double DistanceToOut(const G4ThreeVector& p) const;
150 
152 
154 
155  G4VSolid* Clone() const;
156 
157  std::ostream& StreamInfo(std::ostream& os) const;
158 
159  // Visualisation functions
160 
161  G4VisExtent GetExtent () const;
162  void DescribeYourselfTo(G4VGraphicsScene& scene) const;
164 
165  public: // without description
166 
167  G4Sphere(__void__&);
168  //
169  // Fake default constructor for usage restricted to direct object
170  // persistency for clients requiring preallocation of memory for
171  // persistifiable objects.
172 
173  G4Sphere(const G4Sphere& rhs);
174  G4Sphere& operator=(const G4Sphere& rhs);
175  // Copy constructor and assignment operator.
176 
177  // Old access functions
178 
179  inline G4double GetRmin() const;
180  inline G4double GetRmax() const;
181  inline G4double GetSPhi() const;
182  inline G4double GetDPhi() const;
183  inline G4double GetSTheta() const;
184  inline G4double GetDTheta() const;
185  inline G4double GetInsideRadius() const;
186  inline void SetInsideRadius(G4double newRmin);
187 
188  private:
189 
191  CreateRotatedVertices(const G4AffineTransform& pTransform,
192  G4int& noPolygonVertices) const;
193  //
194  // Creates the List of transformed vertices in the format required
195  // for G4VSolid:: ClipCrossSection and ClipBetweenSections
196 
197  inline void Initialize();
198  //
199  // Reset relevant values to zero
200 
201  inline void CheckThetaAngles(G4double sTheta, G4double dTheta);
202  inline void CheckSPhiAngle(G4double sPhi);
203  inline void CheckDPhiAngle(G4double dPhi);
204  inline void CheckPhiAngles(G4double sPhi, G4double dPhi);
205  //
206  // Reset relevant flags and angle values
207 
208  inline void InitializePhiTrigonometry();
209  inline void InitializeThetaTrigonometry();
210  //
211  // Recompute relevant trigonometric values and cache them
212 
214  //
215  // Algorithm for SurfaceNormal() following the original
216  // specification for points not on the surface
217 
218  private:
219 
220  // Used by distanceToOut
221  //
223 
224  // used by normal
225  //
227 
230  //
231  // Radial and angular tolerances
232 
234  //
235  // Radial and angular dimensions
236 
239  //
240  // Cached trigonometric values for Phi angle
241 
244  //
245  // Cached trigonometric values for Theta angle
246 
248  //
249  // Flags for identification of section, shell or full sphere
250 
252  //
253  // Cached half tolerance values
254 };
255 
256 #include "G4Sphere.icc"
257 
258 #endif
259 
260 #endif
G4double ePhi
Definition: G4Sphere.hh:237
G4bool fFullSphere
Definition: G4Sphere.hh:247
void InitializeThetaTrigonometry()
G4double fDTheta
Definition: G4Sphere.hh:233
G4double fDPhi
Definition: G4Sphere.hh:233
G4double halfCarTolerance
Definition: G4Sphere.hh:251
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Sphere.cc:872
G4double fEpsilon
Definition: G4Sphere.hh:228
void SetStartPhiAngle(G4double newSphi, G4bool trig=true)
CLHEP::Hep3Vector G4ThreeVector
G4double tanETheta2
Definition: G4Sphere.hh:242
G4double fSTheta
Definition: G4Sphere.hh:233
G4double fRmin
Definition: G4Sphere.hh:233
G4double fRminTolerance
Definition: G4Sphere.hh:228
G4double hDPhi
Definition: G4Sphere.hh:237
void CheckSPhiAngle(G4double sPhi)
G4Polyhedron * CreatePolyhedron() const
Definition: G4Sphere.cc:3183
G4double GetInsideRadius() const
G4double fSPhi
Definition: G4Sphere.hh:233
void InitializePhiTrigonometry()
G4double GetDeltaPhiAngle() const
int G4int
Definition: G4Types.hh:78
G4Sphere(const G4String &pName, G4double pRmin, G4double pRmax, G4double pSPhi, G4double pDPhi, G4double pSTheta, G4double pDTheta)
Definition: G4Sphere.cc:90
void CheckDPhiAngle(G4double dPhi)
G4double cosSTheta
Definition: G4Sphere.hh:242
G4double sinSPhi
Definition: G4Sphere.hh:237
G4ThreeVector GetPointOnSurface() const
Definition: G4Sphere.cc:3043
G4double GetSPhi() const
~G4Sphere()
Definition: G4Sphere.cc:145
void SetDeltaThetaAngle(G4double newDTheta)
G4double cosHDPhiOT
Definition: G4Sphere.hh:237
void SetStartThetaAngle(G4double newSTheta)
G4double GetRmin() const
G4double GetDPhi() const
G4double sinEPhi
Definition: G4Sphere.hh:237
G4double GetStartThetaAngle() const
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4Sphere.cc:228
G4double GetDTheta() const
void CheckThetaAngles(G4double sTheta, G4double dTheta)
G4double fRmax
Definition: G4Sphere.hh:233
bool G4bool
Definition: G4Types.hh:79
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=G4bool(false), G4bool *validNorm=0, G4ThreeVector *n=0) const
Definition: G4Sphere.cc:1900
G4double tanSTheta
Definition: G4Sphere.hh:242
G4ThreeVectorList * CreateRotatedVertices(const G4AffineTransform &pTransform, G4int &noPolygonVertices) const
Definition: G4Sphere.cc:2874
EInside Inside(const G4ThreeVector &p) const
Definition: G4Sphere.cc:473
G4double cosSPhi
Definition: G4Sphere.hh:237
void CheckPhiAngles(G4double sPhi, G4double dPhi)
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4VSolid.hh:79
G4GeometryType GetEntityType() const
Definition: G4Sphere.cc:3001
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Sphere.cc:570
const G4int n
G4double sinCPhi
Definition: G4Sphere.hh:237
void SetOuterRadius(G4double newRmax)
void SetInsideRadius(G4double newRmin)
G4double GetInnerRadius() const
G4double cosETheta
Definition: G4Sphere.hh:242
G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector &p) const
Definition: G4Sphere.cc:699
G4VSolid * Clone() const
Definition: G4Sphere.cc:3010
G4Sphere & operator=(const G4Sphere &rhs)
Definition: G4Sphere.cc:179
G4double fRmaxTolerance
Definition: G4Sphere.hh:228
G4double tanSTheta2
Definition: G4Sphere.hh:242
G4double GetRmax() const
void SetDeltaPhiAngle(G4double newDphi)
G4double kRadTolerance
Definition: G4Sphere.hh:228
G4double tanETheta
Definition: G4Sphere.hh:242
G4double eTheta
Definition: G4Sphere.hh:242
EInside
Definition: geomdefs.hh:58
G4bool fFullPhiSphere
Definition: G4Sphere.hh:247
EAxis
Definition: geomdefs.hh:54
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Sphere.cc:3019
G4double sinETheta
Definition: G4Sphere.hh:242
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Sphere.cc:217
G4double cPhi
Definition: G4Sphere.hh:237
G4VisExtent GetExtent() const
Definition: G4Sphere.cc:3172
G4double GetSurfaceArea()
Definition: G4Sphere.cc:3125
G4double kAngTolerance
Definition: G4Sphere.hh:228
void SetInnerRadius(G4double newRMin)
G4double GetSTheta() const
G4double GetOuterRadius() const
G4double GetStartPhiAngle() const
G4double cosEPhi
Definition: G4Sphere.hh:237
void Initialize()
double G4double
Definition: G4Types.hh:76
G4double GetCubicVolume()
G4double sinSTheta
Definition: G4Sphere.hh:242
G4double GetDeltaThetaAngle() const
G4double halfAngTolerance
Definition: G4Sphere.hh:251
G4double cosHDPhiIT
Definition: G4Sphere.hh:237
G4bool fFullThetaSphere
Definition: G4Sphere.hh:247
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Sphere.cc:3178
G4double cosCPhi
Definition: G4Sphere.hh:237