Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4QuadrangularFacet Class Reference

G4QuadrangularFacet defines a facet with 4 vertices, used for the contruction of G4TessellatedSolid. Vertices shall be in the same plane and be supplied in anti-clockwise order looking from the outsider of the solid where it belongs. More...

#include <G4QuadrangularFacet.hh>

Inheritance diagram for G4QuadrangularFacet:

Public Member Functions

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

Additional Inherited Members

Protected Attributes inherited from G4VFacet
G4double kCarTolerance
Static Protected Attributes inherited from G4VFacet
static const G4double dirTolerance = 1.0E-14

Detailed Description

G4QuadrangularFacet defines a facet with 4 vertices, used for the contruction of G4TessellatedSolid. Vertices shall be in the same plane and be supplied in anti-clockwise order looking from the outsider of the solid where it belongs.

Definition at line 68 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 vType )

Constructs a facet with 4 vertices, given its parameters.

Parameters
[in]Pt0The anchor point, first vertex.
[in]vt1Second vertex.
[in]vt2Third vertex.
[in]vt3Fourth vertex.
[in]vTypeThe positioning type for the vertices, either: "ABSOLUTE" - vertices set in anti-clockwise order when looking from the outsider. "RELATIVE" - first vertex is Pt0, second is Pt0+vt1, third vertex is Pt0+vt2 and fourth is Pt0+vt3, still in anti-clockwise order.

Definition at line 50 of file G4QuadrangularFacet.cc.

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

Referenced by G4QuadrangularFacet(), GetClone(), operator=(), and ~G4QuadrangularFacet().

◆ ~G4QuadrangularFacet()

G4QuadrangularFacet::~G4QuadrangularFacet ( )
overridedefault

Default Destructor.

◆ G4QuadrangularFacet() [2/2]

G4QuadrangularFacet::G4QuadrangularFacet ( const G4QuadrangularFacet & right)

Copy constructor and assignment operator.

Definition at line 208 of file G4QuadrangularFacet.cc.

209 : G4VFacet(rhs)
210{
211 fFacet1 = rhs.fFacet1;
212 fFacet2 = rhs.fFacet2;
213 fRadius = 0.0;
214}

Member Function Documentation

◆ Distance() [1/3]

G4ThreeVector G4QuadrangularFacet::Distance ( const G4ThreeVector & p)

Determines the vector between p and the closest point on the facet to p.

Definition at line 241 of file G4QuadrangularFacet.cc.

242{
243 G4ThreeVector v1 = fFacet1.Distance(p);
244 G4ThreeVector v2 = fFacet2.Distance(p);
245
246 if (v1.mag2() < v2.mag2())
247 {
248 return v1;
249 }
250 return v2;
251}

Referenced by Distance(), and Distance().

◆ Distance() [2/3]

G4double G4QuadrangularFacet::Distance ( const G4ThreeVector & p,
G4double minDist )
overridevirtual

Determines the closest distance between point p and the facet.

Implements G4VFacet.

Definition at line 255 of file G4QuadrangularFacet.cc.

257{
258 G4double dist = Distance(p).mag();
259 return dist;
260}
G4ThreeVector Distance(const G4ThreeVector &p)

◆ Distance() [3/3]

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

Determines the distance to point 'p'. kInfinity is returned if either: (1) outgoing is TRUE and the dot product of the normal vector to the facet and the displacement vector from p to the triangle is negative. (2) outgoing is FALSE and the dot product of the normal vector to the facet and the displacement vector from p to the triangle is positive.

Implements G4VFacet.

Definition at line 264 of file G4QuadrangularFacet.cc.

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

◆ Extent()

G4double G4QuadrangularFacet::Extent ( const G4ThreeVector axis)
overridevirtual

Calculates the furthest the quadrangle extends in fA particular direction defined by the vector axis.

Implements G4VFacet.

Definition at line 285 of file G4QuadrangularFacet.cc.

286{
287 G4double ss = 0;
288
289 for (G4int i = 0; i <= 3; ++i)
290 {
292 if (sp > ss) { ss = sp; }
293 }
294 return ss;
295}
int G4int
Definition G4Types.hh:85
const axis_t axis_to_type< N >::axis
Definition pugixml.cc:9668

◆ GetArea()

G4double G4QuadrangularFacet::GetArea ( ) const
overridevirtual

Auxiliary method for returning the surface area.

Implements G4VFacet.

Definition at line 337 of file G4QuadrangularFacet.cc.

338{
339 G4double area = fFacet1.GetArea() + fFacet2.GetArea();
340 return area;
341}

◆ GetCircumcentre()

G4ThreeVector G4QuadrangularFacet::GetCircumcentre ( ) const
inlineoverridevirtual

Returns the circumcentre point of the facet.

Implements G4VFacet.

◆ GetClone()

G4VFacet * G4QuadrangularFacet::GetClone ( )
overridevirtual

Returns a pointer to a newly allocated duplicate copy of the facet.

Implements G4VFacet.

Definition at line 232 of file G4QuadrangularFacet.cc.

233{
234 auto c = new G4QuadrangularFacet (GetVertex(0), GetVertex(1),
236 return c;
237}
G4QuadrangularFacet(const G4ThreeVector &Pt0, const G4ThreeVector &vt1, const G4ThreeVector &vt2, const G4ThreeVector &vt3, G4FacetVertexType vType)

◆ GetEntityType()

G4String G4QuadrangularFacet::GetEntityType ( ) const
overridevirtual

Returns the type ID, "G4QuadrangularFacet" of the facet.

Implements G4VFacet.

Definition at line 345 of file G4QuadrangularFacet.cc.

346{
347 return "G4QuadrangularFacet";
348}

◆ GetNumberOfVertices()

G4int G4QuadrangularFacet::GetNumberOfVertices ( ) const
inlineoverridevirtual

Returns the number of vertices, i.e. 4.

Implements G4VFacet.

◆ GetPointOnFace()

G4ThreeVector G4QuadrangularFacet::GetPointOnFace ( ) const
overridevirtual

Auxiliary method to get a uniform random point on the facet.

Implements G4VFacet.

Definition at line 324 of file G4QuadrangularFacet.cc.

325{
326 G4double s1 = fFacet1.GetArea();
327 G4double s2 = fFacet2.GetArea();
328 return ((s1 + s2)*G4QuickRand() < s1)
329 ? fFacet1.GetPointOnFace()
330 : fFacet2.GetPointOnFace();
331}
G4double G4QuickRand(uint32_t seed=0)

◆ GetRadius()

G4double G4QuadrangularFacet::GetRadius ( ) const
inlineoverridevirtual

Returns the radius to the anchor point and centered on the circumcentre.

Implements G4VFacet.

◆ GetSurfaceNormal()

G4ThreeVector G4QuadrangularFacet::GetSurfaceNormal ( ) const
overridevirtual

Returns the normal vector to the face.

Implements G4VFacet.

Definition at line 352 of file G4QuadrangularFacet.cc.

353{
354 return fFacet1.GetSurfaceNormal();
355}

Referenced by Distance().

◆ GetVertex()

G4ThreeVector G4QuadrangularFacet::GetVertex ( G4int i) const
inlineoverridevirtual

Returns the vertex based on the index 'i'.

Implements G4VFacet.

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

◆ Intersect()

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

Finds the next intersection when going from 'p' in the direction of 'v'. If 'outgoing' is true, only consider the face if we are going out through the face; otherwise, if false, only consider the face if we are going in through the face.

Returns
true if there is an intersection, false otherwise.

Implements G4VFacet.

Definition at line 299 of file G4QuadrangularFacet.cc.

305{
306 G4bool intersect =
307 fFacet1.Intersect(p,v,outgoing,distance,distFromSurface,normal);
308 if (!intersect)
309 {
310 intersect = fFacet2.Intersect(p,v,outgoing,distance,distFromSurface,normal);
311 }
312 if (!intersect)
313 {
314 distance = distFromSurface = kInfinity;
315 normal.set(0,0,0);
316 }
317 return intersect;
318}
bool G4bool
Definition G4Types.hh:86
void set(double x, double y, double z)

◆ IsDefined()

G4bool G4QuadrangularFacet::IsDefined ( ) const
inlineoverridevirtual

Returns true if the facet is defined.

Implements G4VFacet.

◆ operator=()

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

Definition at line 219 of file G4QuadrangularFacet.cc.

220{
221 if (this == &rhs) { return *this; }
222
223 fFacet1 = rhs.fFacet1;
224 fFacet2 = rhs.fFacet2;
225 fRadius = 0.0;
226
227 return *this;
228}

◆ SetVertex()

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

Methods to set the vertices.

Implements G4VFacet.

Referenced by G4QuadrangularFacet().

◆ SetVertices()

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

Implements G4VFacet.


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