Geant4  10.02
G4VSolid.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: G4VSolid.hh 92695 2015-09-14 08:49:44Z gcosmo $
28 //
29 //
30 // class G4VSolid
31 //
32 // Class description:
33 //
34 // Abstract base class for solids, physical shapes that can be tracked through.
35 // Each solid has a name, and the constructors and destructors automatically
36 // add and subtract them from the G4SolidStore, a singleton `master' List
37 // of available solids.
38 //
39 // This class defines, but does not implement, functions to compute
40 // distances to/from the shape. Functions are also defined
41 // to check whether a point is inside the shape, to return the
42 // surface normal of the shape at a given point, and to compute
43 // the extent of the shape. [see descriptions below]
44 //
45 // Some protected/private utility functions are implemented for the
46 // clipping of regions for the computation of a solid's extent. Note that
47 // the clipping mechanism is presently inefficient.
48 //
49 // Some visualization/graphics functions are also defined.
50 //
51 // Member Data:
52 //
53 // G4String fshapeName
54 // - Name for this solid.
55 
56 // History:
57 // 12.04.00 J.Allison Implemented GetExtent() in terms of CalculateExtent()
58 // 17.06.98 J.Apostolakis Added pure virtual function GetEntityType()
59 // 26.07.96 P.Kent Added ComputeDimensions for replication mechanism
60 // 27.03.96 J.Allison Methods for visualisation
61 // 30.06.95 P.Kent Initial version, no scoping or visualisation functions
62 // --------------------------------------------------------------------
63 #ifndef G4VSOLID_HH
64 #define G4VSOLID_HH
65 
66 #include "G4Types.hh"
67 #include "G4String.hh"
68 #include "geomdefs.hh"
69 
70 class G4AffineTransform;
71 class G4VoxelLimits;
72 
74 class G4VPhysicalVolume;
75 
76 class G4VGraphicsScene;
77 class G4Polyhedron;
78 class G4VisExtent;
80 
81 #include "G4ThreeVector.hh"
82 #include <vector>
83 
84 typedef std::vector<G4ThreeVector> G4ThreeVectorList;
86 
87 class G4VSolid
88 {
89  public: // with description
90 
91  G4VSolid(const G4String& name);
92  // Creates a new shape, with the supplied name. No provision is made
93  // for sharing a common name amongst multiple classes.
94  virtual ~G4VSolid();
95  // Default destructor.
96 
97  inline G4bool operator==( const G4VSolid& s ) const;
98  // Return true only if addresses are the same.
99 
100  friend std::ostream& operator<< ( std::ostream& os, const G4VSolid& e );
101  // Streaming operator, using DumpInfo().
102 
103  inline G4String GetName() const;
104  // Returns the current shape's name.
105  inline void SetName(const G4String& name);
106  // Sets the current shape's name.
107 
108  inline G4double GetTolerance() const;
109  // Returns the cached geometrical tolerance.
110 
111  virtual G4bool CalculateExtent(const EAxis pAxis,
112  const G4VoxelLimits& pVoxelLimit,
113  const G4AffineTransform& pTransform,
114  G4double& pMin, G4double& pMax) const = 0;
115  // Calculate the minimum and maximum extent of the solid, when under the
116  // specified transform, and within the specified limits. If the solid
117  // is not intersected by the region, return false, else return true.
118 
119  virtual EInside Inside(const G4ThreeVector& p) const = 0;
120  // Returns kOutside if the point at offset p is outside the shapes
121  // boundaries plus Tolerance/2, kSurface if the point is <= Tolerance/2
122  // from a surface, otherwise kInside.
123 
124  virtual G4ThreeVector SurfaceNormal(const G4ThreeVector& p) const = 0;
125  // Returns the outwards pointing unit normal of the shape for the
126  // surface closest to the point at offset p.
127 
128  virtual G4double DistanceToIn(const G4ThreeVector& p,
129  const G4ThreeVector& v) const = 0;
130  // Return the distance along the normalised vector v to the shape,
131  // from the point at offset p. If there is no intersection, return
132  // kInfinity. The first intersection resulting from `leaving' a
133  // surface/volume is discarded. Hence, it is tolerant of points on
134  // the surface of the shape.
135 
136  virtual G4double DistanceToIn(const G4ThreeVector& p) const = 0;
137  // Calculate the distance to the nearest surface of a shape from an
138  // outside point. The distance can be an underestimate.
139 
140  virtual G4double DistanceToOut(const G4ThreeVector& p,
141  const G4ThreeVector& v,
142  const G4bool calcNorm=false,
143  G4bool *validNorm=0,
144  G4ThreeVector *n=0) const = 0;
145  // Return the distance along the normalised vector v to the shape,
146  // from a point at an offset p inside or on the surface of the shape.
147  // Intersections with surfaces, when the point is < Tolerance/2 from a
148  // surface must be ignored.
149  // If calcNorm==true:
150  // validNorm set true if the solid lies entirely behind or on the
151  // exiting surface.
152  // n set to exiting outwards normal vector (undefined Magnitude).
153  // validNorm set to false if the solid does not lie entirely behind
154  // or on the exiting surface
155  // If calcNorm==false:
156  // validNorm and n are unused.
157  //
158  // Must be called as solid.DistanceToOut(p,v) or by specifying all
159  // the parameters.
160 
161  virtual G4double DistanceToOut(const G4ThreeVector& p) const = 0;
162  // Calculate the distance to the nearest surface of a shape from an
163  // inside point. The distance can be an underestimate.
164 
165 
167  const G4int n,
168  const G4VPhysicalVolume* pRep);
169  // Throw exception if ComputeDimensions called frrom an illegal
170  // derived class.
171 
172  virtual G4double GetCubicVolume();
173  // Returns an estimation of the solid volume in internal units.
174  // This method may be overloaded by derived classes to compute the
175  // exact geometrical quantity for solids where this is possible,
176  // or anyway to cache the computed value.
177  // Note: the computed value is NOT cached.
178 
179  virtual G4double GetSurfaceArea();
180  // Return an estimation of the solid surface area in internal units.
181  // This method may be overloaded by derived classes to compute the
182  // exact geometrical quantity for solids where this is possible,
183  // or anyway to cache the computed value.
184  // Note: the computed value is NOT cached.
185 
186  virtual G4GeometryType GetEntityType() const = 0;
187  // Provide identification of the class of an object.
188  // (required for persistency and STEP interface)
189 
190  virtual G4ThreeVector GetPointOnSurface() const;
191  // Returns a random point located on the surface of the solid.
192  // Points returned are not necessarily uniformly distributed.
193 
194  virtual G4VSolid* Clone() const;
195  // Returns a pointer of a dynamically allocated copy of the solid.
196  // Returns NULL pointer with warning in case the concrete solid does not
197  // implement this method. The caller has responsibility for ownership.
198 
199  virtual std::ostream& StreamInfo(std::ostream& os) const = 0;
200  // Dumps contents of the solid to a stream.
201  inline void DumpInfo() const;
202  // Dumps contents of the solid to the standard output.
203 
204  // Visualization functions
205 
206  virtual void DescribeYourselfTo (G4VGraphicsScene& scene) const = 0;
207  // A "double dispatch" function which identifies the solid
208  // to the graphics scene.
209  virtual G4VisExtent GetExtent () const;
210  // Provide extent (bounding box) as possible hint to the graphics view.
211  virtual G4Polyhedron* CreatePolyhedron () const;
212  // Create a G4Polyhedron. (It is the caller's responsibility
213  // to delete it). A null pointer means "not created".
214  virtual G4Polyhedron* GetPolyhedron () const;
215  // Smart access function - creates on request and stores for future
216  // access. A null pointer means "not available".
217 
218  virtual const G4VSolid* GetConstituentSolid(G4int no) const;
219  virtual G4VSolid* GetConstituentSolid(G4int no);
220  // If the solid is made up from a Boolean operation of two solids,
221  // return the "no" solid. If the solid is not a "Boolean", return 0.
222 
223  virtual const G4DisplacedSolid* GetDisplacedSolidPtr() const;
225  // If the solid is a "G4DisplacedSolid", return a self pointer
226  // else return 0.
227 
228  public: // without description
229 
230  G4VSolid(__void__&);
231  // Fake default constructor for usage restricted to direct object
232  // persistency for clients requiring preallocation of memory for
233  // persistifiable objects.
234 
235  G4VSolid(const G4VSolid& rhs);
236  G4VSolid& operator=(const G4VSolid& rhs);
237  // Copy constructor and assignment operator.
238 
240  // Calculate cubic volume based on Inside() method.
241  // Accuracy is limited by the second argument or the statistics
242  // expressed by the first argument.
243 
244  G4double EstimateSurfaceArea(G4int nStat, G4double ell) const;
245  // Calculate surface area only based on Inside() method.
246  // Accuracy is limited by the second argument or the statistics
247  // expressed by the first argument.
248 
249  protected: // with description
250 
252  const G4VoxelLimits& pVoxelLimit,
253  const EAxis pAxis,
254  G4double& pMin, G4double& pMax) const;
255  // Calculate the maximum and minimum extents of the convex polygon
256  // pPolygon along the axis pAxis, within the limits pVoxelLimit.
257  //
258  // If the minimum is <pMin pMin is set to the new minimum.
259  // If the maximum is >pMax pMax is set to the new maximum.
260  //
261  // Modifications to pPolygon are made - it is left in an undefined state.
262 
263  void ClipCrossSection(G4ThreeVectorList* pVertices,
264  const G4int pSectionIndex,
265  const G4VoxelLimits& pVoxelLimit,
266  const EAxis pAxis,
267  G4double& pMin, G4double& pMax) const;
268  // Calculate the maximum and minimum extents of the polygon described
269  // by the vertices: pSectionIndex->pSectionIndex+1->
270  // pSectionIndex+2->pSectionIndex+3->pSectionIndex
271  // in the List pVertices.
272  //
273  // If the minimum is <pMin pMin is set to the new minimum.
274  // If the maximum is >pMax pMax is set to the new maximum.
275  //
276  // No modifications are made to pVertices.
277 
278  void ClipBetweenSections(G4ThreeVectorList* pVertices,
279  const G4int pSectionIndex,
280  const G4VoxelLimits& pVoxelLimit,
281  const EAxis pAxis,
282  G4double& pMin, G4double& pMax) const;
283  // Calculate the maximum and minimum extents of the polygons
284  // joining the CrossSections at pSectionIndex->pSectionIndex+3 and
285  // pSectionIndex+4->pSectionIndex7
286  // in the List pVertices, within the boundaries of the voxel limits
287  // pVoxelLimit.
288  //
289  // If the minimum is <pMin pMin is set to the new minimum.
290  // If the maximum is >pMax pMax is set to the new maximum.
291  //
292  // No modifications are made to pVertices.
293 
294  void ClipPolygon( G4ThreeVectorList& pPolygon,
295  const G4VoxelLimits& pVoxelLimit,
296  const EAxis pAxis ) const;
297  // Clip the specified convex polygon to the given limits, where
298  // the polygon is described by the vertices at (0),(1),...,(n),(0) in
299  // pPolygon.
300  // If the polygon is completely clipped away, the polygon is cleared.
301 
302  protected:
303 
304  G4double kCarTolerance; // Cached geometrical tolerance
305 
306  private:
307 
309  G4ThreeVectorList& outputPolygon,
310  const G4VoxelLimits& pVoxelLimit ) const;
311  // Clip the specified convex polygon to the given limits, storing the
312  // result in outputPolygon. The voxel limits must be limited in one
313  // *plane* only: This is achieved by having only x or y or z limits,
314  // and either the minimum or maximum limit set to -+kInfinity
315  // respectively.
316 
318 };
319 
320 #include "G4VSolid.icc"
321 
322 #endif
virtual G4Polyhedron * GetPolyhedron() const
Definition: G4VSolid.cc:644
G4String GetName() const
void ClipCrossSection(G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:347
virtual G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const =0
CLHEP::Hep3Vector G4ThreeVector
void ClipPolygonToSimpleLimits(G4ThreeVectorList &pPolygon, G4ThreeVectorList &outputPolygon, const G4VoxelLimits &pVoxelLimit) const
Definition: G4VSolid.cc:569
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4VSolid.cc:137
G4String name
Definition: TRTMaterials.hh:40
G4double EstimateCubicVolume(G4int nStat, G4double epsilon) const
Definition: G4VSolid.cc:203
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:188
G4double GetTolerance() const
virtual G4GeometryType GetEntityType() const =0
int G4int
Definition: G4Types.hh:78
void DumpInfo() const
static const double s
Definition: G4SIunits.hh:168
void SetName(const G4String &name)
virtual std::ostream & StreamInfo(std::ostream &os) const =0
void CalculateClippedPolygonExtent(G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:427
virtual EInside Inside(const G4ThreeVector &p) const =0
G4String G4GeometryType
Definition: G4VSolid.hh:85
bool G4bool
Definition: G4Types.hh:79
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4VSolid.hh:79
friend std::ostream & operator<<(std::ostream &os, const G4VSolid &e)
Definition: G4VSolid.cc:128
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
virtual G4Polyhedron * CreatePolyhedron() const
Definition: G4VSolid.cc:639
const G4int n
G4double EstimateSurfaceArea(G4int nStat, G4double ell) const
Definition: G4VSolid.cc:263
virtual ~G4VSolid()
Definition: G4VSolid.cc:101
virtual void DescribeYourselfTo(G4VGraphicsScene &scene) const =0
G4bool operator==(const G4VSolid &s) const
virtual G4VisExtent GetExtent() const
Definition: G4VSolid.cc:621
EInside
Definition: geomdefs.hh:58
EAxis
Definition: geomdefs.hh:54
virtual G4ThreeVector GetPointOnSurface() const
Definition: G4VSolid.cc:152
virtual const G4VSolid * GetConstituentSolid(G4int no) const
Definition: G4VSolid.cc:167
G4VSolid(const G4String &name)
Definition: G4VSolid.cc:60
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:110
G4double kCarTolerance
Definition: G4VSolid.hh:304
virtual G4VSolid * Clone() const
Definition: G4VSolid.cc:324
double G4double
Definition: G4Types.hh:76
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const =0
G4String fshapeName
Definition: G4VSolid.hh:317
virtual G4double GetSurfaceArea()
Definition: G4VSolid.cc:250
void ClipBetweenSections(G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:378
double epsilon(double density, double temperature)
virtual const G4DisplacedSolid * GetDisplacedSolidPtr() const
Definition: G4VSolid.cc:173
void ClipPolygon(G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis) const
Definition: G4VSolid.cc:494