Loading...
Searching...
No Matches
triangleI.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2018-2024 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
29#include "IOstreams.H"
30#include "pointHit.H"
32
33// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34
36(
37 const point& p0,
38 const point& p1,
39 const point& p2
40)
42 a() = p0;
43 b() = p1;
44 c() = p2;
45}
46
47
50 a() = pts.a();
51 b() = pts.b();
52 c() = pts.c();
53}
54
63(
64 const UList<point>& points,
65 const FixedList<label, 3>& indices
66)
67:
68 FixedList<point, 3>(points, indices)
69{}
70
71
73(
74 const UList<point>& points,
75 const label p0,
76 const label p1,
77 const label p2
78)
79{
80 a() = points[p0];
81 b() = points[p1];
82 c() = points[p2];
83}
84
85
86template<class Point, class PointRef>
88(
89 const Point& p0,
90 const Point& p1,
91 const Point& p2
92)
93:
94 a_(p0),
95 b_(p1),
96 c_(p2)
97{}
98
99
100template<class Point, class PointRef>
102(
104)
105:
106 a_(pts.template get<0>()),
107 b_(pts.template get<1>()),
108 c_(pts.template get<2>())
109{}
110
111
112template<class Point, class PointRef>
114(
115 const UList<Point>& points,
116 const FixedList<label, 3>& indices
117)
118:
119 a_(points[indices.template get<0>()]),
120 b_(points[indices.template get<1>()]),
121 c_(points[indices.template get<2>()])
122{}
123
124
125template<class Point, class PointRef>
127(
128 const UList<Point>& points,
129 const label p0,
130 const label p1,
131 const label p2
132)
133:
134 a_(points[p0]),
135 b_(points[p1]),
136 c_(points[p2])
137{}
138
139
140template<class Point, class PointRef>
142{
143 is >> *this;
144}
145
146
147// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
148
149template<class Point, class PointRef>
151(
152 const Point& p0,
153 const Point& p1,
154 const Point& p2
156{
157 return (1.0/3.0)*(p0 + p1 + p2);
158}
159
160
161template<class Point, class PointRef>
163{
164 return (1.0/3.0)*(a_ + b_ + c_);
165}
166
167
169{
170 return triPointRef::centre(a(), b(), c());
171}
172
173
174template<class Point, class PointRef>
176(
177 const Point& p0,
178 const Point& p1,
179 const Point& p2
181{
182 return 0.5*((p1 - p0)^(p2 - p0));
183}
184
185
186template<class Point, class PointRef>
188{
189 return 0.5*((b_ - a_)^(c_ - a_));
190}
191
192
194{
195 return triPointRef::areaNormal(a(), b(), c());
196}
197
198
199template<class Point, class PointRef>
201(
202 const Point& p0,
203 const Point& p1,
204 const Point& p2
205)
206{
208 (void) n.normalise(ROOTVSMALL);
209 return n;
210}
211
212
213template<class Point, class PointRef>
215{
216 return triangle<Point, PointRef>::unitNormal(a_, b_, c_);
217}
218
219
221{
222 return triPointRef::unitNormal(a(), b(), c());
223}
224
225
226template<class Point, class PointRef>
228(
229 const Point& p0,
230 const Point& p1,
231 const Point& p2
232)
233{
234 return Pair<Point>
235 (
236 min(p0, min(p1, p2)),
237 max(p0, max(p1, p2))
238 );
239}
240
241
242template<class Point, class PointRef>
244{
245 return triangle<Point, PointRef>::box(a_, b_, c_);
246}
247
248
250{
251 return triPointRef::box(a(), b(), c());
252}
253
254
255template<class Point, class PointRef>
257{
258 return Point(c_ - b_);
259}
260
261template<class Point, class PointRef>
263{
264 return Point(a_ - c_);
265}
266
267template<class Point, class PointRef>
269{
270 return Point(b_ - a_);
271}
272
275{
276 return vector(c() - b());
277}
278
281{
282 return vector(a() - c());
283}
284
285
288 return vector(b() - a());
289}
290
291
292// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
295{
296 return triPointRef(a(), b(), c());
298
300inline void Foam::triPoints::flip()
301{
302 // swap pt1 <-> pt2
303 point t(b());
304 b() = c();
305 c() = t;
306}
307
308
309inline Foam::scalar Foam::triPoints::mag() const
310{
311 return ::Foam::mag(areaNormal());
312}
313
314
315template<class Point, class PointRef>
316inline Foam::scalar Foam::triangle<Point, PointRef>::mag() const
317{
318 return ::Foam::mag(areaNormal());
319}
320
321
322inline Foam::scalar Foam::triPoints::magSqr() const
323{
324 return ::Foam::magSqr(areaNormal());
325}
326
327
328template<class Point, class PointRef>
330{
331 return ::Foam::magSqr(areaNormal());
332}
333
334
335template<class Point, class PointRef>
337{
338 scalar d1 = (c_ - a_) & (b_ - a_);
339 scalar d2 = -(c_ - b_) & (b_ - a_);
340 scalar d3 = (c_ - a_) & (c_ - b_);
341
342 scalar c1 = d2*d3;
343 scalar c2 = d3*d1;
344 scalar c3 = d1*d2;
345
346 scalar c = c1 + c2 + c3;
347
348 if (Foam::mag(c) < ROOTVSMALL)
349 {
350 // Degenerate triangle, returning centre instead of circumCentre.
351
352 return centre();
353 }
354
355 return
357 ((c2 + c3)*a_ + (c3 + c1)*b_ + (c1 + c2)*c_)/(2*c)
358 );
359}
360
361
362template<class Point, class PointRef>
363inline Foam::scalar Foam::triangle<Point, PointRef>::circumRadius() const
364{
365 const scalar d1 = (c_ - a_) & (b_ - a_);
366 const scalar d2 = -(c_ - b_) & (b_ - a_);
367 const scalar d3 = (c_ - a_) & (c_ - b_);
368
369 const scalar denom = d2*d3 + d3*d1 + d1*d2;
370
371 if (Foam::mag(denom) < VSMALL)
372 {
373 // Degenerate triangle, returning GREAT for circumRadius.
374
375 return GREAT;
376 }
378 const scalar a = (d1 + d2)*(d2 + d3)*(d3 + d1) / denom;
379 return 0.5*Foam::sqrt(clamp(a, scalar(0), scalar(GREAT)));
380}
381
382
383template<class Point, class PointRef>
384inline Foam::scalar Foam::triangle<Point, PointRef>::quality() const
385{
386 scalar c = circumRadius();
387
388 if (c < ROOTVSMALL)
389 {
390 // zero circumRadius, something has gone wrong.
391 return SMALL;
393
394 return mag()/(Foam::sqr(c)*3.0*sqrt(3.0)/4.0);
395}
396
397
398template<class Point, class PointRef>
400(
401 const triangle& t
402) const
403{
404 return (1.0/12.0)*
405 (
406 ((t.a_ - a_) & ((b_ - a_)^(c_ - a_)))
407 + ((t.b_ - b_) & ((c_ - b_)^(t.a_ - b_)))
408 + ((c_ - t.c_) & ((t.b_ - t.c_)^(t.a_ - t.c_)))
409
410 + ((t.a_ - a_) & ((b_ - a_)^(c_ - a_)))
411 + ((b_ - t.b_) & ((t.a_ - t.b_)^(t.c_ - t.b_)))
412 + ((c_ - t.c_) & ((b_ - t.c_)^(t.a_ - t.c_)))
413 );
414}
415
416
417template<class Point, class PointRef>
419(
420 PointRef refPt,
421 scalar density
422) const
423{
424 Point aRel = a_ - refPt;
425 Point bRel = b_ - refPt;
426 Point cRel = c_ - refPt;
427
428 tensor V
429 (
430 aRel.x(), aRel.y(), aRel.z(),
431 bRel.x(), bRel.y(), bRel.z(),
432 cRel.x(), cRel.y(), cRel.z()
433 );
434
435 scalar a = Foam::mag((b_ - a_)^(c_ - a_));
436
437 tensor S = 1/24.0*(tensor::one + I);
438
439 return
440 (
441 a*I/24.0*
442 (
443 (aRel & aRel)
444 + (bRel & bRel)
445 + (cRel & cRel)
446 + ((aRel + bRel + cRel) & (aRel + bRel + cRel))
447 )
448 - a*(V.T() & S & V)
449 )
450 *density;
451}
452
453
454template<class Point, class PointRef>
456{
458}
459
460
461template<class Point, class PointRef>
463(
464 const barycentric2D& bary
465) const
466{
467 return bary[0]*a_ + bary[1]*b_ + bary[2]*c_;
468}
469
470
471template<class Point, class PointRef>
473(
474 const point& pt
475) const
476{
478 pointToBarycentric(pt, bary);
479 return bary;
480}
481
482
483template<class Point, class PointRef>
485(
486 const point& pt,
487 barycentric2D& bary
488) const
489{
490 // Reference:
491 // Real-time collision detection, Christer Ericson, 2005, p47-48
492
493 vector v0 = b_ - a_;
494 vector v1 = c_ - a_;
495 vector v2 = pt - a_;
496
497 scalar d00 = v0 & v0;
498 scalar d01 = v0 & v1;
499 scalar d11 = v1 & v1;
500 scalar d20 = v2 & v0;
501 scalar d21 = v2 & v1;
502
503 scalar denom = d00*d11 - d01*d01;
504
505 if (Foam::mag(denom) < SMALL)
506 {
507 // Degenerate triangle, returning 1/3 barycentric coordinates.
508
509 bary = barycentric2D(1.0/3.0, 1.0/3.0, 1.0/3.0);
510
511 return denom;
512 }
513
514 bary[1] = (d11*d20 - d01*d21)/denom;
515 bary[2] = (d00*d21 - d01*d20)/denom;
516 bary[0] = 1.0 - bary[1] - bary[2];
517
518 return denom;
519}
520
521
522template<class Point, class PointRef>
524(
525 const point& origin,
526 const vector& normal
527) const
528{
529 // Check points cut below(1) or above(2). Mix of above/below == 3
530 // Cf. plane::whichSide()
531 return
532 (
533 (
534 (((a_ - origin) & normal) < 0 ? 1 : 2)
535 | (((b_ - origin) & normal) < 0 ? 1 : 2)
536 | (((c_ - origin) & normal) < 0 ? 1 : 2)
537 ) == 3
538 );
539}
540
541
542template<class Point, class PointRef>
544(
545 const point& origin,
546 const vector::components axis
547) const
548{
549 // Direct check of points cut below(1) or above(2)
550 // Cf. plane::whichSide()
551 return
552 (
553 (
554 (a_[axis] < origin[axis] ? 1 : 2)
555 | (b_[axis] < origin[axis] ? 1 : 2)
556 | (c_[axis] < origin[axis] ? 1 : 2)
557 ) == 3
558 );
559}
560
561
562template<class Point, class PointRef>
564(
565 const point& p,
566 const vector& q,
567 const intersection::algorithm alg,
569) const
570{
571 // Express triangle in terms of baseVertex (point a_) and
572 // two edge vectors
573 vector E0 = b_ - a_;
574 vector E1 = c_ - a_;
575
576 // Initialize intersection to miss.
577 pointHit inter(p);
578
579 vector n(0.5*(E0 ^ E1));
580 scalar area = Foam::mag(n);
581
582 if (area < VSMALL)
583 {
584 // Ineligible miss.
585 inter.setMiss(false);
586
587 // The miss point is the nearest point on the triangle. Take any one.
588 inter.setPoint(a_);
589
590 // The distance to the miss is the distance between the
591 // original point and plane of intersection. No normal so use
592 // distance from p to a_
593 inter.setDistance(Foam::mag(a_ - p));
594
595 return inter;
596 }
597
598 vector q1 = q/Foam::mag(q);
599
601 {
602 n /= area;
603
604 return ray(p, q1 - n, alg, intersection::VECTOR);
605 }
606
607 // Intersection point with triangle plane
608 point pInter;
609 // Is intersection point inside triangle
610 bool hit;
611 {
612 // Reuse the fast ray intersection routine below in FULL_RAY
613 // mode since the original intersection routine has rounding problems.
615 hit = fastInter.hit();
616
617 if (hit)
618 {
619 pInter = fastInter.point();
620 }
621 else
622 {
623 // Calculate intersection of ray with triangle plane
624 vector v = a_ - p;
625 pInter = p + (q1&v)*q1;
626 }
627 }
628
629 // Distance to intersection point
630 scalar dist = q1 & (pInter - p);
631
632 const scalar planarPointTol =
634 (
636 (
637 Foam::mag(E0),
638 Foam::mag(E1)
639 ),
640 Foam::mag(c_ - b_)
642
643 bool eligible =
645 || (alg == intersection::HALF_RAY && dist > -planarPointTol)
646 || (
648 && ((q1 & areaNormal()) < -VSMALL)
649 );
650
651 if (hit && eligible)
652 {
653 // Hit. Set distance to intersection.
654 inter.hitPoint(pInter);
655 inter.setDistance(dist);
656 }
657 else
658 {
659 // Miss or ineligible hit.
660 inter.setMiss(eligible);
661
662 // The miss point is the nearest point on the triangle
663 inter.setPoint(nearestPoint(p).point());
664
665 // The distance to the miss is the distance between the
666 // original point and plane of intersection
667 inter.setDistance(Foam::mag(pInter - p));
668 }
669
670
671 return inter;
672}
673
674
675// From "Fast, Minimum Storage Ray/Triangle Intersection"
676// Moeller/Trumbore.
677template<class Point, class PointRef>
679(
680 const point& orig,
681 const vector& dir,
682 const intersection::algorithm alg,
683 const scalar tol
684) const
685{
686 const vector edge1 = b_ - a_;
687 const vector edge2 = c_ - a_;
688
689 // begin calculating determinant - also used to calculate U parameter
690 const vector pVec = dir ^ edge2;
691
692 // if determinant is near zero, ray lies in plane of triangle
693 const scalar det = edge1 & pVec;
694
695 // Initialise to miss
696 pointHit intersection(false, Zero, GREAT, false);
697
698 if (alg == intersection::VISIBLE)
699 {
700 // Culling branch
701 if (det < ROOTVSMALL)
702 {
703 // Ray on wrong side of triangle. Return miss
704 return intersection;
705 }
706 }
707 else if (alg == intersection::HALF_RAY || alg == intersection::FULL_RAY)
708 {
709 // Non-culling branch
710 if (det > -ROOTVSMALL && det < ROOTVSMALL)
711 {
712 // Ray parallel to triangle. Return miss
713 return intersection;
714 }
715 }
716
717 const scalar inv_det = 1.0 / det;
718
719 /* calculate distance from a_ to ray origin */
720 const vector tVec = orig-a_;
721
722 /* calculate U parameter and test bounds */
723 const scalar u = (tVec & pVec)*inv_det;
724
725 if (u < -tol || u > 1.0+tol)
726 {
727 // return miss
728 return intersection;
729 }
730
731 /* prepare to test V parameter */
732 const vector qVec = tVec ^ edge1;
733
734 /* calculate V parameter and test bounds */
735 const scalar v = (dir & qVec) * inv_det;
736
737 if (v < -tol || u + v > 1.0+tol)
738 {
739 // return miss
740 return intersection;
741 }
742
743 /* calculate t, scale parameters, ray intersects triangle */
744 const scalar t = (edge2 & qVec) * inv_det;
745
746 if (alg == intersection::HALF_RAY && t < -tol)
747 {
748 // Wrong side of orig. Return miss
749 return intersection;
750 }
751
752 intersection.hitPoint(a_ + u*edge1 + v*edge2);
753 intersection.setDistance(t);
754
755 return intersection;
756}
757
758
759template<class Point, class PointRef>
761(
762 const point& p,
763 label& nearType,
764 label& nearLabel
765) const
766{
767 // Adapted from:
768 // Real-time collision detection, Christer Ericson, 2005, p136-142
769
770 // Check if P in vertex region outside A
771 vector ab = b_ - a_;
772 vector ac = c_ - a_;
773 vector ap = p - a_;
774
775 scalar d1 = ab & ap;
776 scalar d2 = ac & ap;
777
778 if (d1 <= 0.0 && d2 <= 0.0)
779 {
780 // barycentric coordinates (1, 0, 0)
781
782 nearType = POINT;
783 nearLabel = 0;
784 return pointHit(false, a_, Foam::mag(a_ - p), true);
785 }
786
787 // Check if P in vertex region outside B
788 vector bp = p - b_;
789 scalar d3 = ab & bp;
790 scalar d4 = ac & bp;
791
792 if (d3 >= 0.0 && d4 <= d3)
793 {
794 // barycentric coordinates (0, 1, 0)
795
796 nearType = POINT;
797 nearLabel = 1;
798 return pointHit(false, b_, Foam::mag(b_ - p), true);
799 }
800
801 // Check if P in edge region of AB, if so return projection of P onto AB
802 scalar vc = d1*d4 - d3*d2;
803
804 if (vc <= 0.0 && d1 >= 0.0 && d3 <= 0.0)
805 {
806 if ((d1 - d3) < ROOTVSMALL)
807 {
808 // Degenerate triangle, for d1 = d3, a_ and b_ are likely coincident
809 nearType = POINT;
810 nearLabel = 0;
811 return pointHit(false, a_, Foam::mag(a_ - p), true);
812 }
813
814 // barycentric coordinates (1-v, v, 0)
815 scalar v = d1/(d1 - d3);
816
817 point nearPt = a_ + v*ab;
818 nearType = EDGE;
819 nearLabel = 0;
820 return pointHit(false, nearPt, Foam::mag(nearPt - p), true);
821 }
822
823 // Check if P in vertex region outside C
824 vector cp = p - c_;
825 scalar d5 = ab & cp;
826 scalar d6 = ac & cp;
827
828 if (d6 >= 0.0 && d5 <= d6)
829 {
830 // barycentric coordinates (0, 0, 1)
831
832 nearType = POINT;
833 nearLabel = 2;
834 return pointHit(false, c_, Foam::mag(c_ - p), true);
835 }
836
837 // Check if P in edge region of AC, if so return projection of P onto AC
838 scalar vb = d5*d2 - d1*d6;
839
840 if (vb <= 0.0 && d2 >= 0.0 && d6 <= 0.0)
841 {
842 if ((d2 - d6) < ROOTVSMALL)
843 {
844 // Degenerate triangle, for d2 = d6, a_ and c_ are likely coincident
845 nearType = POINT;
846 nearLabel = 0;
847 return pointHit(false, a_, Foam::mag(a_ - p), true);
848 }
849
850 // barycentric coordinates (1-w, 0, w)
851 scalar w = d2/(d2 - d6);
852
853 point nearPt = a_ + w*ac;
854 nearType = EDGE;
855 nearLabel = 2;
856 return pointHit(false, nearPt, Foam::mag(nearPt - p), true);
857 }
858
859 // Check if P in edge region of BC, if so return projection of P onto BC
860 scalar va = d3*d6 - d5*d4;
861
862 if (va <= 0.0 && (d4 - d3) >= 0.0 && (d5 - d6) >= 0.0)
863 {
864 if (((d4 - d3) + (d5 - d6)) < ROOTVSMALL)
865 {
866 // Degenerate triangle, for (d4 - d3) = (d6 - d5), b_ and c_ are
867 // likely coincident
868 nearType = POINT;
869 nearLabel = 1;
870 return pointHit(false, b_, Foam::mag(b_ - p), true);
871 }
872
873 // barycentric coordinates (0, 1-w, w)
874 scalar w = (d4 - d3)/((d4 - d3) + (d5 - d6));
875
876 point nearPt = b_ + w*(c_ - b_);
877 nearType = EDGE;
878 nearLabel = 1;
879 return pointHit(false, nearPt, Foam::mag(nearPt - p), true);
880 }
881
882 // P inside face region. Compute Q through its barycentric
883 // coordinates (u, v, w)
884
885 if ((va + vb + vc) < ROOTVSMALL)
886 {
887 // Degenerate triangle, return the centre because no edge or points are
888 // closest
889 point nearPt = centre();
890 nearType = NONE,
891 nearLabel = -1;
892 return pointHit(true, nearPt, Foam::mag(nearPt - p), false);
893 }
894
895 scalar denom = 1.0/(va + vb + vc);
896 scalar v = vb * denom;
897 scalar w = vc * denom;
898
899 // = u*a + v*b + w*c, u = va*denom = 1.0 - v - w
900
901 point nearPt = a_ + ab*v + ac*w;
902 nearType = NONE,
903 nearLabel = -1;
904 return pointHit(true, nearPt, Foam::mag(nearPt - p), false);
905}
906
907
908template<class Point, class PointRef>
910(
911 const point& p
912) const
913{
914 // Dummy labels
915 label nearType = -1;
916 label nearLabel = -1;
917
918 return nearestPointClassify(p, nearType, nearLabel);
919}
920
921
922template<class Point, class PointRef>
924(
925 const point& p,
926 label& nearType,
927 label& nearLabel
928) const
929{
930 return nearestPointClassify(p, nearType, nearLabel).hit();
931}
932
933
934template<class Point, class PointRef>
936(
937 const linePointRef& ln,
938 pointHit& lnInfo
939) const
940{
941 vector q = ln.vec();
942 pointHit triInfo
943 (
945 (
946 ln.start(),
947 q,
949 )
950 );
951
952 if (triInfo.hit())
953 {
954 // Line hits triangle. Find point on line.
955 if (triInfo.distance() > 1)
956 {
957 // Hit beyond endpoint
958 lnInfo.setMiss(true);
959 lnInfo.setPoint(ln.end());
960 scalar dist = triInfo.point().dist(lnInfo.point());
961 lnInfo.setDistance(dist);
962 triInfo.setMiss(true);
963 triInfo.setDistance(dist);
964 }
965 else if (triInfo.distance() < 0)
966 {
967 // Hit beyond startpoint
968 lnInfo.setMiss(true);
969 lnInfo.setPoint(ln.start());
970 scalar dist = triInfo.point().dist(lnInfo.point());
971 lnInfo.setDistance(dist);
972 triInfo.setMiss(true);
973 triInfo.setDistance(dist);
974 }
975 else
976 {
977 // Hit on line
978 lnInfo.hitPoint(triInfo.point());
979 lnInfo.setDistance(0);
980 triInfo.setDistance(0);
981 }
982 }
983 else
984 {
985 // Line skips triangle. See which triangle edge it gets closest to
986
987 point nearestEdgePoint;
988 point nearestLinePoint;
989 //label minEdgeIndex = 0;
990 scalar minDist = ln.nearestDist
991 (
992 linePointRef(a_, b_),
993 nearestLinePoint,
994 nearestEdgePoint
995 );
996
997 {
998 point linePoint;
999 point triEdgePoint;
1000 scalar dist = ln.nearestDist
1001 (
1002 linePointRef(b_, c_),
1003 linePoint,
1004 triEdgePoint
1005 );
1006 if (dist < minDist)
1007 {
1008 minDist = dist;
1009 nearestEdgePoint = triEdgePoint;
1010 nearestLinePoint = linePoint;
1011 //minEdgeIndex = 1;
1012 }
1013 }
1014
1015 {
1016 point linePoint;
1017 point triEdgePoint;
1018 scalar dist = ln.nearestDist
1019 (
1020 linePointRef(c_, a_),
1021 linePoint,
1022 triEdgePoint
1023 );
1024 if (dist < minDist)
1025 {
1026 minDist = dist;
1027 nearestEdgePoint = triEdgePoint;
1028 nearestLinePoint = linePoint;
1029 //minEdgeIndex = 2;
1030 }
1031 }
1032
1033 lnInfo.setDistance(minDist);
1034 triInfo.setDistance(minDist);
1035 triInfo.setMiss(false);
1036 triInfo.setPoint(nearestEdgePoint);
1037
1038 // Convert point on line to pointHit
1039 if (Foam::mag(nearestLinePoint-ln.start()) < SMALL)
1040 {
1041 lnInfo.setMiss(true);
1042 lnInfo.setPoint(ln.start());
1043 }
1044 else if (Foam::mag(nearestLinePoint-ln.end()) < SMALL)
1045 {
1046 lnInfo.setMiss(true);
1047 lnInfo.setPoint(ln.end());
1048 }
1049 else
1050 {
1051 lnInfo.hitPoint(nearestLinePoint);
1053 }
1054 return triInfo;
1055}
1056
1057
1058template<class Point, class PointRef>
1060(
1061 const point& p,
1062 const scalar tol
1063) const
1064{
1065 const scalar dist = ((p - a_) & unitNormal());
1066
1067 return ((dist < -tol) ? -1 : (dist > tol) ? +1 : 0);
1068}
1069
1070
1071template<class Point, class PointRef>
1073(
1074 const triPoints& tri
1076{
1077 area_ += triangle<Point, const Point&>(tri).mag();
1078}
1079
1080
1081template<class Point, class PointRef>
1083(
1084 triIntersectionList& tris,
1085 label& nTris
1086)
1088 tris_(tris),
1089 nTris_(nTris)
1090{}
1091
1092
1093template<class Point, class PointRef>
1095(
1096 const triPoints& tri
1097)
1098{
1099 tris_[nTris_++] = tri;
1100}
1101
1102
1103template<class Point, class PointRef>
1104inline Foam::point Foam::triangle<Point, PointRef>::planeIntersection
1105(
1106 const FixedList<scalar, 3>& d,
1107 const triPoints& t,
1108 const label negI,
1109 const label posI
1110)
1111{
1112 return (d[posI]*t[negI] - d[negI]*t[posI])/(-d[negI] + d[posI]);
1113}
1114
1115
1116// * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
1117
1118template<class Point, class PointRef>
1119inline Foam::Istream& Foam::operator>>
1120(
1121 Istream& is,
1122 triangle<Point, PointRef>& t
1123)
1124{
1125 is.readBegin("triangle");
1126 is >> t.a_ >> t.b_ >> t.c_;
1127 is.readEnd("triangle");
1128
1129 is.check(FUNCTION_NAME);
1130 return is;
1131}
1132
1133
1134template<class Point, class PointRef>
1135inline Foam::Ostream& Foam::operator<<
1136(
1137 Ostream& os,
1139)
1140{
1141 // Format like FixedList
1142 os << nl
1144 << t.a_ << token::SPACE
1145 << t.b_ << token::SPACE
1146 << t.c_
1147 << token::END_LIST;
1148
1149 return os;
1150}
1151
1152
1153// ************************************************************************* //
CGAL::Point_3< K > Point
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
label n
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
An ordered pair of two objects of type <T> with first() and second() elements.
Definition Pair.H:66
scalar distance() const noexcept
Return distance to hit.
Definition pointHit.H:169
void setPoint(const point_type &p)
Set the point.
Definition pointHit.H:235
void setDistance(const scalar d) noexcept
Set the distance.
Definition pointHit.H:243
bool hit() const noexcept
Is there a hit.
Definition pointHit.H:145
void setMiss(const bool eligible) noexcept
Set the hit status off and set the eligible miss status.
Definition pointHit.H:226
const point_type & point() const noexcept
Return the point, no checks.
Definition pointHit.H:161
const point_type & hitPoint() const
Return the hit point. Fatal if not hit.
Definition pointHit.H:177
Random number generator.
Definition Random.H:56
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
scalar dist(const Vector< Cmpt > &v2) const
The L2-norm distance from another vector. The mag() of the difference.
Definition VectorI.H:107
Foam::intersection.
static scalar planarTol()
Return planar tolerance.
@ BEGIN_LIST
Begin list [isseparator].
Definition token.H:174
@ END_LIST
End list [isseparator].
Definition token.H:175
@ SPACE
Space [isspace].
Definition token.H:144
Triangle point storage. Default constructable (triangle is not).
Definition triangle.H:77
void flip()
Flip triangle orientation by swapping second and third vertices.
Definition triangleI.H:293
triPointRef tri() const
Return as triangle reference.
Definition triangleI.H:287
const point & c() const noexcept
The third vertex.
Definition triangle.H:151
vector vecA() const
Edge vector opposite point a(): from b() to c().
Definition triangleI.H:267
vector areaNormal() const
The area normal - with magnitude equal to area of triangle.
Definition triangleI.H:186
const point & b() const noexcept
The second vertex.
Definition triangle.H:146
Pair< point > box() const
The enclosing (bounding) box for the triangle.
Definition triangleI.H:242
vector unitNormal() const
Return unit normal.
Definition triangleI.H:213
scalar mag() const
The magnitude of the triangle area.
Definition triangleI.H:302
vector vecC() const
Edge vector opposite point c(): from a() to b().
Definition triangleI.H:279
triPoints()=default
Default construct.
const point & a() const noexcept
The first vertex.
Definition triangle.H:141
scalar magSqr() const
The magnitude squared of the triangle area.
Definition triangleI.H:315
point centre() const
Return centre (centroid).
Definition triangleI.H:161
vector vecB() const
Edge vector opposite point b(): from c() to a().
Definition triangleI.H:273
A triangle primitive used to calculate face normals and swept volumes. Uses referred points.
Definition triangle.H:234
Point circumCentre() const
Return circum-centre.
Definition triangleI.H:329
barycentric2D pointToBarycentric(const point &pt) const
Calculate the barycentric coordinates from the given point.
Definition triangleI.H:466
static vector areaNormal(const Point &p0, const Point &p1, const Point &p2)
The area normal for a triangle defined by three points (right-hand rule). Magnitude equal to area of ...
Definition triangleI.H:169
scalar circumRadius() const
Return circum-radius.
Definition triangleI.H:356
Point barycentricToPoint(const barycentric2D &bary) const
Calculate the point from the given barycentric coordinates.
Definition triangleI.H:456
scalar sweptVol(const triangle &t) const
Return swept-volume.
Definition triangleI.H:393
vector areaNormal() const
The area normal - with magnitude equal to area of triangle.
Definition triangleI.H:180
Point centre() const
Return centre (centroid).
Definition triangleI.H:155
Point vecB() const
Edge vector opposite point b(): from c() to a().
Definition triangleI.H:255
FixedList< triPoints, 27 > triIntersectionList
Storage type for triangles originating from intersecting triangle with another triangle.
Definition triangle.H:248
pointHit nearestPointClassify(const point &p, label &nearType, label &nearLabel) const
Find the nearest point to p on the triangle and classify it:
Definition triangleI.H:754
@ POINT
Close to point.
Definition triangle.H:256
@ EDGE
Close to edge.
Definition triangle.H:257
@ NONE
Unknown proximity.
Definition triangle.H:255
triangle(const Point &p0, const Point &p1, const Point &p2)
Construct from three points.
Definition triangleI.H:81
static Point centre(const Point &p0, const Point &p1, const Point &p2)
The centre (centroid) of three points.
Definition triangleI.H:144
pointHit ray(const point &p, const vector &q, const intersection::algorithm=intersection::FULL_RAY, const intersection::direction dir=intersection::VECTOR) const
Return point intersection with a ray.
Definition triangleI.H:557
vector unitNormal() const
Return unit normal.
Definition triangleI.H:207
tensor inertia(PointRef refPt=Zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
Definition triangleI.H:412
scalar mag() const
The magnitude of the triangle area.
Definition triangleI.H:309
bool intersects(const point &origin, const vector &normal) const
Fast intersection detection with a plane.
Definition triangleI.H:517
pointHit nearestPoint(const point &p) const
Return nearest point to p on triangle.
Definition triangleI.H:903
Point randomPoint(Random &rndGen) const
Return a random point on the triangle from a uniform distribution.
Definition triangleI.H:448
Pair< Point > box() const
The enclosing (bounding) box for the triangle.
Definition triangleI.H:236
pointHit intersection(const point &p, const vector &q, const intersection::algorithm alg, const scalar tol=0.0) const
Fast intersection with a ray.
Definition triangleI.H:672
bool classify(const point &p, label &nearType, label &nearLabel) const
Classify nearest point to p in triangle plane.
Definition triangleI.H:917
scalar magSqr() const
The magnitude squared of the triangle area.
Definition triangleI.H:322
int sign(const point &p, const scalar tol=SMALL) const
The sign for which side of the face plane the point is on.
Definition triangleI.H:1053
static vector unitNormal(const Point &p0, const Point &p1, const Point &p2)
The unit normal for a triangle defined by three points (right-hand rule).
Definition triangleI.H:194
scalar quality() const
Return quality: Ratio of triangle and circum-circle.
Definition triangleI.H:377
const point & c() const noexcept
Definition triangle.H:408
Point vecA() const
Edge vector opposite point a(): from b() to c().
Definition triangleI.H:249
const Point & a() const noexcept
The first vertex.
Definition triangle.H:398
Point vecC() const
Edge vector opposite point c(): from a() to b().
Definition triangleI.H:261
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
volScalarField & p
const volScalarField & p0
Definition EEqn.H:36
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
#define FUNCTION_NAME
@ NONE
No type, or default initialized type.
Definition exprTraits.H:83
dimensionedScalar det(const dimensionedSphericalTensor &dt)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
barycentric2D barycentric2D01(Random &rndGen)
Generate a random barycentric coordinate within the unit triangle.
PointHit< point > pointHit
A PointHit with a 3D point.
Definition pointHit.H:51
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Barycentric2D< scalar > barycentric2D
A scalar version of the templated Barycentric2D.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
static const Identity< scalar > I
Definition Identity.H:100
Tensor< scalar > tensor
Definition symmTensor.H:57
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
line< point, const point & > linePointRef
A line using referred points.
Definition line.H:66
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition triangleFwd.H:46
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy the source to the destination (recursively if necessary).
Definition POSIX.C:1065
Vector< scalar > vector
Definition vector.H:57
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition POSIX.C:1239
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
const pointField & pts
volScalarField & b
storeOp(triIntersectionList &, label &)
Definition triangleI.H:1076
triIntersectionList & tris_
Definition triangle.H:294
Random rndGen