Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TriangularFacet.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. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// G4TriangularFacet implementation
27//
28// 31.10.2004, P R Truscott, QinetiQ Ltd, UK - Created.
29// 01.08.2007, P R Truscott, QinetiQ Ltd, UK
30// Significant modification to correct for errors and enhance
31// based on patches/observations kindly provided by Rickard
32// Holmberg.
33// 12.10.2012, M Gayer, CERN
34// New implementation reducing memory requirements by 50%,
35// and considerable CPU speedup together with the new
36// implementation of G4TessellatedSolid.
37// 23.02.2016, E Tcherniaev, CERN
38// Improved test to detect degenerate (too small or
39// too narrow) triangles.
40// --------------------------------------------------------------------
41
42#include "G4TriangularFacet.hh"
43#include "G4QuickRand.hh"
45
46using namespace std;
47
48///////////////////////////////////////////////////////////////////////////////
49//
50// Definition of triangular facet using absolute vectors to fVertices.
51// From this for first vector is retained to define the facet location and
52// two relative vectors (E0 and E1) define the sides and orientation of
53// the outward surface normal.
54//
56 const G4ThreeVector& vt1,
57 const G4ThreeVector& vt2,
58 G4FacetVertexType vertexType)
59{
60 fVertices = new vector<G4ThreeVector>(3);
61
62 SetVertex(0, vt0);
63 if (vertexType == ABSOLUTE)
64 {
65 SetVertex(1, vt1);
66 SetVertex(2, vt2);
67 fE1 = vt1 - vt0;
68 fE2 = vt2 - vt0;
69 }
70 else
71 {
72 SetVertex(1, vt0 + vt1);
73 SetVertex(2, vt0 + vt2);
74 fE1 = vt1;
75 fE2 = vt2;
76 }
77
78 G4ThreeVector E1xE2 = fE1.cross(fE2);
79 fArea = 0.5 * E1xE2.mag();
80 for (G4int i = 0; i < 3; ++i)
81 {
82 fIndices[i] = -1;
83 }
84
85 fIsDefined = true;
86 G4double delta = kCarTolerance; // Set tolerance for checking
87
88 // Check length of edges
89 //
90 G4double leng1 = fE1.mag();
91 G4double leng2 = (fE2-fE1).mag();
92 G4double leng3 = fE2.mag();
93 if (leng1 <= delta || leng2 <= delta || leng3 <= delta)
94 {
95 fIsDefined = false;
96 }
97
98 // Check min height of triangle
99 //
100 if (fIsDefined)
101 {
102 if (2.*fArea/std::max(std::max(leng1,leng2),leng3) <= delta)
103 {
104 fIsDefined = false;
105 }
106 }
107
108 // Define facet
109 //
110 if (!fIsDefined)
111 {
112 ostringstream message;
113 message << "Facet is too small or too narrow." << G4endl
114 << "Triangle area = " << fArea << G4endl
115 << "P0 = " << GetVertex(0) << G4endl
116 << "P1 = " << GetVertex(1) << G4endl
117 << "P2 = " << GetVertex(2) << G4endl
118 << "Side1 length (P0->P1) = " << leng1 << G4endl
119 << "Side2 length (P1->P2) = " << leng2 << G4endl
120 << "Side3 length (P2->P0) = " << leng3;
121 G4Exception("G4TriangularFacet::G4TriangularFacet()",
122 "GeomSolids1001", JustWarning, message);
123 fSurfaceNormal.set(0,0,0);
124 fA = fB = fC = 0.0;
125 fDet = 0.0;
126 fCircumcentre = vt0 + 0.5*fE1 + 0.5*fE2;
127 fArea = fRadius = 0.0;
128 }
129 else
130 {
131 fSurfaceNormal = E1xE2.unit();
132 fA = fE1.mag2();
133 fB = fE1.dot(fE2);
134 fC = fE2.mag2();
135 fDet = std::fabs(fA*fC - fB*fB);
136
137 fCircumcentre =
138 vt0 + (E1xE2.cross(fE1)*fC + fE2.cross(E1xE2)*fA) / (2.*E1xE2.mag2());
139 fRadius = (fCircumcentre - vt0).mag();
140 }
141}
142
143///////////////////////////////////////////////////////////////////////////////
144//
146{
147 fVertices = new vector<G4ThreeVector>(3);
148 G4ThreeVector zero(0,0,0);
149 SetVertex(0, zero);
150 SetVertex(1, zero);
151 SetVertex(2, zero);
152 for (G4int i = 0; i < 3; ++i)
153 {
154 fIndices[i] = -1;
155 }
156 fIsDefined = false;
157 fSurfaceNormal.set(0,0,0);
158 fA = fB = fC = 0;
159 fE1 = zero;
160 fE2 = zero;
161 fDet = 0.0;
162 fArea = fRadius = 0.0;
163}
164
165///////////////////////////////////////////////////////////////////////////////
166//
171
172///////////////////////////////////////////////////////////////////////////////
173//
174void G4TriangularFacet::CopyFrom (const G4TriangularFacet& rhs)
175{
176 auto p = (char *) &rhs;
177 copy(p, p + sizeof(*this), (char *)this);
178
179 if (fIndices[0] < 0 && fVertices == nullptr)
180 {
181 fVertices = new vector<G4ThreeVector>(3);
182 for (G4int i = 0; i < 3; ++i)
183 {
184 (*fVertices)[i] = (*rhs.fVertices)[i];
185 }
186 }
187}
188
189///////////////////////////////////////////////////////////////////////////////
190//
191void G4TriangularFacet::MoveFrom (G4TriangularFacet& rhs)
192{
193 fSurfaceNormal = std::move(rhs.fSurfaceNormal);
194 fArea = rhs.fArea;
195 fCircumcentre = std::move(rhs.fCircumcentre);
196 fRadius = rhs.fRadius;
197 fIndices = rhs.fIndices;
198 fA = rhs.fA; fB = rhs.fB; fC = rhs.fC;
199 fDet = rhs.fDet;
200 fSqrDist = rhs.fSqrDist;
201 fE1 = std::move(rhs.fE1); fE2 = std::move(rhs.fE2);
202 fIsDefined = rhs.fIsDefined;
203 fVertices = rhs.fVertices;
204 rhs.fVertices = nullptr;
205}
206
207///////////////////////////////////////////////////////////////////////////////
208//
210 : G4VFacet(rhs)
211{
212 CopyFrom(rhs);
213}
214
215///////////////////////////////////////////////////////////////////////////////
216//
218 : G4VFacet(rhs)
219{
220 MoveFrom(rhs);
221}
222
223///////////////////////////////////////////////////////////////////////////////
224//
227{
228 SetVertices(nullptr);
229
230 if (this != &rhs)
231 {
232 delete fVertices;
233 CopyFrom(rhs);
234 }
235
236 return *this;
237}
238
239///////////////////////////////////////////////////////////////////////////////
240//
243{
244 SetVertices(nullptr);
245
246 if (this != &rhs)
247 {
248 delete fVertices;
249 MoveFrom(rhs);
250 }
251
252 return *this;
253}
254
255///////////////////////////////////////////////////////////////////////////////
256//
257// GetClone
258//
259// Simple member function to generate fA duplicate of the triangular facet.
260//
262{
263 auto fc = new G4TriangularFacet (GetVertex(0), GetVertex(1),
264 GetVertex(2), ABSOLUTE);
265 return fc;
266}
267
268///////////////////////////////////////////////////////////////////////////////
269//
270// GetFlippedFacet
271//
272// Member function to generate an identical facet, but with the normal vector
273// pointing at 180 degrees.
274//
276{
277 auto flipped = new G4TriangularFacet (GetVertex(0), GetVertex(1),
278 GetVertex(2), ABSOLUTE);
279 return flipped;
280}
281
282///////////////////////////////////////////////////////////////////////////////
283//
284// Distance (G4ThreeVector)
285//
286// Determines the vector between p and the closest point on the facet to p.
287// This is based on the algorithm published in "Geometric Tools for Computer
288// Graphics," Philip J Scheider and David H Eberly, Elsevier Science (USA),
289// 2003. at the time of writing, the algorithm is also available in fA
290// technical note "Distance between point and triangle in 3D," by David Eberly
291// at http://www.geometrictools.com/Documentation/DistancePoint3Triangle3.pdf
292//
293// The by-product is the square-distance fSqrDist, which is retained
294// in case needed by the other "Distance" member functions.
295//
297{
298 G4ThreeVector D = GetVertex(0) - p;
299 G4double d = fE1.dot(D);
300 G4double e = fE2.dot(D);
301 G4double f = D.mag2();
302 G4double q = fB*e - fC*d;
303 G4double t = fB*d - fA*e;
304 fSqrDist = 0.;
305
306 if (q+t <= fDet)
307 {
308 if (q < 0.0)
309 {
310 if (t < 0.0)
311 {
312 //
313 // We are in region 4.
314 //
315 if (d < 0.0)
316 {
317 t = 0.0;
318 if (-d >= fA) {q = 1.0; fSqrDist = fA + 2.0*d + f;}
319 else {q = -d/fA; fSqrDist = d*q + f;}
320 }
321 else
322 {
323 q = 0.0;
324 if (e >= 0.0) {t = 0.0; fSqrDist = f;}
325 else if (-e >= fC) {t = 1.0; fSqrDist = fC + 2.0*e + f;}
326 else {t = -e/fC; fSqrDist = e*t + f;}
327 }
328 }
329 else
330 {
331 //
332 // We are in region 3.
333 //
334 q = 0.0;
335 if (e >= 0.0) {t = 0.0; fSqrDist = f;}
336 else if (-e >= fC) {t = 1.0; fSqrDist = fC + 2.0*e + f;}
337 else {t = -e/fC; fSqrDist = e*t + f;}
338 }
339 }
340 else if (t < 0.0)
341 {
342 //
343 // We are in region 5.
344 //
345 t = 0.0;
346 if (d >= 0.0) {q = 0.0; fSqrDist = f;}
347 else if (-d >= fA) {q = 1.0; fSqrDist = fA + 2.0*d + f;}
348 else {q = -d/fA; fSqrDist = d*q + f;}
349 }
350 else
351 {
352 //
353 // We are in region 0.
354 //
355 G4double dist = fSurfaceNormal.dot(D);
356 fSqrDist = dist*dist;
357 return fSurfaceNormal*dist;
358 }
359 }
360 else
361 {
362 if (q < 0.0)
363 {
364 //
365 // We are in region 2.
366 //
367 G4double tmp0 = fB + d;
368 G4double tmp1 = fC + e;
369 if (tmp1 > tmp0)
370 {
371 G4double numer = tmp1 - tmp0;
372 G4double denom = fA - 2.0*fB + fC;
373 if (numer >= denom) {q = 1.0; t = 0.0; fSqrDist = fA + 2.0*d + f;}
374 else
375 {
376 q = numer/denom;
377 t = 1.0 - q;
378 fSqrDist = q*(fA*q + fB*t +2.0*d) + t*(fB*q + fC*t + 2.0*e) + f;
379 }
380 }
381 else
382 {
383 q = 0.0;
384 if (tmp1 <= 0.0) {t = 1.0; fSqrDist = fC + 2.0*e + f;}
385 else if (e >= 0.0) {t = 0.0; fSqrDist = f;}
386 else {t = -e/fC; fSqrDist = e*t + f;}
387 }
388 }
389 else if (t < 0.0)
390 {
391 //
392 // We are in region 6.
393 //
394 G4double tmp0 = fB + e;
395 G4double tmp1 = fA + d;
396 if (tmp1 > tmp0)
397 {
398 G4double numer = tmp1 - tmp0;
399 G4double denom = fA - 2.0*fB + fC;
400 if (numer >= denom) {t = 1.0; q = 0.0; fSqrDist = fC + 2.0*e + f;}
401 else
402 {
403 t = numer/denom;
404 q = 1.0 - t;
405 fSqrDist = q*(fA*q + fB*t +2.0*d) + t*(fB*q + fC*t + 2.0*e) + f;
406 }
407 }
408 else
409 {
410 t = 0.0;
411 if (tmp1 <= 0.0) {q = 1.0; fSqrDist = fA + 2.0*d + f;}
412 else if (d >= 0.0) {q = 0.0; fSqrDist = f;}
413 else {q = -d/fA; fSqrDist = d*q + f;}
414 }
415 }
416 else
417 //
418 // We are in region 1.
419 //
420 {
421 G4double numer = fC + e - fB - d;
422 if (numer <= 0.0)
423 {
424 q = 0.0;
425 t = 1.0;
426 fSqrDist = fC + 2.0*e + f;
427 }
428 else
429 {
430 G4double denom = fA - 2.0*fB + fC;
431 if (numer >= denom) {q = 1.0; t = 0.0; fSqrDist = fA + 2.0*d + f;}
432 else
433 {
434 q = numer/denom;
435 t = 1.0 - q;
436 fSqrDist = q*(fA*q + fB*t + 2.0*d) + t*(fB*q + fC*t + 2.0*e) + f;
437 }
438 }
439 }
440 }
441 //
442 //
443 // Do fA check for rounding errors in the distance-squared. It appears that
444 // the conventional methods for calculating fSqrDist breaks down when very
445 // near to or at the surface (as required by transport).
446 // We'll therefore also use the magnitude-squared of the vector displacement.
447 // (Note that I've also tried to get around this problem by using the
448 // existing equations for
449 //
450 // fSqrDist = function(fA,fB,fC,d,q,t)
451 //
452 // and use fA more accurate addition process which minimises errors and
453 // breakdown of cummutitivity [where (A+B)+C != A+(B+C)] but this still
454 // doesn't work.
455 // Calculation from u = D + q*fE1 + t*fE2 is less efficient, but appears
456 // more robust.
457 //
458 if (fSqrDist < 0.0) { fSqrDist = 0.; }
459 G4ThreeVector u = D + q*fE1 + t*fE2;
460 G4double u2 = u.mag2();
461 //
462 // The following (part of the roundoff error check) is from Oliver Merle'q
463 // updates.
464 //
465 if (fSqrDist > u2) { fSqrDist = u2; }
466
467 return u;
468}
469
470///////////////////////////////////////////////////////////////////////////////
471//
472// Distance (G4ThreeVector, G4double)
473//
474// Determines the closest distance between point p and the facet. This makes
475// use of G4ThreeVector G4TriangularFacet::Distance, which stores the
476// square of the distance in variable fSqrDist. If approximate methods show
477// the distance is to be greater than minDist, then forget about further
478// computation and return fA very large number.
479//
481 G4double minDist)
482{
483 //
484 // Start with quicky test to determine if the surface of the sphere enclosing
485 // the triangle is any closer to p than minDist. If not, then don't bother
486 // about more accurate test.
487 //
488 G4double dist = kInfinity;
489 if ((p-fCircumcentre).mag()-fRadius < minDist)
490 {
491 //
492 // It's possible that the triangle is closer than minDist,
493 // so do more accurate assessment.
494 //
495 dist = Distance(p).mag();
496 }
497 return dist;
498}
499
500///////////////////////////////////////////////////////////////////////////////
501//
502// Distance (G4ThreeVector, G4double, G4bool)
503//
504// Determine the distance to point p. kInfinity is returned if either:
505// (1) outgoing is TRUE and the dot product of the normal vector to the facet
506// and the displacement vector from p to the triangle is negative.
507// (2) outgoing is FALSE and the dot product of the normal vector to the facet
508// and the displacement vector from p to the triangle is positive.
509// If approximate methods show the distance is to be greater than minDist, then
510// forget about further computation and return fA very large number.
511//
512// This method has been heavily modified thanks to the valuable comments and
513// corrections of Rickard Holmberg.
514//
516 G4double minDist,
517 const G4bool outgoing)
518{
519 //
520 // Start with quicky test to determine if the surface of the sphere enclosing
521 // the triangle is any closer to p than minDist. If not, then don't bother
522 // about more accurate test.
523 //
524 G4double dist = kInfinity;
525 if ((p-fCircumcentre).mag()-fRadius < minDist)
526 {
527 //
528 // It's possible that the triangle is closer than minDist,
529 // so do more accurate assessment.
530 //
531 G4ThreeVector v = Distance(p);
532 G4double dist1 = sqrt(fSqrDist);
533 G4double dir = v.dot(fSurfaceNormal);
534 G4bool wrongSide = (dir > 0.0 && !outgoing) || (dir < 0.0 && outgoing);
535 if (dist1 <= kCarTolerance)
536 {
537 //
538 // Point p is very close to triangle. Check if it's on the wrong side,
539 // in which case return distance of 0.0 otherwise .
540 //
541 if (wrongSide)
542 {
543 dist = 0.0;
544 }
545 else
546 {
547 dist = dist1;
548 }
549 }
550 else if (!wrongSide)
551 {
552 dist = dist1;
553 }
554 }
555 return dist;
556}
557
558///////////////////////////////////////////////////////////////////////////////
559//
560// Extent
561//
562// Calculates the furthest the triangle extends in fA particular direction
563// defined by the vector axis.
564//
566{
567 G4double ss = GetVertex(0).dot(axis);
568 G4double sp = GetVertex(1).dot(axis);
569 if (sp > ss) { ss = sp; }
570 sp = GetVertex(2).dot(axis);
571 if (sp > ss) { ss = sp; }
572 return ss;
573}
574
575///////////////////////////////////////////////////////////////////////////////
576//
577// Intersect
578//
579// Member function to find the next intersection when going from p in the
580// direction of v. If:
581// (1) "outgoing" is TRUE, only consider the face if we are going out through
582// the face.
583// (2) "outgoing" is FALSE, only consider the face if we are going in through
584// the face.
585// Member functions returns TRUE if there is an intersection, FALSE otherwise.
586// Sets the distance (distance along w), distFromSurface (orthogonal distance)
587// and normal.
588//
589// Also considers intersections that happen with negative distance for small
590// distances of distFromSurface = 0.5*kCarTolerance in the wrong direction.
591// This is to detect kSurface without doing fA full Inside(p) in
592// G4TessellatedSolid::Distance(p,v) calculation.
593//
594// This member function is thanks the valuable work of Rickard Holmberg. PT.
595// However, "gotos" are the Work of the Devil have been exorcised with
596// extreme prejudice!!
597//
598// IMPORTANT NOTE: These calculations are predicated on v being fA unit
599// vector. If G4TessellatedSolid or other classes call this member function
600// with |v| != 1 then there will be errors.
601//
603 const G4ThreeVector& v,
604 G4bool outgoing,
605 G4double& distance,
606 G4double& distFromSurface,
607 G4ThreeVector& normal)
608{
609 //
610 // Check whether the direction of the facet is consistent with the vector v
611 // and the need to be outgoing or ingoing. If inconsistent, disregard and
612 // return false.
613 //
614 G4double w = v.dot(fSurfaceNormal);
615 if ((outgoing && w < -dirTolerance) || (!outgoing && w > dirTolerance))
616 {
617 distance = kInfinity;
618 distFromSurface = kInfinity;
619 normal.set(0,0,0);
620 return false;
621 }
622 //
623 // Calculate the orthogonal distance from p to the surface containing the
624 // triangle. Then determine if we're on the right or wrong side of the
625 // surface (at fA distance greater than kCarTolerance to be consistent with
626 // "outgoing".
627 //
628 const G4ThreeVector& p0 = GetVertex(0);
629 G4ThreeVector D = p0 - p;
630 distFromSurface = D.dot(fSurfaceNormal);
631 G4bool wrongSide = (outgoing && distFromSurface < -0.5*kCarTolerance) ||
632 (!outgoing && distFromSurface > 0.5*kCarTolerance);
633
634 if (wrongSide)
635 {
636 distance = kInfinity;
637 distFromSurface = kInfinity;
638 normal.set(0,0,0);
639 return false;
640 }
641
642 wrongSide = (outgoing && distFromSurface < 0.0)
643 || (!outgoing && distFromSurface > 0.0);
644 if (wrongSide)
645 {
646 //
647 // We're slightly on the wrong side of the surface. Check if we're close
648 // enough using fA precise distance calculation.
649 //
650 G4ThreeVector u = Distance(p);
651 if (fSqrDist <= kCarTolerance*kCarTolerance)
652 {
653 //
654 // We're very close. Therefore return fA small negative number
655 // to pretend we intersect.
656 //
657 // distance = -0.5*kCarTolerance
658 distance = 0.0;
659 normal = fSurfaceNormal;
660 return true;
661 }
662
663 //
664 // We're close to the surface containing the triangle, but sufficiently
665 // far from the triangle, and on the wrong side compared to the directions
666 // of the surface normal and v. There is no intersection.
667 //
668 distance = kInfinity;
669 distFromSurface = kInfinity;
670 normal.set(0,0,0);
671 return false;
672 }
673 if (w < dirTolerance && w > -dirTolerance)
674 {
675 //
676 // The ray is within the plane of the triangle. Project the problem into 2D
677 // in the plane of the triangle. First try to create orthogonal unit vectors
678 // mu and nu, where mu is fE1/|fE1|. This is kinda like
679 // the original algorithm due to Rickard Holmberg, but with better
680 // mathematical justification than the original method ... however,
681 // beware Rickard's was less time-consuming.
682 //
683 // Note that vprime is not fA unit vector. We need to keep it unnormalised
684 // since the values of distance along vprime (s0 and s1) for intersection
685 // with the triangle will be used to determine if we cut the plane at the
686 // same time.
687 //
688 G4ThreeVector mu = fE1.unit();
689 G4ThreeVector nu = fSurfaceNormal.cross(mu);
690 G4TwoVector pprime(p.dot(mu), p.dot(nu));
691 G4TwoVector vprime(v.dot(mu), v.dot(nu));
692 G4TwoVector P0prime(p0.dot(mu), p0.dot(nu));
693 G4TwoVector E0prime(fE1.mag(), 0.0);
694 G4TwoVector E1prime(fE2.dot(mu), fE2.dot(nu));
695 G4TwoVector loc[2];
697 vprime, P0prime, E0prime, E1prime, loc))
698 {
699 //
700 // There is an intersection between the line and triangle in 2D.
701 // Now check which part of the line intersects with the plane
702 // containing the triangle in 3D.
703 //
704 G4double vprimemag = vprime.mag();
705 G4double s0 = (loc[0] - pprime).mag()/vprimemag;
706 G4double s1 = (loc[1] - pprime).mag()/vprimemag;
707 G4double normDist0 = fSurfaceNormal.dot(s0*v) - distFromSurface;
708 G4double normDist1 = fSurfaceNormal.dot(s1*v) - distFromSurface;
709
710 if ((normDist0 < 0.0 && normDist1 < 0.0)
711 || (normDist0 > 0.0 && normDist1 > 0.0)
712 || (normDist0 == 0.0 && normDist1 == 0.0) )
713 {
714 distance = kInfinity;
715 distFromSurface = kInfinity;
716 normal.set(0,0,0);
717 return false;
718 }
719
720 G4double dnormDist = normDist1 - normDist0;
721 if (fabs(dnormDist) < DBL_EPSILON)
722 {
723 distance = s0;
724 normal = fSurfaceNormal;
725 if (!outgoing) { distFromSurface = -distFromSurface; }
726 return true;
727 }
728
729 distance = s0 - normDist0*(s1-s0)/dnormDist;
730 normal = fSurfaceNormal;
731 if (!outgoing) { distFromSurface = -distFromSurface; }
732 return true;
733 }
734 distance = kInfinity;
735 distFromSurface = kInfinity;
736 normal.set(0,0,0);
737 return false;
738 }
739 //
740 //
741 // Use conventional algorithm to determine the whether there is an
742 // intersection. This involves determining the point of intersection of the
743 // line with the plane containing the triangle, and then calculating if the
744 // point is within the triangle.
745 //
746 distance = distFromSurface / w;
747 G4ThreeVector pp = p + v*distance;
748 G4ThreeVector DD = p0 - pp;
749 G4double d = fE1.dot(DD);
750 G4double e = fE2.dot(DD);
751 G4double ss = fB*e - fC*d;
752 G4double t = fB*d - fA*e;
753
754 G4double sTolerance = (fabs(fB)+ fabs(fC) + fabs(d) + fabs(e))*kCarTolerance;
755 G4double tTolerance = (fabs(fA)+ fabs(fB) + fabs(d) + fabs(e))*kCarTolerance;
756 G4double detTolerance = (fabs(fA)+ fabs(fC) + 2*fabs(fB) )*kCarTolerance;
757
758 //if (ss < 0.0 || t < 0.0 || ss+t > fDet)
759 if (ss < -sTolerance || t < -tTolerance || ( ss+t - fDet ) > detTolerance)
760 {
761 //
762 // The intersection is outside of the triangle.
763 //
764 distance = distFromSurface = kInfinity;
765 normal.set(0,0,0);
766 return false;
767 }
768
769 //
770 // There is an intersection. Now we only need to set the surface normal.
771 //
772 normal = fSurfaceNormal;
773 if (!outgoing) { distFromSurface = -distFromSurface; }
774 return true;
775}
776
777////////////////////////////////////////////////////////////////////////
778//
779// GetPointOnFace
780//
781// Auxiliary method, returns a uniform random point on the facet
782//
784{
785 G4double u = G4QuickRand();
786 G4double v = G4QuickRand();
787 if (u + v > 1.)
788 {
789 u = 1. - u;
790 v = 1. - v;
791 }
792 return GetVertex(0) + u*fE1 + v*fE2;
793}
794
795////////////////////////////////////////////////////////////////////////
796//
797// GetArea
798//
799// Auxiliary method for returning the surface fArea
800//
802{
803 return fArea;
804}
805
806////////////////////////////////////////////////////////////////////////
807//
809{
810 return "G4TriangularFacet";
811}
812
813////////////////////////////////////////////////////////////////////////
814//
816{
817 return fSurfaceNormal;
818}
819
820////////////////////////////////////////////////////////////////////////
821//
823{
824 fSurfaceNormal = normal;
825}
G4double D(G4double temp)
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4double G4QuickRand(uint32_t seed=0)
CLHEP::Hep3Vector G4ThreeVector
CLHEP::Hep2Vector G4TwoVector
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
G4String G4GeometryType
Definition G4VSolid.hh:70
#define G4endl
Definition G4ios.hh:67
double mag() const
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)
static G4bool IntersectLineAndTriangle2D(const G4TwoVector &p, const G4TwoVector &v, const G4TwoVector &p0, const G4TwoVector &e0, const G4TwoVector &e1, G4TwoVector location[2])
G4TriangularFacet defines a facet with 3 vertices, used for the contruction of G4TessellatedSolid....
G4GeometryType GetEntityType() const override
G4double Extent(const G4ThreeVector axis) override
void SetVertex(G4int i, const G4ThreeVector &val) override
G4ThreeVector GetSurfaceNormal() const override
G4ThreeVector GetVertex(G4int i) const override
G4bool Intersect(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool outgoing, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal) override
G4TriangularFacet & operator=(const G4TriangularFacet &right)
void SetVertices(std::vector< G4ThreeVector > *v) override
G4ThreeVector GetPointOnFace() const override
G4TriangularFacet * GetFlippedFacet()
G4ThreeVector Distance(const G4ThreeVector &p)
G4VFacet * GetClone() override
void SetSurfaceNormal(const G4ThreeVector &normal)
G4double GetArea() const override
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
#define DBL_EPSILON
Definition templates.hh:66