Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 100820 2016-11-02 15:18:48Z 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 
82 #include "G4CSGSolid.hh"
83 #include "G4Polyhedron.hh"
84 
85 class G4VisExtent;
86 
87 class G4Sphere : public G4CSGSolid
88 {
89  public: // with description
90 
91  G4Sphere(const G4String& pName,
92  G4double pRmin, G4double pRmax,
93  G4double pSPhi, G4double pDPhi,
94  G4double pSTheta, G4double pDTheta);
95  //
96  // Constructs a sphere or sphere shell section
97  // with the given name and dimensions
98 
99  ~G4Sphere();
100  //
101  // Destructor
102 
103  // Accessors
104 
105  inline G4double GetInnerRadius () const;
106  inline G4double GetOuterRadius () const;
107  inline G4double GetStartPhiAngle () const;
108  inline G4double GetDeltaPhiAngle () const;
109  inline G4double GetStartThetaAngle() const;
110  inline G4double GetDeltaThetaAngle() const;
111  inline G4double GetSinStartPhi () const;
112  inline G4double GetCosStartPhi () const;
113  inline G4double GetSinEndPhi () const;
114  inline G4double GetCosEndPhi () const;
115  inline G4double GetSinStartTheta () const;
116  inline G4double GetCosStartTheta () const;
117  inline G4double GetSinEndTheta () const;
118  inline G4double GetCosEndTheta () const;
119 
120  // Modifiers
121 
122  inline void SetInnerRadius (G4double newRMin);
123  inline void SetOuterRadius (G4double newRmax);
124  inline void SetStartPhiAngle (G4double newSphi, G4bool trig=true);
125  inline void SetDeltaPhiAngle (G4double newDphi);
126  inline void SetStartThetaAngle(G4double newSTheta);
127  inline void SetDeltaThetaAngle(G4double newDTheta);
128 
129  // Methods for solid
130 
131  inline G4double GetCubicVolume();
133 
135  const G4int n,
136  const G4VPhysicalVolume* pRep);
137 
138  void Extent(G4ThreeVector& pMin, G4ThreeVector& pMax) const;
139 
140  G4bool CalculateExtent(const EAxis pAxis,
141  const G4VoxelLimits& pVoxelLimit,
142  const G4AffineTransform& pTransform,
143  G4double& pmin, G4double& pmax) const;
144 
145  EInside Inside(const G4ThreeVector& p) const;
146 
147  G4ThreeVector SurfaceNormal( const G4ThreeVector& p) const;
148 
150  const G4ThreeVector& v) const;
151 
152  G4double DistanceToIn(const G4ThreeVector& p) const;
153 
155  const G4ThreeVector& v,
156  const G4bool calcNorm=G4bool(false),
157  G4bool *validNorm=0,
158  G4ThreeVector *n=0) const;
159 
160  G4double DistanceToOut(const G4ThreeVector& p) const;
161 
163 
165 
166  G4VSolid* Clone() const;
167 
168  std::ostream& StreamInfo(std::ostream& os) const;
169 
170  // Visualisation functions
171 
172  G4VisExtent GetExtent () const;
173  void DescribeYourselfTo(G4VGraphicsScene& scene) const;
175 
176  public: // without description
177 
178  G4Sphere(__void__&);
179  //
180  // Fake default constructor for usage restricted to direct object
181  // persistency for clients requiring preallocation of memory for
182  // persistifiable objects.
183 
184  G4Sphere(const G4Sphere& rhs);
185  G4Sphere& operator=(const G4Sphere& rhs);
186  // Copy constructor and assignment operator.
187 
188  // Old access functions
189 
190  inline G4double GetRmin() const;
191  inline G4double GetRmax() const;
192  inline G4double GetSPhi() const;
193  inline G4double GetDPhi() const;
194  inline G4double GetSTheta() const;
195  inline G4double GetDTheta() const;
196  inline G4double GetInsideRadius() const;
197  inline void SetInsideRadius(G4double newRmin);
198 
199  private:
200 
201  inline void Initialize();
202  //
203  // Reset relevant values to zero
204 
205  inline void CheckThetaAngles(G4double sTheta, G4double dTheta);
206  inline void CheckSPhiAngle(G4double sPhi);
207  inline void CheckDPhiAngle(G4double dPhi);
208  inline void CheckPhiAngles(G4double sPhi, G4double dPhi);
209  //
210  // Reset relevant flags and angle values
211 
212  inline void InitializePhiTrigonometry();
213  inline void InitializeThetaTrigonometry();
214  //
215  // Recompute relevant trigonometric values and cache them
216 
217  G4ThreeVector ApproxSurfaceNormal( const G4ThreeVector& p) const;
218  //
219  // Algorithm for SurfaceNormal() following the original
220  // specification for points not on the surface
221 
222  private:
223 
224  // Used by distanceToOut
225  //
226  enum ESide {kNull,kRMin,kRMax,kSPhi,kEPhi,kSTheta,kETheta};
227 
228  // used by normal
229  //
230  enum ENorm {kNRMin,kNRMax,kNSPhi,kNEPhi,kNSTheta,kNETheta};
231 
232  G4double fRminTolerance, fRmaxTolerance, kAngTolerance,
233  kRadTolerance, fEpsilon;
234  //
235  // Radial and angular tolerances
236 
237  G4double fRmin, fRmax, fSPhi, fDPhi, fSTheta, fDTheta;
238  //
239  // Radial and angular dimensions
240 
241  G4double sinCPhi, cosCPhi, cosHDPhiOT, cosHDPhiIT,
242  sinSPhi, cosSPhi, sinEPhi, cosEPhi, hDPhi, cPhi, ePhi;
243  //
244  // Cached trigonometric values for Phi angle
245 
246  G4double sinSTheta, cosSTheta, sinETheta, cosETheta,
247  tanSTheta, tanSTheta2, tanETheta, tanETheta2, eTheta;
248  //
249  // Cached trigonometric values for Theta angle
250 
251  G4bool fFullPhiSphere, fFullThetaSphere, fFullSphere;
252  //
253  // Flags for identification of section, shell or full sphere
254 
255  G4double halfCarTolerance, halfAngTolerance;
256  //
257  // Cached half tolerance values
258 };
259 
260 #include "G4Sphere.icc"
261 
262 #endif
263 
264 #endif
G4double GetCosStartTheta() const
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Sphere.cc:726
void SetStartPhiAngle(G4double newSphi, G4bool trig=true)
G4double GetSinStartPhi() const
G4double GetCosEndTheta() const
const char * p
Definition: xmltok.h:285
G4double GetSinEndTheta() const
G4Polyhedron * CreatePolyhedron() const
Definition: G4Sphere.cc:2937
G4double GetInsideRadius() const
G4double GetSinStartTheta() const
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:95
G4ThreeVector GetPointOnSurface() const
Definition: G4Sphere.cc:2797
G4double GetSPhi() const
~G4Sphere()
Definition: G4Sphere.cc:150
void SetDeltaThetaAngle(G4double newDTheta)
void SetStartThetaAngle(G4double newSTheta)
G4double GetRmin() const
G4double GetDPhi() const
G4double GetStartThetaAngle() const
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4Sphere.cc:289
G4double GetDTheta() const
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:1771
EInside Inside(const G4ThreeVector &p) const
Definition: G4Sphere.cc:310
G4GeometryType GetEntityType() const
Definition: G4Sphere.cc:2755
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Sphere.cc:424
void SetOuterRadius(G4double newRmax)
void SetInsideRadius(G4double newRmin)
G4double GetInnerRadius() const
G4VSolid * Clone() const
Definition: G4Sphere.cc:2764
G4Sphere & operator=(const G4Sphere &rhs)
Definition: G4Sphere.cc:184
G4double GetRmax() const
void SetDeltaPhiAngle(G4double newDphi)
EInside
Definition: geomdefs.hh:58
EAxis
Definition: geomdefs.hh:54
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Sphere.cc:2773
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Sphere.cc:222
G4double GetCosStartPhi() const
G4VisExtent GetExtent() const
Definition: G4Sphere.cc:2926
G4double GetSurfaceArea()
Definition: G4Sphere.cc:2879
void SetInnerRadius(G4double newRMin)
G4double GetSTheta() const
G4double GetSinEndPhi() const
void Extent(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4Sphere.cc:233
G4double GetOuterRadius() const
G4double GetStartPhiAngle() const
double G4double
Definition: G4Types.hh:76
G4double GetCubicVolume()
G4double GetDeltaThetaAngle() const
G4double GetCosEndPhi() const
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Sphere.cc:2932