Geant4  10.03
G4VCSGfaceted.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: G4VCSGfaceted.hh 83572 2014-09-01 15:23:27Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------------
31 // GEANT 4 class header file
32 //
33 //
34 // G4VCSGfaceted
35 //
36 // Class description:
37 //
38 // Virtual class defining CSG-like type shape that is built entire
39 // of G4CSGface faces.
40 
41 // Author:
42 // David C. Williams (davidw@scipp.ucsc.edu)
43 // --------------------------------------------------------------------
44 
45 #ifndef G4VCSGfaceted_hh
46 #define G4VCSGfaceted_hh
47 
48 #include "G4VSolid.hh"
49 
50 class G4VCSGface;
51 class G4VisExtent;
52 
53 class G4VCSGfaceted : public G4VSolid
54 {
55  public: // with description
56 
57  G4VCSGfaceted( const G4String& name );
58  virtual ~G4VCSGfaceted();
59 
60  G4VCSGfaceted( const G4VCSGfaceted &source );
61  G4VCSGfaceted &operator=( const G4VCSGfaceted &source );
62 
63  virtual G4bool CalculateExtent( const EAxis pAxis,
64  const G4VoxelLimits& pVoxelLimit,
65  const G4AffineTransform& pTransform,
66  G4double& pmin,G4double& pmax ) const;
67 
68  virtual EInside Inside( const G4ThreeVector& p ) const;
69 
70  virtual G4ThreeVector SurfaceNormal( const G4ThreeVector& p ) const;
71 
72  virtual G4double DistanceToIn( const G4ThreeVector& p,
73  const G4ThreeVector& v ) const;
74  virtual G4double DistanceToIn( const G4ThreeVector& p ) const;
75  virtual G4double DistanceToOut( const G4ThreeVector& p,
76  const G4ThreeVector& v,
77  const G4bool calcNorm=false,
78  G4bool *validNorm=0,
79  G4ThreeVector *n=0 ) const;
80  virtual G4double DistanceToOut( const G4ThreeVector& p ) const;
81 
82  virtual G4GeometryType GetEntityType() const;
83 
84  virtual std::ostream& StreamInfo(std::ostream& os) const;
85 
86  virtual G4Polyhedron* CreatePolyhedron() const = 0;
87 
88  virtual void DescribeYourselfTo( G4VGraphicsScene& scene ) const;
89 
90  virtual G4VisExtent GetExtent() const;
91 
92  virtual G4Polyhedron* GetPolyhedron () const;
93 
94  G4int GetCubVolStatistics() const;
95  G4double GetCubVolEpsilon() const;
96  void SetCubVolStatistics(G4int st);
97  void SetCubVolEpsilon(G4double ep);
98  G4int GetAreaStatistics() const;
99  G4double GetAreaAccuracy() const;
100  void SetAreaStatistics(G4int st);
101  void SetAreaAccuracy(G4double ep);
102 
103  virtual G4double GetCubicVolume();
104  // Returns an estimation of the geometrical cubic volume of the
105  // solid. Caches the computed value once computed the first time.
106  virtual G4double GetSurfaceArea();
107  // Returns an estimation of the geometrical surface area of the
108  // solid. Caches the computed value once computed the first time.
109 
110  public: // without description
111 
112  G4VCSGfaceted(__void__&);
113  // Fake default constructor for usage restricted to direct object
114  // persistency for clients requiring preallocation of memory for
115  // persistifiable objects.
116 
117  protected: // without description
118 
125 
126  virtual G4double DistanceTo( const G4ThreeVector &p,
127  const G4bool outgoing ) const;
128 
130  // Returns a random point located on the surface of the solid
131  // in case of generic Polycone or generic Polyhedra.
132 
133  void CopyStuff( const G4VCSGfaceted &source );
134  void DeleteStuff();
135 
136  private:
137 
141  // Statistics, error accuracy for volume estimation.
142 
143 };
144 
145 #endif
virtual void DescribeYourselfTo(G4VGraphicsScene &scene) const
G4double fSurfaceArea
void SetAreaStatistics(G4int st)
G4VCSGfaceted & operator=(const G4VCSGfaceted &source)
G4bool fRebuildPolyhedron
G4Polyhedron * fpPolyhedron
CLHEP::Hep3Vector G4ThreeVector
G4int GetCubVolStatistics() const
void SetCubVolStatistics(G4int st)
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
virtual ~G4VCSGfaceted()
G4VCSGfaceted(const G4String &name)
void CopyStuff(const G4VCSGfaceted &source)
const char * name(G4int ptype)
int G4int
Definition: G4Types.hh:78
G4double fCubVolEpsilon
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
G4ThreeVector GetPointOnSurfaceGeneric() const
virtual G4double GetSurfaceArea()
G4double fCubicVolume
virtual G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
virtual G4Polyhedron * CreatePolyhedron() const =0
G4double GetAreaAccuracy() const
bool G4bool
Definition: G4Types.hh:79
virtual G4VisExtent GetExtent() const
const G4int n
virtual G4GeometryType GetEntityType() const
G4VCSGface ** faces
virtual std::ostream & StreamInfo(std::ostream &os) const
EInside
Definition: geomdefs.hh:58
EAxis
Definition: geomdefs.hh:54
G4double fAreaAccuracy
virtual G4double DistanceTo(const G4ThreeVector &p, const G4bool outgoing) const
void SetAreaAccuracy(G4double ep)
G4int GetAreaStatistics() const
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
virtual EInside Inside(const G4ThreeVector &p) const
virtual G4double GetCubicVolume()
G4double GetCubVolEpsilon() const
void SetCubVolEpsilon(G4double ep)
virtual G4Polyhedron * GetPolyhedron() const