Geant4-11
Public Member Functions | Protected Attributes | Static Protected Attributes | Private Member Functions | Private Attributes
G4QuadrangularFacet Class Reference

#include <G4QuadrangularFacet.hh>

Inheritance diagram for G4QuadrangularFacet:
G4VFacet

Public Member Functions

void ApplyTranslation (const G4ThreeVector v)
 
G4ThreeVector Distance (const G4ThreeVector &p)
 
G4double Distance (const G4ThreeVector &p, G4double minDist)
 
G4double Distance (const G4ThreeVector &p, G4double minDist, const G4bool outgoing)
 
G4double Extent (const G4ThreeVector axis)
 
 G4QuadrangularFacet (const G4QuadrangularFacet &right)
 
 G4QuadrangularFacet (const G4ThreeVector &Pt0, const G4ThreeVector &vt1, const G4ThreeVector &vt2, const G4ThreeVector &vt3, G4FacetVertexType)
 
G4double GetArea () const
 
G4ThreeVector GetCircumcentre () const
 
G4VFacetGetClone ()
 
G4GeometryType GetEntityType () const
 
G4int GetNumberOfVertices () const
 
G4ThreeVector GetPointOnFace () const
 
G4double GetRadius () const
 
G4ThreeVector GetSurfaceNormal () const
 
G4ThreeVector GetVertex (G4int i) const
 
G4bool Intersect (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool outgoing, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal)
 
G4bool IsDefined () const
 
G4bool IsInside (const G4ThreeVector &p) const
 
G4QuadrangularFacetoperator= (const G4QuadrangularFacet &right)
 
G4bool operator== (const G4VFacet &right) const
 
void SetVertex (G4int i, const G4ThreeVector &val)
 
void SetVertices (std::vector< G4ThreeVector > *v)
 
std::ostream & StreamInfo (std::ostream &os) const
 
 ~G4QuadrangularFacet ()
 

Protected Attributes

G4double kCarTolerance
 

Static Protected Attributes

static const G4double dirTolerance = 1.0E-14
 

Private Member Functions

G4int AllocatedMemory ()
 
G4int GetVertexIndex (G4int i) const
 
void SetVertexIndex (G4int i, G4int val)
 

Private Attributes

G4ThreeVector fCircumcentre
 
G4TriangularFacet fFacet1
 
G4TriangularFacet fFacet2
 
G4double fRadius = 0.0
 

Detailed Description

Definition at line 61 of file G4QuadrangularFacet.hh.

Constructor & Destructor Documentation

◆ G4QuadrangularFacet() [1/2]

G4QuadrangularFacet::G4QuadrangularFacet ( const G4ThreeVector Pt0,
const G4ThreeVector vt1,
const G4ThreeVector vt2,
const G4ThreeVector vt3,
G4FacetVertexType  vertexType 
)

Definition at line 51 of file G4QuadrangularFacet.cc.

56 : G4VFacet()
57{
58 G4double delta = 1.0 * kCarTolerance; // dimension tolerance
59 G4double epsilon = 0.01 * kCarTolerance; // planarity tolerance
60
62 SetVertex(0, vt0);
63 if (vertexType == ABSOLUTE)
64 {
65 SetVertex(1, vt1);
66 SetVertex(2, vt2);
67 SetVertex(3, vt3);
68
69 e1 = vt1 - vt0;
70 e2 = vt2 - vt0;
71 e3 = vt3 - vt0;
72 }
73 else
74 {
75 SetVertex(1, vt0 + vt1);
76 SetVertex(2, vt0 + vt2);
77 SetVertex(3, vt0 + vt3);
78
79 e1 = vt1;
80 e2 = vt2;
81 e3 = vt3;
82 }
83
84 // Check length of sides and diagonals
85 //
86 G4double leng1 = e1.mag();
87 G4double leng2 = (e2-e1).mag();
88 G4double leng3 = (e3-e2).mag();
89 G4double leng4 = e3.mag();
90
91 G4double diag1 = e2.mag();
92 G4double diag2 = (e3-e1).mag();
93
94 if (leng1 <= delta || leng2 <= delta || leng3 <= delta || leng4 <= delta ||
95 diag1 <= delta || diag2 <= delta)
96 {
97 ostringstream message;
98 message << "Sides/diagonals of facet are too small." << G4endl
99 << "P0 = " << GetVertex(0) << G4endl
100 << "P1 = " << GetVertex(1) << G4endl
101 << "P2 = " << GetVertex(2) << G4endl
102 << "P3 = " << GetVertex(3) << G4endl
103 << "Side1 length (P0->P1) = " << leng1 << G4endl
104 << "Side2 length (P1->P2) = " << leng2 << G4endl
105 << "Side3 length (P2->P3) = " << leng3 << G4endl
106 << "Side4 length (P3->P0) = " << leng4 << G4endl
107 << "Diagonal1 length (P0->P2) = " << diag1 << G4endl
108 << "Diagonal2 length (P1->P3) = " << diag2;
109 G4Exception("G4QuadrangularFacet::G4QuadrangularFacet()",
110 "GeomSolids1001", JustWarning, message);
111 return;
112 }
113
114 // Check that vertices are not collinear
115 //
116 G4double s1 = (e1.cross(e2)).mag()*0.5;
117 G4double s2 = ((e2-e1).cross(e3-e2)).mag()*0.5;
118 G4double s3 = (e2.cross(e3)).mag()*0.5;
119 G4double s4 = (e1.cross(e3)).mag()*0.5;
120
121 G4double h1 = 2.*s1 / std::max(std::max(leng1,leng2),diag1);
122 G4double h2 = 2.*s2 / std::max(std::max(leng2,leng3),diag2);
123 G4double h3 = 2.*s3 / std::max(std::max(leng3,leng4),diag1);
124 G4double h4 = 2.*s4 / std::max(std::max(leng4,leng1),diag2);
125
126 if (h1 <= delta || h2 <= delta || h3 <= delta || h4 <= delta )
127 {
128 ostringstream message;
129 message << "Facet has three or more collinear vertices." << G4endl
130 << "P0 = " << GetVertex(0) << G4endl
131 << "P1 = " << GetVertex(1) << G4endl
132 << "P2 = " << GetVertex(2) << G4endl
133 << "P3 = " << GetVertex(3) << G4endl
134 << "Height in P0-P1-P2 = " << h1 << G4endl
135 << "Height in P1-P2-P3 = " << h2 << G4endl
136 << "Height in P2-P3-P4 = " << h3 << G4endl
137 << "Height in P4-P0-P1 = " << h4;
138 G4Exception("G4QuadrangularFacet::G4QuadrangularFacet()",
139 "GeomSolids1001", JustWarning, message);
140 return;
141 }
142
143 // Check that vertices are coplanar by computing minimal
144 // height of tetrahedron comprising of vertices
145 //
146 G4double smax = std::max( std::max(s1,s2), std::max(s3,s4) );
147 G4double hmin = 0.5 * std::fabs( e1.dot(e2.cross(e3)) ) / smax;
148 if (hmin >= epsilon)
149 {
150 ostringstream message;
151 message << "Facet is not planar." << G4endl
152 << "Disrepancy = " << hmin << G4endl
153 << "P0 = " << GetVertex(0) << G4endl
154 << "P1 = " << GetVertex(1) << G4endl
155 << "P2 = " << GetVertex(2) << G4endl
156 << "P3 = " << GetVertex(3);
157 G4Exception("G4QuadrangularFacet::G4QuadrangularFacet()",
158 "GeomSolids1001", JustWarning, message);
159 return;
160 }
161
162 // Check that facet is convex by computing crosspoint
163 // of diagonals
164 //
165 G4ThreeVector normal = e2.cross(e3-e1);
166 G4double s = kInfinity, t = kInfinity, magnitude2 = normal.mag2();
167 if (magnitude2 > delta*delta) // check: magnitude2 != 0.
168 {
169 s = normal.dot(e1.cross(e3-e1)) / magnitude2;
170 t = normal.dot(e1.cross(e2)) / magnitude2;
171 }
172 if (s <= 0. || s >= 1. || t <= 0. || t >= 1.)
173 {
174 ostringstream message;
175 message << "Facet is not convex." << G4endl
176 << "Parameters of crosspoint of diagonals: "
177 << s << " and " << t << G4endl
178 << "should both be within (0,1) range" << G4endl
179 << "P0 = " << GetVertex(0) << G4endl
180 << "P1 = " << GetVertex(1) << G4endl
181 << "P2 = " << GetVertex(2) << G4endl
182 << "P3 = " << GetVertex(3);
183 G4Exception("G4QuadrangularFacet::G4QuadrangularFacet()",
184 "GeomSolids1001", JustWarning, message);
185 return;
186 }
187
188 // Define facet
189 //
192
193 normal = normal.unit();
196
197 G4ThreeVector vtmp = 0.5 * (e1 + e2);
198 fCircumcentre = GetVertex(0) + vtmp;
199 G4double radiusSqr = vtmp.mag2();
200 fRadius = std::sqrt(radiusSqr);
201 // 29.02.2016 Remark by E.Tcherniaev: computation
202 // of fCircumcenter and fRadius is wrong, however
203 // it did not create any problem till now.
204 // Bizarre! Need to investigate!
205}
static const G4double e1[44]
static const G4double e2[44]
static const G4double e3[45]
G4double epsilon(G4double density, G4double temperature)
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
static constexpr double s
Definition: G4SIunits.hh:154
double G4double
Definition: G4Types.hh:83
@ ABSOLUTE
Definition: G4VFacet.hh:48
#define G4endl
Definition: G4ios.hh:57
double mag2() const
G4TriangularFacet fFacet1
void SetVertex(G4int i, const G4ThreeVector &val)
G4ThreeVector GetVertex(G4int i) const
G4TriangularFacet fFacet2
void SetSurfaceNormal(G4ThreeVector normal)
G4VFacet()
Definition: G4VFacet.cc:44
G4double kCarTolerance
Definition: G4VFacet.hh:93
static const G4double kInfinity
Definition: geomdefs.hh:41
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:79
T max(const T t1, const T t2)
brief Return the largest of the two arguments

References ABSOLUTE, e1, e2, e3, epsilon(), fCircumcentre, fFacet1, fFacet2, fRadius, G4endl, G4Exception(), GetVertex(), JustWarning, G4VFacet::kCarTolerance, kInfinity, CLHEP::Hep3Vector::mag2(), G4INCL::Math::max(), CLHEP::normal(), s, G4TriangularFacet::SetSurfaceNormal(), and SetVertex().

Referenced by GetClone().

◆ G4QuadrangularFacet() [2/2]

G4QuadrangularFacet::G4QuadrangularFacet ( const G4QuadrangularFacet right)

Definition at line 215 of file G4QuadrangularFacet.cc.

216 : G4VFacet(rhs)
217{
218 fFacet1 = rhs.fFacet1;
219 fFacet2 = rhs.fFacet2;
220 fRadius = 0.0;
221}

References fFacet1, fFacet2, and fRadius.

◆ ~G4QuadrangularFacet()

G4QuadrangularFacet::~G4QuadrangularFacet ( )

Definition at line 209 of file G4QuadrangularFacet.cc.

210{
211}

Member Function Documentation

◆ AllocatedMemory()

G4int G4QuadrangularFacet::AllocatedMemory ( )
inlineprivatevirtual

Implements G4VFacet.

Definition at line 199 of file G4QuadrangularFacet.hh.

200{
201 return sizeof(*this) + fFacet1.AllocatedMemory() + fFacet2.AllocatedMemory();
202}

References G4TriangularFacet::AllocatedMemory(), fFacet1, and fFacet2.

◆ ApplyTranslation()

void G4VFacet::ApplyTranslation ( const G4ThreeVector  v)
inherited

Definition at line 85 of file G4VFacet.cc.

86{
88 for (G4int i = 0; i < n; ++i)
89 {
90 SetVertex(i, GetVertex(i) + v);
91 }
92}
int G4int
Definition: G4Types.hh:85
virtual G4ThreeVector GetVertex(G4int i) const =0
virtual G4int GetNumberOfVertices() const =0
virtual void SetVertex(G4int i, const G4ThreeVector &val)=0

References G4VFacet::GetNumberOfVertices(), G4VFacet::GetVertex(), CLHEP::detail::n, and G4VFacet::SetVertex().

◆ Distance() [1/3]

G4ThreeVector G4QuadrangularFacet::Distance ( const G4ThreeVector p)

Definition at line 249 of file G4QuadrangularFacet.cc.

250{
253
254 if (v1.mag2() < v2.mag2()) return v1;
255 else return v2;
256}
G4ThreeVector Distance(const G4ThreeVector &p)

References G4TriangularFacet::Distance(), fFacet1, fFacet2, and CLHEP::Hep3Vector::mag2().

Referenced by Distance().

◆ Distance() [2/3]

G4double G4QuadrangularFacet::Distance ( const G4ThreeVector p,
G4double  minDist 
)
virtual

Implements G4VFacet.

Definition at line 260 of file G4QuadrangularFacet.cc.

262{
263 G4double dist = Distance(p).mag();
264 return dist;
265}
double mag() const
G4ThreeVector Distance(const G4ThreeVector &p)

References Distance(), and CLHEP::Hep3Vector::mag().

◆ Distance() [3/3]

G4double G4QuadrangularFacet::Distance ( const G4ThreeVector p,
G4double  minDist,
const G4bool  outgoing 
)
virtual

Implements G4VFacet.

Definition at line 269 of file G4QuadrangularFacet.cc.

271{
272 G4double dist;
273
274 G4ThreeVector v = Distance(p);
275 G4double dir = v.dot(GetSurfaceNormal());
276 if ( ((dir > dirTolerance) && (!outgoing))
277 || ((dir < -dirTolerance) && outgoing))
278 dist = kInfinity;
279 else
280 dist = v.mag();
281 return dist;
282}
double dot(const Hep3Vector &) const
G4ThreeVector GetSurfaceNormal() const
static const G4double dirTolerance
Definition: G4VFacet.hh:92

References G4VFacet::dirTolerance, Distance(), CLHEP::Hep3Vector::dot(), GetSurfaceNormal(), kInfinity, and CLHEP::Hep3Vector::mag().

◆ Extent()

G4double G4QuadrangularFacet::Extent ( const G4ThreeVector  axis)
virtual

Implements G4VFacet.

Definition at line 286 of file G4QuadrangularFacet.cc.

287{
288 G4double ss = 0;
289
290 for (G4int i = 0; i <= 3; ++i)
291 {
292 G4double sp = GetVertex(i).dot(axis);
293 if (sp > ss) ss = sp;
294 }
295 return ss;
296}

References CLHEP::Hep3Vector::dot(), GetVertex(), and G4InuclParticleNames::sp.

◆ GetArea()

G4double G4QuadrangularFacet::GetArea ( ) const
virtual

Implements G4VFacet.

Definition at line 335 of file G4QuadrangularFacet.cc.

336{
338 return area;
339}
G4double GetArea() const

References fFacet1, fFacet2, and G4TriangularFacet::GetArea().

◆ GetCircumcentre()

G4ThreeVector G4QuadrangularFacet::GetCircumcentre ( ) const
inlinevirtual

Implements G4VFacet.

Definition at line 135 of file G4QuadrangularFacet.hh.

136{
137 return fCircumcentre;
138}

References fCircumcentre.

◆ GetClone()

G4VFacet * G4QuadrangularFacet::GetClone ( )
virtual

Implements G4VFacet.

Definition at line 239 of file G4QuadrangularFacet.cc.

240{
242 GetVertex(2), GetVertex(3),
243 ABSOLUTE);
244 return c;
245}
G4QuadrangularFacet(const G4ThreeVector &Pt0, const G4ThreeVector &vt1, const G4ThreeVector &vt2, const G4ThreeVector &vt3, G4FacetVertexType)

References ABSOLUTE, G4QuadrangularFacet(), and GetVertex().

◆ GetEntityType()

G4String G4QuadrangularFacet::GetEntityType ( ) const
virtual

Implements G4VFacet.

Definition at line 343 of file G4QuadrangularFacet.cc.

344{
345 return "G4QuadrangularFacet";
346}

◆ GetNumberOfVertices()

G4int G4QuadrangularFacet::GetNumberOfVertices ( ) const
inlinevirtual

Implements G4VFacet.

Definition at line 119 of file G4QuadrangularFacet.hh.

120{
121 return 4;
122}

◆ GetPointOnFace()

G4ThreeVector G4QuadrangularFacet::GetPointOnFace ( ) const
virtual

Implements G4VFacet.

Definition at line 323 of file G4QuadrangularFacet.cc.

324{
325 G4double s1 = fFacet1.GetArea();
326 G4double s2 = fFacet2.GetArea();
327 return ((s1+s2)*G4UniformRand() < s1) ?
329}
#define G4UniformRand()
Definition: Randomize.hh:52
G4ThreeVector GetPointOnFace() const

References fFacet1, fFacet2, G4UniformRand, G4TriangularFacet::GetArea(), and G4TriangularFacet::GetPointOnFace().

◆ GetRadius()

G4double G4QuadrangularFacet::GetRadius ( ) const
inlinevirtual

Implements G4VFacet.

Definition at line 130 of file G4QuadrangularFacet.hh.

131{
132 return fRadius;
133}

References fRadius.

◆ GetSurfaceNormal()

G4ThreeVector G4QuadrangularFacet::GetSurfaceNormal ( ) const
virtual

Implements G4VFacet.

Definition at line 350 of file G4QuadrangularFacet.cc.

351{
352 return fFacet1.GetSurfaceNormal();
353}
G4ThreeVector GetSurfaceNormal() const

References fFacet1, and G4TriangularFacet::GetSurfaceNormal().

Referenced by Distance().

◆ GetVertex()

G4ThreeVector G4QuadrangularFacet::GetVertex ( G4int  i) const
inlinevirtual

Implements G4VFacet.

Definition at line 124 of file G4QuadrangularFacet.hh.

125{
126 return i == 3 ? fFacet2.GetVertex(2) : fFacet1.GetVertex(i);
127}
G4ThreeVector GetVertex(G4int i) const

References fFacet1, fFacet2, and G4TriangularFacet::GetVertex().

Referenced by Extent(), G4QuadrangularFacet(), and GetClone().

◆ GetVertexIndex()

G4int G4QuadrangularFacet::GetVertexIndex ( G4int  i) const
inlineprivatevirtual

Implements G4VFacet.

Definition at line 172 of file G4QuadrangularFacet.hh.

173{
174 return i == 3 ? fFacet2.GetVertexIndex(2) : fFacet1.GetVertexIndex(i);
175}
G4int GetVertexIndex(G4int i) const

References fFacet1, fFacet2, and G4TriangularFacet::GetVertexIndex().

◆ Intersect()

G4bool G4QuadrangularFacet::Intersect ( const G4ThreeVector p,
const G4ThreeVector v,
const G4bool  outgoing,
G4double distance,
G4double distFromSurface,
G4ThreeVector normal 
)
virtual

Implements G4VFacet.

Definition at line 300 of file G4QuadrangularFacet.cc.

306{
307 G4bool intersect =
308 fFacet1.Intersect(p,v,outgoing,distance,distFromSurface,normal);
309 if (!intersect) intersect =
310 fFacet2.Intersect(p,v,outgoing,distance,distFromSurface,normal);
311 if (!intersect)
312 {
313 distance = distFromSurface = kInfinity;
314 normal.set(0,0,0);
315 }
316 return intersect;
317}
bool G4bool
Definition: G4Types.hh:86
G4bool Intersect(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool outgoing, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal)

References fFacet1, fFacet2, G4TriangularFacet::Intersect(), kInfinity, and CLHEP::normal().

◆ IsDefined()

G4bool G4QuadrangularFacet::IsDefined ( ) const
inlinevirtual

Implements G4VFacet.

Definition at line 167 of file G4QuadrangularFacet.hh.

168{
169 return fFacet1.IsDefined();
170}
G4bool IsDefined() const

References fFacet1, and G4TriangularFacet::IsDefined().

◆ IsInside()

G4bool G4VFacet::IsInside ( const G4ThreeVector p) const
inherited

Definition at line 112 of file G4VFacet.cc.

113{
114 G4ThreeVector d = p-GetVertex(0);
115 G4double displacement = d.dot(GetSurfaceNormal());
116 return displacement <= 0.0;
117}
virtual G4ThreeVector GetSurfaceNormal() const =0

References CLHEP::Hep3Vector::dot(), G4VFacet::GetSurfaceNormal(), and G4VFacet::GetVertex().

Referenced by G4TessellatedSolid::SetExtremeFacets().

◆ operator=()

G4QuadrangularFacet & G4QuadrangularFacet::operator= ( const G4QuadrangularFacet right)

Definition at line 226 of file G4QuadrangularFacet.cc.

227{
228 if (this == &rhs) return *this;
229
230 fFacet1 = rhs.fFacet1;
231 fFacet2 = rhs.fFacet2;
232 fRadius = 0.0;
233
234 return *this;
235}

References fFacet1, fFacet2, and fRadius.

◆ operator==()

G4bool G4VFacet::operator== ( const G4VFacet right) const
inherited

Definition at line 57 of file G4VFacet.cc.

58{
59 G4double tolerance = kCarTolerance*kCarTolerance/4.0;
60
62 return false;
63 else if ((GetCircumcentre()-right.GetCircumcentre()).mag2() > tolerance)
64 return false;
65 else if (std::fabs((right.GetSurfaceNormal()).dot(GetSurfaceNormal())) < 0.9999999999)
66 return false;
67
68 G4bool coincident = true;
69 G4int i = 0;
70 do // Loop checking, 13.08.2015, G.Cosmo
71 {
72 coincident = false;
73 G4int j = 0;
74 do // Loop checking, 13.08.2015, G.Cosmo
75 {
76 coincident = (GetVertex(i)-right.GetVertex(j)).mag2() < tolerance;
77 } while (!coincident && ++j < GetNumberOfVertices());
78 } while (coincident && ++i < GetNumberOfVertices());
79
80 return coincident;
81}
virtual G4ThreeVector GetCircumcentre() const =0

References G4VFacet::GetCircumcentre(), G4VFacet::GetNumberOfVertices(), G4VFacet::GetSurfaceNormal(), G4VFacet::GetVertex(), and G4VFacet::kCarTolerance.

◆ SetVertex()

void G4QuadrangularFacet::SetVertex ( G4int  i,
const G4ThreeVector val 
)
inlinevirtual

Implements G4VFacet.

Definition at line 140 of file G4QuadrangularFacet.hh.

141{
142 switch (i)
143 {
144 case 0:
145 fFacet1.SetVertex(0, val);
146 fFacet2.SetVertex(0, val);
147 break;
148 case 1:
149 fFacet1.SetVertex(1, val);
150 break;
151 case 2:
152 fFacet1.SetVertex(2, val);
153 fFacet2.SetVertex(1, val);
154 break;
155 case 3:
156 fFacet2.SetVertex(2, val);
157 break;
158 }
159}
void SetVertex(G4int i, const G4ThreeVector &val)

References fFacet1, fFacet2, and G4TriangularFacet::SetVertex().

Referenced by G4QuadrangularFacet().

◆ SetVertexIndex()

void G4QuadrangularFacet::SetVertexIndex ( G4int  i,
G4int  val 
)
inlineprivatevirtual

Implements G4VFacet.

Definition at line 178 of file G4QuadrangularFacet.hh.

179{
180 switch (i)
181 {
182 case 0:
183 fFacet1.SetVertexIndex(0, val);
184 fFacet2.SetVertexIndex(0, val);
185 break;
186 case 1:
187 fFacet1.SetVertexIndex(1, val);
188 break;
189 case 2:
190 fFacet1.SetVertexIndex(2, val);
191 fFacet2.SetVertexIndex(1, val);
192 break;
193 case 3:
194 fFacet2.SetVertexIndex(2, val);
195 break;
196 }
197}
void SetVertexIndex(G4int i, G4int j)

References fFacet1, fFacet2, and G4TriangularFacet::SetVertexIndex().

◆ SetVertices()

void G4QuadrangularFacet::SetVertices ( std::vector< G4ThreeVector > *  v)
inlinevirtual

Implements G4VFacet.

Definition at line 161 of file G4QuadrangularFacet.hh.

162{
165}
void SetVertices(std::vector< G4ThreeVector > *v)

References fFacet1, fFacet2, and G4TriangularFacet::SetVertices().

◆ StreamInfo()

std::ostream & G4VFacet::StreamInfo ( std::ostream &  os) const
inherited

Definition at line 96 of file G4VFacet.cc.

97{
98 os << G4endl;
99 os << "*********************************************************************"
100 << G4endl;
101 os << "FACET TYPE = " << GetEntityType() << G4endl;
102 os << "ABSOLUTE VECTORS = " << G4endl;
104 for (G4int i = 0; i < n; ++i)
105 os << "P[" << i << "] = " << GetVertex(i) << G4endl;
106 os << "*********************************************************************"
107 << G4endl;
108
109 return os;
110}
virtual G4GeometryType GetEntityType() const =0

References G4endl, G4VFacet::GetEntityType(), G4VFacet::GetNumberOfVertices(), G4VFacet::GetVertex(), and CLHEP::detail::n.

Referenced by G4TessellatedSolid::AddFacet(), and G4TessellatedSolid::StreamInfo().

Field Documentation

◆ dirTolerance

const G4double G4VFacet::dirTolerance = 1.0E-14
staticprotectedinherited

Definition at line 92 of file G4VFacet.hh.

Referenced by Distance(), and G4TriangularFacet::Intersect().

◆ fCircumcentre

G4ThreeVector G4QuadrangularFacet::fCircumcentre
private

Definition at line 110 of file G4QuadrangularFacet.hh.

Referenced by G4QuadrangularFacet(), and GetCircumcentre().

◆ fFacet1

G4TriangularFacet G4QuadrangularFacet::fFacet1
private

◆ fFacet2

G4TriangularFacet G4QuadrangularFacet::fFacet2
private

◆ fRadius

G4double G4QuadrangularFacet::fRadius = 0.0
private

Definition at line 108 of file G4QuadrangularFacet.hh.

Referenced by G4QuadrangularFacet(), GetRadius(), and operator=().

◆ kCarTolerance

G4double G4VFacet::kCarTolerance
protectedinherited

The documentation for this class was generated from the following files: