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