source: trunk/source/geometry/solids/CSG/include/G4Sphere.hh@ 1232

Last change on this file since 1232 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

File size: 8.1 KB
Line 
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,v 1.24 2009/03/31 07:51:49 gcosmo Exp $
28// GEANT4 tag $Name: geant4-09-03 $
29//
30//
31// --------------------------------------------------------------------
32// GEANT 4 class header file
33//
34// G4Sphere
35//
36// Class description:
37//
38// A G4Sphere is, in the general case, a section of a spherical shell,
39// between specified phi and theta angles
40//
41// The phi and theta segments are described by a starting angle,
42// and the +ve delta angle for the shape.
43// If the delta angle is >=2*pi, or >=pi the shape is treated as
44// continuous in phi or theta respectively.
45//
46// Theta must lie between 0-pi (incl).
47//
48// Member Data:
49//
50// fRmin inner radius
51// fRmax outer radius
52//
53// fSPhi starting angle of the segment in radians
54// fDPhi delta angle of the segment in radians
55//
56// fSTheta starting angle of the segment in radians
57// fDTheta delta angle of the segment in radians
58//
59//
60// Note:
61// Internally fSPhi & fDPhi are adjusted so that fDPhi<=2PI,
62// and fDPhi+fSPhi<=2PI. This enables simpler comparisons to be
63// made with (say) Phi of a point.
64
65// History:
66// 28.3.94 P.Kent: old C++ code converted to tolerant geometry
67// 17.9.96 V.Grichine: final modifications to commit
68// --------------------------------------------------------------------
69
70#ifndef G4Sphere_HH
71#define G4Sphere_HH
72
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;
99 inline G4double GetStartThetaAngle() const;
100 inline G4double GetDeltaThetaAngle() 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
113 inline G4double GetCubicVolume();
114 G4double GetSurfaceArea();
115
116 void ComputeDimensions( G4VPVParameterisation* p,
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
127 G4ThreeVector SurfaceNormal( const G4ThreeVector& p) const;
128
129 G4double DistanceToIn(const G4ThreeVector& p,
130 const G4ThreeVector& v) const;
131
132 G4double DistanceToIn(const G4ThreeVector& p) const;
133
134 G4double DistanceToOut(const G4ThreeVector& p,
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
142 G4GeometryType GetEntityType() const;
143
144 G4ThreeVector GetPointOnSurface() const;
145
146 std::ostream& StreamInfo(std::ostream& os) const;
147
148 // Visualisation functions
149
150 G4VisExtent GetExtent () const;
151 void DescribeYourselfTo(G4VGraphicsScene& scene) const;
152 G4Polyhedron* CreatePolyhedron() const;
153 G4NURBS* CreateNURBS() const;
154
155 public: // without description
156
157 G4Sphere(__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 // Old access functions
164
165 inline G4double GetRmin() const;
166 inline G4double GetRmax() const;
167 inline G4double GetSPhi() const;
168 inline G4double GetDPhi() const;
169 inline G4double GetSTheta() const;
170 inline G4double GetDTheta() const;
171 inline G4double GetInsideRadius() const;
172 inline void SetInsideRadius(G4double newRmin);
173
174 private:
175
176 G4ThreeVectorList*
177 CreateRotatedVertices(const G4AffineTransform& pTransform,
178 G4int& noPolygonVertices) 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 CheckThetaAngles(G4double sTheta, G4double dTheta);
188 inline void CheckSPhiAngle(G4double sPhi);
189 inline void CheckDPhiAngle(G4double dPhi);
190 inline void CheckPhiAngles(G4double sPhi, G4double dPhi);
191 //
192 // Reset relevant flags and angle values
193
194 inline void InitializePhiTrigonometry();
195 inline void InitializeThetaTrigonometry();
196 //
197 // Recompute relevant trigonometric values and cache them
198
199 G4ThreeVector ApproxSurfaceNormal( const G4ThreeVector& p) const;
200 //
201 // Algorithm for SurfaceNormal() following the original
202 // specification for points not on the surface
203
204 private:
205
206 // Used by distanceToOut
207 //
208 enum ESide {kNull,kRMin,kRMax,kSPhi,kEPhi,kSTheta,kETheta};
209
210 // used by normal
211 //
212 enum ENorm {kNRMin,kNRMax,kNSPhi,kNEPhi,kNSTheta,kNETheta};
213
214 G4double fEpsilon, fRminTolerance, fRmaxTolerance, kAngTolerance;
215 //
216 // Radial and angular tolerances
217
218 G4double fRmin, fRmax, fSPhi, fDPhi, fSTheta, fDTheta;
219 //
220 // Radial and angular dimensions
221
222 G4double sinCPhi, cosCPhi, cosHDPhiOT, cosHDPhiIT,
223 sinSPhi, cosSPhi, sinEPhi, cosEPhi, hDPhi, cPhi, ePhi;
224 //
225 // Cached trigonometric values for Phi angle
226
227 G4double sinSTheta, cosSTheta, sinETheta, cosETheta,
228 tanSTheta, tanSTheta2, tanETheta, tanETheta2, eTheta;
229 //
230 // Cached trigonometric values for Theta angle
231
232 G4bool fFullPhiSphere, fFullThetaSphere, fFullSphere;
233 //
234 // Flags for identification of section, shell or full sphere
235};
236
237#include "G4Sphere.icc"
238
239#endif
Note: See TracBrowser for help on using the repository browser.