Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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//
31// $Id$
32//
33//
34// --------------------------------------------------------------------
35// GEANT 4 class header file
36//
37// G4Tet
38//
39// Class description:
40//
41// A G4Tet is a tetrahedrasolid.
42//
43
44// History:
45// -------
46// 03.09.2004 - M.H.Mendenhall & R.A.Weller (Vanderbilt University, USA)
47// 10.02.2005 - D.Anninos (CERN) - Added GetPointOnSurface() method.
48// 12.11.2006 - M.H.Mendenhall - Added GetSurfaceArea() concrete implementation.
49// 20.09.2010 - G.Cosmo (CERN) - Added copy-ctor and operator=().
50// --------------------------------------------------------------------
51#ifndef G4TET_HH
52#define G4TET_HH
53
54#include "G4VSolid.hh"
55
56class G4Tet : public G4VSolid
57{
58
59 public: // with description
60
61 G4Tet(const G4String& pName,
62 G4ThreeVector anchor,
65 G4ThreeVector p4,
66 G4bool *degeneracyFlag=0);
67
68 virtual ~G4Tet();
69
71 const G4int n,
72 const G4VPhysicalVolume* pRep);
73
74 G4bool CalculateExtent(const EAxis pAxis,
75 const G4VoxelLimits& pVoxelLimit,
76 const G4AffineTransform& pTransform,
77 G4double& pmin, G4double& pmax) const;
78 // Methods for solid
79
80 EInside Inside(const G4ThreeVector& p) const;
81
83
84 G4double DistanceToIn(const G4ThreeVector& p, const G4ThreeVector& v) const;
85
86 G4double DistanceToIn(const G4ThreeVector& p) const;
87
89 const G4bool calcNorm=false,
90 G4bool *validNorm=0, G4ThreeVector *n=0) const;
91
92 G4double DistanceToOut(const G4ThreeVector& p) const;
93
96
98
99 G4VSolid* Clone() const;
100
101 std::ostream& StreamInfo(std::ostream& os) const;
102
104
105 // Functions for visualization
106
107 void DescribeYourselfTo (G4VGraphicsScene& scene) const;
108 G4VisExtent GetExtent () const;
110 G4NURBS* CreateNURBS () const;
111 G4Polyhedron* GetPolyhedron () const;
112
113 public: // without description
114
115 G4Tet(__void__&);
116 // Fake default constructor for usage restricted to direct object
117 // persistency for clients requiring preallocation of memory for
118 // persistifiable objects.
119
120 G4Tet(const G4Tet& rhs);
121 G4Tet& operator=(const G4Tet& rhs);
122 // Copy constructor and assignment operator.
123
124 const char* CVSHeaderVers()
125 { return "$Id: G4Tet.hh,v 1.11 2010-10-20 08:54:18 gcosmo Exp $"; }
126 const char* CVSFileVers()
127 { return CVSVers; }
129 { warningFlag=flag; }
131 G4ThreeVector p2,
132 G4ThreeVector p3,
133 G4ThreeVector p4);
134 std::vector<G4ThreeVector> GetVertices() const;
135 // Return the four vertices of the shape.
136
137 protected: // with description
138
140 CreateRotatedVertices(const G4AffineTransform& pTransform) const;
141 // Create the List of transformed vertices in the format required
142 // for G4VSolid:: ClipCrossSection and ClipBetweenSections.
143
144 private:
145
146 G4double fCubicVolume, fSurfaceArea;
147
148 mutable G4Polyhedron* fpPolyhedron;
149
150 G4ThreeVector GetPointOnFace(G4ThreeVector p1, G4ThreeVector p2,
151 G4ThreeVector p3, G4double& area) const;
152 static const char CVSVers[];
153
154 private:
155
156 G4ThreeVector fAnchor, fP2, fP3, fP4, fMiddle;
157 G4ThreeVector fNormal123, fNormal142, fNormal134, fNormal234;
158
159 G4bool warningFlag;
160
161 G4double fCdotN123, fCdotN142, fCdotN134, fCdotN234;
162 G4double fXMin, fXMax, fYMin, fYMax, fZMin, fZMax;
163 G4double fDx, fDy, fDz, fTol, fMaxSize;
164};
165
166#endif
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4VSolid.hh:85
Definition: G4Tet.hh:57
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Tet.cc:436
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Tet.cc:281
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4Tet.cc:291
void PrintWarnings(G4bool flag)
Definition: G4Tet.hh:128
G4Tet & operator=(const G4Tet &rhs)
Definition: G4Tet.cc:234
virtual ~G4Tet()
Definition: G4Tet.cc:204
G4Polyhedron * GetPolyhedron() const
Definition: G4Tet.cc:820
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Tet.cc:417
G4GeometryType GetEntityType() const
Definition: G4Tet.cc:650
G4ThreeVector GetPointOnSurface() const
Definition: G4Tet.cc:718
G4double GetSurfaceArea()
Definition: G4Tet.cc:763
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Tet.cc:668
static G4bool CheckDegeneracy(G4ThreeVector anchor, G4ThreeVector p2, G4ThreeVector p3, G4ThreeVector p4)
Definition: G4Tet.cc:265
const char * CVSFileVers()
Definition: G4Tet.hh:126
G4VSolid * Clone() const
Definition: G4Tet.cc:659
G4double GetCubicVolume()
Definition: G4Tet.cc:754
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
Definition: G4Tet.cc:529
G4Polyhedron * CreatePolyhedron() const
Definition: G4Tet.cc:792
std::vector< G4ThreeVector > GetVertices() const
Definition: G4Tet.cc:739
const char * CVSHeaderVers()
Definition: G4Tet.hh:124
G4ThreeVectorList * CreateRotatedVertices(const G4AffineTransform &pTransform) const
Definition: G4Tet.cc:619
G4VisExtent GetExtent() const
Definition: G4Tet.cc:783
EInside Inside(const G4ThreeVector &p) const
Definition: G4Tet.cc:386
G4NURBS * CreateNURBS() const
Definition: G4Tet.cc:811
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Tet.cc:774
EAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:58