Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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$
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
73#include "G4CSGSolid.hh"
74
75class G4VisExtent;
76
77class G4Sphere : public G4CSGSolid
78{
79 public: // with description
80
81 G4Sphere(const G4String& pName,
82 G4double pRmin, G4double pRmax,
83 G4double pSPhi, G4double pDPhi,
84 G4double pSTheta, G4double pDTheta);
85 //
86 // Constructs a sphere or sphere shell section
87 // with the given name and dimensions
88
89 ~G4Sphere();
90 //
91 // Destructor
92
93 // Accessors
94
95 inline G4double GetInnerRadius () const;
96 inline G4double GetOuterRadius () const;
97 inline G4double GetStartPhiAngle () const;
98 inline G4double GetDeltaPhiAngle () const;
101
102 // Modifiers
103
104 inline void SetInnerRadius (G4double newRMin);
105 inline void SetOuterRadius (G4double newRmax);
106 inline void SetStartPhiAngle (G4double newSphi, G4bool trig=true);
107 inline void SetDeltaPhiAngle (G4double newDphi);
108 inline void SetStartThetaAngle(G4double newSTheta);
109 inline void SetDeltaThetaAngle(G4double newDTheta);
110
111 // Methods for solid
112
115
117 const G4int n,
118 const G4VPhysicalVolume* pRep);
119
120 G4bool CalculateExtent(const EAxis pAxis,
121 const G4VoxelLimits& pVoxelLimit,
122 const G4AffineTransform& pTransform,
123 G4double& pmin, G4double& pmax) const;
124
125 EInside Inside(const G4ThreeVector& p) const;
126
128
130 const G4ThreeVector& v) const;
131
132 G4double DistanceToIn(const G4ThreeVector& p) const;
133
135 const G4ThreeVector& v,
136 const G4bool calcNorm=G4bool(false),
137 G4bool *validNorm=0,
138 G4ThreeVector *n=0) const;
139
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 G4VisExtent GetExtent () const;
153 void DescribeYourselfTo(G4VGraphicsScene& scene) const;
155 G4NURBS* CreateNURBS() const;
156
157 public: // without description
158
159 G4Sphere(__void__&);
160 //
161 // Fake default constructor for usage restricted to direct object
162 // persistency for clients requiring preallocation of memory for
163 // persistifiable objects.
164
165 G4Sphere(const G4Sphere& rhs);
166 G4Sphere& operator=(const G4Sphere& rhs);
167 // Copy constructor and assignment operator.
168
169 // Old access functions
170
171 inline G4double GetRmin() const;
172 inline G4double GetRmax() const;
173 inline G4double GetSPhi() const;
174 inline G4double GetDPhi() const;
175 inline G4double GetSTheta() const;
176 inline G4double GetDTheta() const;
177 inline G4double GetInsideRadius() const;
178 inline void SetInsideRadius(G4double newRmin);
179
180 private:
181
183 CreateRotatedVertices(const G4AffineTransform& pTransform,
184 G4int& noPolygonVertices) const;
185 //
186 // Creates the List of transformed vertices in the format required
187 // for G4VSolid:: ClipCrossSection and ClipBetweenSections
188
189 inline void Initialize();
190 //
191 // Reset relevant values to zero
192
193 inline void CheckThetaAngles(G4double sTheta, G4double dTheta);
194 inline void CheckSPhiAngle(G4double sPhi);
195 inline void CheckDPhiAngle(G4double dPhi);
196 inline void CheckPhiAngles(G4double sPhi, G4double dPhi);
197 //
198 // Reset relevant flags and angle values
199
200 inline void InitializePhiTrigonometry();
201 inline void InitializeThetaTrigonometry();
202 //
203 // Recompute relevant trigonometric values and cache them
204
205 G4ThreeVector ApproxSurfaceNormal( const G4ThreeVector& p) const;
206 //
207 // Algorithm for SurfaceNormal() following the original
208 // specification for points not on the surface
209
210 private:
211
212 // Used by distanceToOut
213 //
214 enum ESide {kNull,kRMin,kRMax,kSPhi,kEPhi,kSTheta,kETheta};
215
216 // used by normal
217 //
218 enum ENorm {kNRMin,kNRMax,kNSPhi,kNEPhi,kNSTheta,kNETheta};
219
220 G4double fRminTolerance, fRmaxTolerance, kAngTolerance,
221 kRadTolerance, fEpsilon;
222 //
223 // Radial and angular tolerances
224
225 G4double fRmin, fRmax, fSPhi, fDPhi, fSTheta, fDTheta;
226 //
227 // Radial and angular dimensions
228
229 G4double sinCPhi, cosCPhi, cosHDPhiOT, cosHDPhiIT,
230 sinSPhi, cosSPhi, sinEPhi, cosEPhi, hDPhi, cPhi, ePhi;
231 //
232 // Cached trigonometric values for Phi angle
233
234 G4double sinSTheta, cosSTheta, sinETheta, cosETheta,
235 tanSTheta, tanSTheta2, tanETheta, tanETheta2, eTheta;
236 //
237 // Cached trigonometric values for Theta angle
238
239 G4bool fFullPhiSphere, fFullThetaSphere, fFullSphere;
240 //
241 // Flags for identification of section, shell or full sphere
242};
243
244#include "G4Sphere.icc"
245
246#endif
ENorm
Definition: G4Cons.cc:72
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4VSolid.hh:85
G4ThreeVector GetPointOnSurface() const
Definition: G4Sphere.cc:3042
~G4Sphere()
Definition: G4Sphere.cc:141
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=G4bool(false), G4bool *validNorm=0, G4ThreeVector *n=0) const
Definition: G4Sphere.cc:1897
G4double GetStartPhiAngle() const
void SetDeltaPhiAngle(G4double newDphi)
G4VSolid * Clone() const
Definition: G4Sphere.cc:3009
G4NURBS * CreateNURBS() const
Definition: G4Sphere.cc:3187
void SetStartThetaAngle(G4double newSTheta)
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Sphere.cc:3018
G4double GetRmax() const
G4double GetSPhi() const
G4double GetDeltaPhiAngle() const
G4VisExtent GetExtent() const
Definition: G4Sphere.cc:3171
void SetOuterRadius(G4double newRmax)
G4double GetSTheta() const
G4GeometryType GetEntityType() const
Definition: G4Sphere.cc:3000
void SetDeltaThetaAngle(G4double newDTheta)
void SetInnerRadius(G4double newRMin)
G4Polyhedron * CreatePolyhedron() const
Definition: G4Sphere.cc:3182
G4double GetInnerRadius() const
EInside Inside(const G4ThreeVector &p) const
Definition: G4Sphere.cc:465
G4double GetOuterRadius() const
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4Sphere.cc:220
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Sphere.cc:562
void SetInsideRadius(G4double newRmin)
G4double GetCubicVolume()
G4double GetDeltaThetaAngle() const
G4double GetDPhi() const
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Sphere.cc:867
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Sphere.cc:209
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Sphere.cc:3177
G4Sphere & operator=(const G4Sphere &rhs)
Definition: G4Sphere.cc:173
G4double GetSurfaceArea()
Definition: G4Sphere.cc:3124
G4double GetStartThetaAngle() const
G4double GetDTheta() const
G4double GetRmin() const
void SetStartPhiAngle(G4double newSphi, G4bool trig=true)
G4double GetInsideRadius() const
EAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:58