Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4Tet.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 intellectual property of the *
19// * Vanderbilt University Free Electron Laser Center *
20// * Vanderbilt University, Nashville, TN, USA *
21// * Development supported by: *
22// * United States MFEL program under grant FA9550-04-1-0045 *
23// * and NASA under contract number NNG04CT05P. *
24// * Written by Marcus H. Mendenhall and Robert A. Weller. *
25// * *
26// * Contributed to the Geant4 Core, January, 2005. *
27// * *
28// ********************************************************************
29//
30// G4Tet
31//
32// Class description:
33//
34// A G4Tet is a tetrahedra solid.
35
36// 03.09.2004 - M.H.Mendenhall & R.A.Weller (Vanderbilt University, USA)
37// 08.01.2020 - E.Tcherniaev, complete revision, speed up
38// --------------------------------------------------------------------
39#ifndef G4TET_HH
40#define G4TET_HH
41
42#include "G4GeomTypes.hh"
43
44#if defined(G4GEOM_USE_USOLIDS)
45#define G4GEOM_USE_UTET 1
46#endif
47
48#if defined(G4GEOM_USE_UTET)
49 #define G4UTet G4Tet
50 #include "G4UTet.hh"
51#else
52
53#include "G4VSolid.hh"
54
55class G4Tet : public G4VSolid
56{
57
58 public: // with description
59
60 // Constructor
61 G4Tet(const G4String& pName,
62 const G4ThreeVector& anchor,
63 const G4ThreeVector& p1,
64 const G4ThreeVector& p2,
65 const G4ThreeVector& p3,
66 G4bool* degeneracyFlag = nullptr);
67
68 // Destructor
69 virtual ~G4Tet();
70
71 // Modifier
72 void SetVertices(const G4ThreeVector& anchor,
73 const G4ThreeVector& p1,
74 const G4ThreeVector& p2,
75 const G4ThreeVector& p3);
76
77 // Accessors, return the four vertices of the shape
78 void GetVertices(G4ThreeVector& anchor,
79 G4ThreeVector& p1,
80 G4ThreeVector& p2,
81 G4ThreeVector& p3) const;
82 std::vector<G4ThreeVector> GetVertices() const;
83
84 // Set warning flag - depricated (dummy)
86
87 // Return true if the tetrahedron is degenerate
89 const G4ThreeVector& p1,
90 const G4ThreeVector& p2,
91 const G4ThreeVector& p3) const;
92
93 // Standard methods
95 const G4int n,
96 const G4VPhysicalVolume* pRep);
97
98 void BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const;
99 G4bool CalculateExtent(const EAxis pAxis,
100 const G4VoxelLimits& pVoxelLimit,
101 const G4AffineTransform& pTransform,
102 G4double& pmin, G4double& pmax) const;
103
104 EInside Inside(const G4ThreeVector& p) const;
107 const G4ThreeVector& v) const;
108 G4double DistanceToIn(const G4ThreeVector& p) const;
110 const G4ThreeVector& v,
111 const G4bool calcNorm = false,
112 G4bool* validNorm = nullptr,
113 G4ThreeVector* n = nullptr) const;
114 G4double DistanceToOut(const G4ThreeVector& p) const;
115
117
118 G4VSolid* Clone() const;
119
120 std::ostream& StreamInfo(std::ostream& os) const;
121
124
126
127 // Methods for visualization
128 void DescribeYourselfTo (G4VGraphicsScene& scene) const;
129 G4VisExtent GetExtent () const;
131 G4Polyhedron* GetPolyhedron () const;
132
133 public: // without description
134
135 // Fake default constructor for usage restricted to direct object
136 // persistency for clients requiring preallocation of memory for
137 // persistifiable objects
138 G4Tet(__void__&);
139
140 // Copy constructor
141 G4Tet(const G4Tet& rhs);
142
143 // Assignment operator
144 G4Tet& operator=(const G4Tet& rhs);
145
146 private:
147
148 // Set data members
149 void Initialize(const G4ThreeVector& p0,
150 const G4ThreeVector& p1,
151 const G4ThreeVector& p2,
152 const G4ThreeVector& p3);
153
154 // Return normal to surface closest to p
155 G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector& p) const;
156
157 private:
158
159 G4double halfTolerance = 0;
160 G4double fCubicVolume = 0; // Volume
161 G4double fSurfaceArea = 0; // Surface area
162 mutable G4bool fRebuildPolyhedron = false;
163 mutable G4Polyhedron* fpPolyhedron = nullptr;
164
165 G4ThreeVector fVertex[4]; // thetrahedron vertices
166 G4ThreeVector fNormal[4]; // normals to faces
167 G4double fDist[4] = {0}; // distances from origin to faces
168 G4double fArea[4] = {0}; // face areas
169 G4ThreeVector fBmin, fBmax; // bounding box
170};
171
172#endif
173
174#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
Definition: G4Tet.hh:56
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Tet.cc:447
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Tet.cc:309
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4Tet.cc:329
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const
Definition: G4Tet.cc:487
G4Tet & operator=(const G4Tet &rhs)
Definition: G4Tet.cc:140
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4Tet.cc:319
virtual ~G4Tet()
Definition: G4Tet.cc:113
G4Polyhedron * GetPolyhedron() const
Definition: G4Tet.cc:678
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Tet.cc:393
G4GeometryType GetEntityType() const
Definition: G4Tet.cc:544
void SetVertices(const G4ThreeVector &anchor, const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3)
Definition: G4Tet.cc:250
G4ThreeVector GetPointOnSurface() const
Definition: G4Tet.cc:583
void PrintWarnings(G4bool)
Definition: G4Tet.hh:85
G4double GetSurfaceArea()
Definition: G4Tet.cc:617
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Tet.cc:562
G4VSolid * Clone() const
Definition: G4Tet.cc:553
G4double GetCubicVolume()
Definition: G4Tet.cc:608
G4bool CheckDegeneracy(const G4ThreeVector &p0, const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3) const
Definition: G4Tet.cc:173
G4Polyhedron * CreatePolyhedron() const
Definition: G4Tet.cc:646
std::vector< G4ThreeVector > GetVertices() const
Definition: G4Tet.cc:297
G4VisExtent GetExtent() const
Definition: G4Tet.cc:635
EInside Inside(const G4ThreeVector &p) const
Definition: G4Tet.cc:379
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Tet.cc:626
EAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:67