Geant4 10.7.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 "geomdefs.hh"
42#include "Randomize.hh"
43
44using namespace std;
45
46///////////////////////////////////////////////////////////////////////////////
47//
48// Constructing two adjacent G4TriangularFacet
49// Not efficient, but practical...
50//
52 const G4ThreeVector& vt1,
53 const G4ThreeVector& vt2,
54 const G4ThreeVector& vt3,
55 G4FacetVertexType vertexType)
56 : G4VFacet()
57{
58 G4double delta = 1.0 * kCarTolerance; // dimension tolerance
59 G4double epsilon = 0.01 * kCarTolerance; // planarity tolerance
60
61 G4ThreeVector e1, e2, e3;
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();
194 fFacet1.SetSurfaceNormal(normal);
195 fFacet2.SetSurfaceNormal(normal);
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}
206
207///////////////////////////////////////////////////////////////////////////////
208//
210{
211}
212
213///////////////////////////////////////////////////////////////////////////////
214//
216 : G4VFacet(rhs)
217{
218 fFacet1 = rhs.fFacet1;
219 fFacet2 = rhs.fFacet2;
220 fRadius = 0.0;
221}
222
223///////////////////////////////////////////////////////////////////////////////
224//
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}
236
237///////////////////////////////////////////////////////////////////////////////
238//
240{
242 GetVertex(2), GetVertex(3),
243 ABSOLUTE);
244 return c;
245}
246
247///////////////////////////////////////////////////////////////////////////////
248//
250{
251 G4ThreeVector v1 = fFacet1.Distance(p);
252 G4ThreeVector v2 = fFacet2.Distance(p);
253
254 if (v1.mag2() < v2.mag2()) return v1;
255 else return v2;
256}
257
258///////////////////////////////////////////////////////////////////////////////
259//
261 G4double)
262{
263 G4double dist = Distance(p).mag();
264 return dist;
265}
266
267///////////////////////////////////////////////////////////////////////////////
268//
270 const G4bool outgoing)
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}
283
284///////////////////////////////////////////////////////////////////////////////
285//
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}
297
298///////////////////////////////////////////////////////////////////////////////
299//
301 const G4ThreeVector& v,
302 G4bool outgoing,
303 G4double& distance,
304 G4double& distFromSurface,
305 G4ThreeVector& normal)
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}
318
319///////////////////////////////////////////////////////////////////////////////
320//
321// Auxiliary method to get a uniform random point on the facet
322//
324{
325 G4double s1 = fFacet1.GetArea();
326 G4double s2 = fFacet2.GetArea();
327 return ((s1+s2)*G4UniformRand() < s1) ?
328 fFacet1.GetPointOnFace() : fFacet2.GetPointOnFace();
329}
330
331///////////////////////////////////////////////////////////////////////////////
332//
333// Auxiliary method for returning the surface area
334//
336{
337 G4double area = fFacet1.GetArea() + fFacet2.GetArea();
338 return area;
339}
340
341///////////////////////////////////////////////////////////////////////////////
342//
344{
345 return "G4QuadrangularFacet";
346}
347
348///////////////////////////////////////////////////////////////////////////////
349//
351{
352 return fFacet1.GetSurfaceNormal();
353}
double epsilon(double density, double temperature)
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
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:57
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
double mag2() const
Hep3Vector cross(const Hep3Vector &) const
double dot(const Hep3Vector &) const
double mag() const
void set(double x, double y, double z)
G4ThreeVector Distance(const G4ThreeVector &p)
G4bool Intersect(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool outgoing, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal)
G4QuadrangularFacet(const G4ThreeVector &Pt0, const G4ThreeVector &vt1, const G4ThreeVector &vt2, const G4ThreeVector &vt3, G4FacetVertexType)
void SetVertex(G4int i, const G4ThreeVector &val)
G4ThreeVector GetSurfaceNormal() const
G4QuadrangularFacet & operator=(const G4QuadrangularFacet &right)
G4ThreeVector GetVertex(G4int i) const
G4double GetArea() const
G4GeometryType GetEntityType() const
G4double Extent(const G4ThreeVector axis)
G4ThreeVector GetPointOnFace() const
G4double GetArea() const
void SetSurfaceNormal(G4ThreeVector normal)
G4ThreeVector GetPointOnFace() const
G4bool Intersect(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool outgoing, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal)
G4ThreeVector Distance(const G4ThreeVector &p)
G4ThreeVector GetSurfaceNormal() const
static const G4double dirTolerance
Definition: G4VFacet.hh:92
G4double kCarTolerance
Definition: G4VFacet.hh:93