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

Last change on this file since 1086 was 1058, checked in by garnier, 17 years ago

file release beta

File size: 6.8 KB
RevLine 
[831]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//
[921]27// $Id: G4Sphere.hh,v 1.21 2008/11/21 09:50:05 gcosmo Exp $
[1058]28// GEANT4 tag $Name: geant4-09-02-ref-02 $
[831]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, 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 virtual ~G4Sphere() ;
87
88 // Accessors
89
[921]90 inline G4double GetInnerRadius () const;
[831]91 inline G4double GetOuterRadius () const;
92 inline G4double GetStartPhiAngle () const;
93 inline G4double GetDeltaPhiAngle () const;
94 inline G4double GetStartThetaAngle() const;
95 inline G4double GetDeltaThetaAngle() const;
96
97 // Modifiers
98
[921]99 inline void SetInnerRadius (G4double newRMin);
[831]100 inline void SetOuterRadius (G4double newRmax);
101 inline void SetStartPhiAngle (G4double newSphi);
102 inline void SetDeltaPhiAngle (G4double newDphi);
103 inline void SetStartThetaAngle(G4double newSTheta);
104 inline void SetDeltaThetaAngle(G4double newDTheta);
105
106 // Methods for solid
107
108 inline G4double GetCubicVolume();
109 inline G4double GetSurfaceArea();
110
111 void ComputeDimensions( G4VPVParameterisation* p,
112 const G4int n,
113 const G4VPhysicalVolume* pRep);
114
115 G4bool CalculateExtent(const EAxis pAxis,
116 const G4VoxelLimits& pVoxelLimit,
117 const G4AffineTransform& pTransform,
118 G4double& pmin, G4double& pmax) const;
119
120 EInside Inside(const G4ThreeVector& p) const;
121
122 G4ThreeVector SurfaceNormal( const G4ThreeVector& p) const;
123
124 G4double DistanceToIn(const G4ThreeVector& p,
125 const G4ThreeVector& v) const;
126
127 G4double DistanceToIn(const G4ThreeVector& p) const;
128
129 G4double DistanceToOut(const G4ThreeVector& p,
130 const G4ThreeVector& v,
131 const G4bool calcNorm=G4bool(false),
132 G4bool *validNorm=0,
133 G4ThreeVector *n=0) const;
134
135 G4double DistanceToOut(const G4ThreeVector& p) const;
136
137 G4GeometryType GetEntityType() const;
138
139 G4ThreeVector GetPointOnSurface() const;
140
141 std::ostream& StreamInfo(std::ostream& os) const;
142
143 // Visualisation functions
144 G4VisExtent GetExtent () const;
145 void DescribeYourselfTo(G4VGraphicsScene& scene) const;
146 G4Polyhedron* CreatePolyhedron() const;
147 G4NURBS* CreateNURBS() const;
148
149 public: // without description
150
151 G4Sphere(__void__&);
152 // Fake default constructor for usage restricted to direct object
153 // persistency for clients requiring preallocation of memory for
154 // persistifiable objects.
155
156 // Old access functions
157
158 inline G4double GetRmin() const;
159 inline G4double GetRmax() const;
160 inline G4double GetSPhi() const;
161 inline G4double GetDPhi() const;
162 inline G4double GetSTheta() const;
163 inline G4double GetDTheta() const;
[921]164 inline G4double GetInsideRadius() const;
165 inline void SetInsideRadius(G4double newRmin);
[831]166
167 protected:
168
169 G4ThreeVectorList*
170 CreateRotatedVertices(const G4AffineTransform& pTransform,
171 G4int& noPolygonVertices) const;
172
173 // Used by distanceToOut
174
175 enum ESide {kNull,kRMin,kRMax,kSPhi,kEPhi,kSTheta,kETheta};
176
177 // used by normal
178
179 enum ENorm {kNRMin,kNRMax,kNSPhi,kNEPhi,kNSTheta,kNETheta};
180
181 private:
182
183 G4ThreeVector ApproxSurfaceNormal( const G4ThreeVector& p) const;
184 // Algorithm for SurfaceNormal() following the original
185 // specification for points not on the surface
186
187 private:
188
189 G4double kAngTolerance, kRadTolerance;
190
191 G4double fRmin,fRmax,
192 fSPhi,fDPhi,
193 fSTheta,fDTheta;
194 G4double fEpsilon;
195};
196
197#include "G4Sphere.icc"
198
199#endif
Note: See TracBrowser for help on using the repository browser.