Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4QuadrangularFacet.cc
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 and of QinetiQ Ltd, *
20// * subject to DEFCON 705 IPR conditions. *
21// * By using, copying, modifying or distributing the software (or *
22// * any work based on the software) you agree to acknowledge its *
23// * use in resulting scientific publications, and indicate your *
24// * acceptance of all terms of the Geant4 Software license. *
25// ********************************************************************
26//
27// G4QuadrangularFacet class implementation.
28//
29// 31 October 2004 P R Truscott, QinetiQ Ltd, UK - Created.
30// 12 October 2012 M Gayer, CERN
31// New implementation reducing memory requirements by 50%,
32// and considerable CPU speedup together with the new
33// implementation of G4TessellatedSolid.
34// 29 February 2016 E Tcherniaev, CERN
35// Added exhaustive tests to catch various problems with a
36// quadrangular facet: collinear vertices, non planar surface,
37// degenerate, concave or self intersecting quadrilateral.
38// --------------------------------------------------------------------
39
41#include "G4QuickRand.hh"
42
43using namespace std;
44
45///////////////////////////////////////////////////////////////////////////////
46//
47// Constructing two adjacent G4TriangularFacet
48// Not efficient, but practical...
49//
51 const G4ThreeVector& vt1,
52 const G4ThreeVector& vt2,
53 const G4ThreeVector& vt3,
54 G4FacetVertexType vertexType)
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 //
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}
205
206///////////////////////////////////////////////////////////////////////////////
207//
209 : G4VFacet(rhs)
210{
211 fFacet1 = rhs.fFacet1;
212 fFacet2 = rhs.fFacet2;
213 fRadius = 0.0;
214}
215
216///////////////////////////////////////////////////////////////////////////////
217//
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}
229
230///////////////////////////////////////////////////////////////////////////////
231//
233{
234 auto c = new G4QuadrangularFacet (GetVertex(0), GetVertex(1),
236 return c;
237}
238
239///////////////////////////////////////////////////////////////////////////////
240//
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}
252
253///////////////////////////////////////////////////////////////////////////////
254//
256 G4double)
257{
258 G4double dist = Distance(p).mag();
259 return dist;
260}
261
262///////////////////////////////////////////////////////////////////////////////
263//
265 const G4bool outgoing)
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}
282
283///////////////////////////////////////////////////////////////////////////////
284//
286{
287 G4double ss = 0;
288
289 for (G4int i = 0; i <= 3; ++i)
290 {
291 G4double sp = GetVertex(i).dot(axis);
292 if (sp > ss) { ss = sp; }
293 }
294 return ss;
295}
296
297///////////////////////////////////////////////////////////////////////////////
298//
300 const G4ThreeVector& v,
301 G4bool outgoing,
302 G4double& distance,
303 G4double& distFromSurface,
304 G4ThreeVector& normal)
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}
319
320///////////////////////////////////////////////////////////////////////////////
321//
322// Auxiliary method to get a uniform random point on the facet
323//
325{
326 G4double s1 = fFacet1.GetArea();
327 G4double s2 = fFacet2.GetArea();
328 return ((s1 + s2)*G4QuickRand() < s1)
329 ? fFacet1.GetPointOnFace()
330 : fFacet2.GetPointOnFace();
331}
332
333///////////////////////////////////////////////////////////////////////////////
334//
335// Auxiliary method for returning the surface area
336//
338{
339 G4double area = fFacet1.GetArea() + fFacet2.GetArea();
340 return area;
341}
342
343///////////////////////////////////////////////////////////////////////////////
344//
346{
347 return "G4QuadrangularFacet";
348}
349
350///////////////////////////////////////////////////////////////////////////////
351//
353{
354 return fFacet1.GetSurfaceNormal();
355}
G4double epsilon(G4double density, G4double temperature)
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4double G4QuickRand(uint32_t seed=0)
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4FacetVertexType
Definition G4VFacet.hh:48
@ ABSOLUTE
Definition G4VFacet.hh:48
#define G4endl
Definition G4ios.hh:67
Hep3Vector unit() const
double mag2() const
double dot(const Hep3Vector &) const
double mag() const
void set(double x, double y, double z)
G4QuadrangularFacet defines a facet with 4 vertices, used for the contruction of G4TessellatedSolid....
G4double GetArea() const override
G4ThreeVector Distance(const G4ThreeVector &p)
G4ThreeVector GetSurfaceNormal() const override
G4QuadrangularFacet & operator=(const G4QuadrangularFacet &right)
G4ThreeVector GetVertex(G4int i) const override
G4QuadrangularFacet(const G4ThreeVector &Pt0, const G4ThreeVector &vt1, const G4ThreeVector &vt2, const G4ThreeVector &vt3, G4FacetVertexType vType)
G4VFacet * GetClone() override
G4GeometryType GetEntityType() const override
void SetVertex(G4int i, const G4ThreeVector &val) override
G4double Extent(const G4ThreeVector axis) override
G4ThreeVector GetPointOnFace() const override
G4bool Intersect(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool outgoing, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal) override
G4TriangularFacet defines a facet with 3 vertices, used for the contruction of G4TessellatedSolid....
static const G4double dirTolerance
Definition G4VFacet.hh:184
G4double kCarTolerance
Definition G4VFacet.hh:185
const axis_t axis_to_type< N >::axis
Definition pugixml.cc:9668