Geant4-11
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// G4Sphere
27//
28// Class description:
29//
30// A G4Sphere is, in the general case, a section of a spherical shell,
31// between specified phi and theta angles
32//
33// The phi and theta segments are described by a starting angle,
34// and the +ve delta angle for the shape.
35// If the delta angle is >=2*pi, or >=pi the shape is treated as
36// continuous in phi or theta respectively.
37//
38// Theta must lie between 0-pi (incl).
39//
40// Member Data:
41//
42// fRmin inner radius
43// fRmax outer radius
44//
45// fSPhi starting angle of the segment in radians
46// fDPhi delta angle of the segment in radians
47//
48// fSTheta starting angle of the segment in radians
49// fDTheta delta angle of the segment in radians
50//
51//
52// Note:
53// Internally fSPhi & fDPhi are adjusted so that fDPhi<=2PI,
54// and fDPhi+fSPhi<=2PI. This enables simpler comparisons to be
55// made with (say) Phi of a point.
56
57// 28.3.94 P.Kent: old C++ code converted to tolerant geometry
58// 17.9.96 V.Grichine: final modifications to commit
59// --------------------------------------------------------------------
60#ifndef G4SPHERE_HH
61#define G4SPHERE_HH
62
63#include "G4GeomTypes.hh"
64
65#if defined(G4GEOM_USE_USOLIDS)
66#define G4GEOM_USE_USPHERE 1
67#endif
68
69#if defined(G4GEOM_USE_USPHERE)
70 #define G4USphere G4Sphere
71 #include "G4USphere.hh"
72#else
73
75#include "G4CSGSolid.hh"
76#include "G4Polyhedron.hh"
77
78class G4VisExtent;
79
80class G4Sphere : public G4CSGSolid
81{
82 public: // with description
83
84 G4Sphere(const G4String& pName,
85 G4double pRmin, G4double pRmax,
86 G4double pSPhi, G4double pDPhi,
87 G4double pSTheta, G4double pDTheta);
88 //
89 // Constructs a sphere or sphere shell section
90 // with the given name and dimensions
91
92 ~G4Sphere();
93 //
94 // Destructor
95
96 // Accessors
97
98 inline G4double GetInnerRadius () const;
99 inline G4double GetOuterRadius () const;
100 inline G4double GetStartPhiAngle () const;
101 inline G4double GetDeltaPhiAngle () const;
104 inline G4double GetSinStartPhi () const;
105 inline G4double GetCosStartPhi () const;
106 inline G4double GetSinEndPhi () const;
107 inline G4double GetCosEndPhi () const;
108 inline G4double GetSinStartTheta () const;
109 inline G4double GetCosStartTheta () const;
110 inline G4double GetSinEndTheta () const;
111 inline G4double GetCosEndTheta () const;
112
113 // Modifiers
114
115 inline void SetInnerRadius (G4double newRMin);
116 inline void SetOuterRadius (G4double newRmax);
117 inline void SetStartPhiAngle (G4double newSphi, G4bool trig = true);
118 inline void SetDeltaPhiAngle (G4double newDphi);
119 inline void SetStartThetaAngle(G4double newSTheta);
120 inline void SetDeltaThetaAngle(G4double newDTheta);
121
122 // Methods for solid
123
126
128 const G4int n,
129 const G4VPhysicalVolume* pRep);
130
132
133 G4bool CalculateExtent(const EAxis pAxis,
134 const G4VoxelLimits& pVoxelLimit,
135 const G4AffineTransform& pTransform,
136 G4double& pmin, G4double& pmax) const;
137
138 EInside Inside(const G4ThreeVector& p) const;
139
141
143 const G4ThreeVector& v) const;
144
145 G4double DistanceToIn(const G4ThreeVector& p) const;
146
148 const G4ThreeVector& v,
149 const G4bool calcNorm = false,
150 G4bool* validNorm = nullptr,
151 G4ThreeVector* n = nullptr) const;
152
153 G4double DistanceToOut(const G4ThreeVector& p) const;
154
156
158
159 G4VSolid* Clone() const;
160
161 std::ostream& StreamInfo(std::ostream& os) const;
162
163 // Visualisation functions
164
165 G4VisExtent GetExtent () const;
166 void DescribeYourselfTo(G4VGraphicsScene& scene) const;
168
169 public: // without description
170
171 G4Sphere(__void__&);
172 //
173 // Fake default constructor for usage restricted to direct object
174 // persistency for clients requiring preallocation of memory for
175 // persistifiable objects.
176
177 G4Sphere(const G4Sphere& rhs);
178 G4Sphere& operator=(const G4Sphere& rhs);
179 // Copy constructor and assignment operator.
180
181
182 private:
183
184 inline void Initialize();
185 //
186 // Reset relevant values to zero
187
188 inline void CheckThetaAngles(G4double sTheta, G4double dTheta);
189 inline void CheckSPhiAngle(G4double sPhi);
190 inline void CheckDPhiAngle(G4double dPhi);
191 inline void CheckPhiAngles(G4double sPhi, G4double dPhi);
192 //
193 // Reset relevant flags and angle values
194
197 //
198 // Recompute relevant trigonometric values and cache them
199
201 //
202 // Algorithm for SurfaceNormal() following the original
203 // specification for points not on the surface
204
205 private:
206
207 // Used by distanceToOut
208 //
210
211 // used by normal
212 //
214
217 //
218 // Radial and angular tolerances
219
221 //
222 // Radial and angular dimensions
223
226 //
227 // Cached trigonometric values for Phi angle
228
231 //
232 // Cached trigonometric values for Theta angle
233
235 //
236 // Flags for identification of section, shell or full sphere
237
239 //
240 // Cached half tolerance values
241};
242
243#include "G4Sphere.icc"
244
245#endif
246
247#endif
static const G4double pMax
static const G4double pMin
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4ThreeVector GetPointOnSurface() const
Definition: G4Sphere.cc:2808
void InitializePhiTrigonometry()
G4double fRmin
Definition: G4Sphere.hh:220
~G4Sphere()
Definition: G4Sphere.cc:126
G4double cosHDPhiIT
Definition: G4Sphere.hh:224
G4double cPhi
Definition: G4Sphere.hh:225
G4double GetStartPhiAngle() const
G4double tanSTheta2
Definition: G4Sphere.hh:230
void SetDeltaPhiAngle(G4double newDphi)
void Initialize()
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4Sphere.cc:210
G4bool fFullSphere
Definition: G4Sphere.hh:234
G4Sphere(const G4String &pName, G4double pRmin, G4double pRmax, G4double pSPhi, G4double pDPhi, G4double pSTheta, G4double pDTheta)
Definition: G4Sphere.cc:72
G4double tanETheta2
Definition: G4Sphere.hh:230
G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector &p) const
Definition: G4Sphere.cc:531
G4VSolid * Clone() const
Definition: G4Sphere.cc:2742
G4double GetSinStartTheta() const
G4double cosSTheta
Definition: G4Sphere.hh:229
G4double GetCosStartPhi() const
void SetStartThetaAngle(G4double newSTheta)
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Sphere.cc:2751
G4double fEpsilon
Definition: G4Sphere.hh:216
G4bool fFullPhiSphere
Definition: G4Sphere.hh:234
void CheckPhiAngles(G4double sPhi, G4double dPhi)
G4double cosHDPhiOT
Definition: G4Sphere.hh:224
G4double fSPhi
Definition: G4Sphere.hh:220
G4double sinSTheta
Definition: G4Sphere.hh:229
G4double halfAngTolerance
Definition: G4Sphere.hh:238
G4double GetDeltaPhiAngle() const
G4VisExtent GetExtent() const
Definition: G4Sphere.cc:2859
G4double fDPhi
Definition: G4Sphere.hh:220
G4double GetCosEndTheta() const
G4double halfCarTolerance
Definition: G4Sphere.hh:238
void SetOuterRadius(G4double newRmax)
G4double fSTheta
Definition: G4Sphere.hh:220
G4double fRmaxTolerance
Definition: G4Sphere.hh:215
G4GeometryType GetEntityType() const
Definition: G4Sphere.cc:2733
G4double hDPhi
Definition: G4Sphere.hh:225
void SetDeltaThetaAngle(G4double newDTheta)
G4double kAngTolerance
Definition: G4Sphere.hh:215
void SetInnerRadius(G4double newRMin)
G4Polyhedron * CreatePolyhedron() const
Definition: G4Sphere.cc:2870
G4double cosHDPhi
Definition: G4Sphere.hh:224
G4double cosEPhi
Definition: G4Sphere.hh:225
G4double GetInnerRadius() const
EInside Inside(const G4ThreeVector &p) const
Definition: G4Sphere.cc:288
G4double GetOuterRadius() const
G4double fRmax
Definition: G4Sphere.hh:220
G4double tanETheta
Definition: G4Sphere.hh:230
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4Sphere.cc:267
G4double sinEPhi
Definition: G4Sphere.hh:225
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Sphere.cc:402
G4double GetCosEndPhi() const
void CheckSPhiAngle(G4double sPhi)
G4double GetSinEndTheta() const
G4double GetCubicVolume()
Definition: G4Sphere.cc:2775
void CheckDPhiAngle(G4double dPhi)
G4double GetDeltaThetaAngle() const
G4double cosSPhi
Definition: G4Sphere.hh:225
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Sphere.cc:704
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Sphere.cc:199
G4double GetSinEndPhi() const
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Sphere.cc:2865
G4double GetSinStartPhi() const
@ kNSTheta
Definition: G4Sphere.hh:213
@ kNETheta
Definition: G4Sphere.hh:213
G4Sphere & operator=(const G4Sphere &rhs)
Definition: G4Sphere.cc:160
G4double GetSurfaceArea()
Definition: G4Sphere.cc:2790
G4double GetStartThetaAngle() const
G4double GetCosStartTheta() const
G4double sinCPhi
Definition: G4Sphere.hh:224
G4bool fFullThetaSphere
Definition: G4Sphere.hh:234
G4double tanSTheta
Definition: G4Sphere.hh:230
G4double kRadTolerance
Definition: G4Sphere.hh:216
G4double ePhi
Definition: G4Sphere.hh:225
G4double cosCPhi
Definition: G4Sphere.hh:224
G4double cosETheta
Definition: G4Sphere.hh:229
void CheckThetaAngles(G4double sTheta, G4double dTheta)
G4double fDTheta
Definition: G4Sphere.hh:220
G4double sinETheta
Definition: G4Sphere.hh:229
G4double sinSPhi
Definition: G4Sphere.hh:225
void InitializeThetaTrigonometry()
void SetStartPhiAngle(G4double newSphi, G4bool trig=true)
G4double eTheta
Definition: G4Sphere.hh:230
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const
Definition: G4Sphere.cc:1749
G4double fRminTolerance
Definition: G4Sphere.hh:215
EAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:67