Loading...
Searching...
No Matches
tetrahedronI.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-2022 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 "triangle.H"
30#include "IOstreams.H"
31#include "plane.H"
32
33// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34
36(
37 const point& p0,
38 const point& p1,
39 const point& p2,
40 const point& p3
41)
42{
43 a() = p0;
44 b() = p1;
45 c() = p2;
46 d() = p3;
47}
48
49
51{
52 a() = pts.a();
53 b() = pts.b();
54 c() = pts.c();
55 d() = pts.d();
56}
57
66(
67 const UList<point>& points,
68 const FixedList<label, 4>& indices
70:
71 FixedList<point, 4>(points, indices)
72{}
73
74
75template<class Point, class PointRef>
77(
78 const Point& p0,
79 const Point& p1,
80 const Point& p2,
81 const Point& p3
82)
83:
84 a_(p0),
85 b_(p1),
86 c_(p2),
87 d_(p3)
88{}
89
90
91template<class Point, class PointRef>
93(
95)
96:
97 a_(pts.template get<0>()),
98 b_(pts.template get<1>()),
99 c_(pts.template get<2>()),
100 d_(pts.template get<3>())
101{}
102
103
104template<class Point, class PointRef>
106(
107 const UList<Point>& points,
108 const FixedList<label, 4>& indices
109)
110:
111 a_(points[indices.template get<0>()]),
112 b_(points[indices.template get<1>()]),
113 c_(points[indices.template get<2>()]),
114 d_(points[indices.template get<3>()])
115{}
116
117
118template<class Point, class PointRef>
121 is >> *this;
122}
123
124
125// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
128{
129 return tetPointRef(a(), b(), c(), d());
130}
131
132
133inline void Foam::tetPoints::flip()
134{
135 // swap pt2 <-> pt3
136 point t(c());
137 c() = d();
138 d() = t;
139}
140
141
142// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
143
144template<class Point, class PointRef>
146(
147 const Point& p0,
148 const Point& p1,
149 const Point& p2,
150 const Point& p3
151)
152{
153 return Pair<Point>
154 (
155 min(p0, min(p1, min(p2, p3))),
156 max(p0, max(p1, max(p2, p3)))
157 );
158}
159
160
161template<class Point, class PointRef>
163{
164 return tetrahedron<Point, PointRef>::box(a_, b_, c_, d_);
165}
166
167
169{
170 return tetPointRef::box(a(), b(), c(), d());
171}
172
173
174template<class Point, class PointRef>
176{
177 return treeBoundBox(box());
178}
179
180
182{
183 return treeBoundBox(box());
185
186
187// Warning. Ordering of faces needs to be the same for a tetrahedron class,
188// tetrahedron cell shape model and a tetCell
189
190template<class Point, class PointRef>
192(
193 const label facei
194) const
195{
196 // Warning. Ordering of faces needs to be the same for a tetrahedron
197 // class, a tetrahedron cell shape model and a tetCell
198
199 if (facei == 0)
200 {
201 return triPointRef(b_, c_, d_);
202 }
203 else if (facei == 1)
204 {
205 return triPointRef(a_, d_, c_);
206 }
207 else if (facei == 2)
208 {
209 return triPointRef(a_, b_, d_);
210 }
211 else if (facei == 3)
212 {
213 return triPointRef(a_, c_, b_);
214 }
215
217 << "Face index (" << facei << ") out of range 0..3\n"
218 << abort(FatalError);
219 return triPointRef(b_, c_, d_);
220}
221
222
223template<class Point, class PointRef>
228
229
230template<class Point, class PointRef>
235
236
237template<class Point, class PointRef>
242
243
244template<class Point, class PointRef>
249
251template<class Point, class PointRef>
253{
254 return 0.25*(a_ + b_ + c_ + d_);
255}
256
257
258template<class Point, class PointRef>
260{
261 return (1.0/6.0)*(((b_ - a_) ^ (c_ - a_)) & (d_ - a_));
262}
264
265template<class Point, class PointRef>
267{
268 vector a = b_ - a_;
269 vector b = c_ - a_;
270 vector c = d_ - a_;
271
272 scalar lambda = magSqr(c) - (a & c);
273 scalar mu = magSqr(b) - (a & b);
274
275 vector ba = b ^ a;
276 vector ca = c ^ a;
277
278 vector num = lambda*ba - mu*ca;
279 scalar denom = (c & ba);
280
281 if (Foam::mag(denom) < ROOTVSMALL)
282 {
283 // Degenerate tetrahedron, returning centre instead of circumCentre.
284
285 return centre();
287
288 return a_ + 0.5*(a + num/denom);
289}
290
291
292template<class Point, class PointRef>
293inline Foam::scalar Foam::tetrahedron<Point, PointRef>::circumRadius() const
294{
295 vector a = b_ - a_;
296 vector b = c_ - a_;
297 vector c = d_ - a_;
298
299 scalar lambda = magSqr(c) - (a & c);
300 scalar mu = magSqr(b) - (a & b);
301
302 vector ba = b ^ a;
303 vector ca = c ^ a;
304
305 vector num = lambda*ba - mu*ca;
306 scalar denom = (c & ba);
307
308 if (Foam::mag(denom) < ROOTVSMALL)
309 {
310 // Degenerate tetrahedron, returning GREAT for circumRadius.
311 return GREAT;
313
314 return Foam::mag(0.5*(a + num/denom));
315}
316
317
318template<class Point, class PointRef>
319inline Foam::scalar Foam::tetrahedron<Point, PointRef>::quality() const
320{
321 return
322 mag()
323 /(
324 8.0/(9.0*sqrt(3.0))
326 + ROOTVSMALL
327 );
328}
329
330
331template<class Point, class PointRef>
333(
335) const
336{
338}
339
340
341template<class Point, class PointRef>
343(
344 const barycentric& bary
345) const
346{
347 return Point(bary.a()*a_ + bary.b()*b_ + bary.c()*c_ + bary.d()*d_);
348}
349
350
351template<class Point, class PointRef>
353(
354 const point& pt
355) const
356{
358 pointToBarycentric(pt, bary);
359 return bary;
360}
361
362
363template<class Point, class PointRef>
365(
366 const point& pt,
367 barycentric& bary
368) const
369{
370 // Reference:
371 // http://en.wikipedia.org/wiki/Barycentric_coordinate_system_(mathematics)
372
373 vector e0(a_ - d_);
374 vector e1(b_ - d_);
375 vector e2(c_ - d_);
376
377 tensor t
378 (
379 e0.x(), e1.x(), e2.x(),
380 e0.y(), e1.y(), e2.y(),
381 e0.z(), e1.z(), e2.z()
382 );
383
384 scalar detT = det(t);
385
386 if (Foam::mag(detT) < SMALL)
387 {
388 // Degenerate tetrahedron, returning 1/4 barycentric coordinates
389
390 bary = barycentric(0.25, 0.25, 0.25, 0.25);
391
392 return detT;
393 }
394
395 vector res = inv(t, detT) & (pt - d_);
396
397 bary[0] = res.x();
398 bary[1] = res.y();
399 bary[2] = res.z();
400 bary[3] = 1 - bary[0] - bary[1] - bary[2];
401
402 return detT;
403}
404
405
406template<class Point, class PointRef>
408(
409 const point& p
410) const
411{
412 // Adapted from:
413 // Real-time collision detection, Christer Ericson, 2005, p142-144
414
415 // Assuming initially that the point is inside all of the
416 // halfspaces, so closest to itself.
417
418 point closestPt = p;
419
420 scalar minOutsideDistance = VGREAT;
421
422 bool inside = true;
423
424 // Side a
425 {
426 const triangle<Point, PointRef> tria(b_, c_, d_);
427
428 if (((p - b_) & tria.areaNormal()) >= 0)
429 {
430 // p is outside halfspace plane of tri
431 pointHit info = tria.nearestPoint(p);
432
433 inside = false;
434
435 if (info.distance() < minOutsideDistance)
436 {
437 closestPt = info.point();
438
439 minOutsideDistance = info.distance();
440 }
441 }
442 }
443
444 // Side b
445 {
446 const triangle<Point, PointRef> tria(a_, d_, c_);
447
448 if (((p - a_) & tria.areaNormal()) >= 0)
449 {
450 // p is outside halfspace plane of tri
451 pointHit info = tria.nearestPoint(p);
452
453 inside = false;
454
455 if (info.distance() < minOutsideDistance)
456 {
457 closestPt = info.point();
458
459 minOutsideDistance = info.distance();
460 }
461 }
462 }
463
464 // Side c
465 {
466 const triangle<Point, PointRef> tria(a_, b_, d_);
467
468 if (((p - a_) & tria.areaNormal()) >= 0)
469 {
470 // p is outside halfspace plane of tri
471 pointHit info = tria.nearestPoint(p);
472
473 inside = false;
474
475 if (info.distance() < minOutsideDistance)
476 {
477 closestPt = info.point();
478
479 minOutsideDistance = info.distance();
480 }
481 }
482 }
483
484 // Side c
485 {
486 const triangle<Point, PointRef> tria(a_, c_, b_);
487
488 if (((p - a_) & tria.areaNormal()) >= 0)
489 {
490 // p is outside halfspace plane of tri
491 pointHit info = tria.nearestPoint(p);
492
493 inside = false;
494
495 if (info.distance() < minOutsideDistance)
496 {
497 closestPt = info.point();
498
499 minOutsideDistance = info.distance();
500 }
501 }
502 }
503
504 // If the point is inside, then the distance to the closest point
505 // is zero
506 if (inside)
507 {
508 minOutsideDistance = 0;
509 }
510
511 return pointHit
512 (
513 inside,
514 closestPt,
515 minOutsideDistance,
516 !inside
517 );
518}
519
520
521template<class Point, class PointRef>
523{
524 // For robustness, assuming that the point is in the tet unless
525 // "definitively" shown otherwise by obtaining a positive dot
526 // product greater than a tolerance of SMALL.
527
528 // The tet is defined: tet(Cc, tetBasePt, pA, pB) where the area normal
529 // vectors and base points for the half-space planes are:
530 // area[0] = Sa();
531 // area[1] = Sb();
532 // area[2] = Sc();
533 // area[3] = Sd();
534 // planeBase[0] = tetBasePt = b_
535 // planeBase[1] = ptA = c_
536 // planeBase[2] = tetBasePt = b_
537 // planeBase[3] = tetBasePt = b_
538
539 // Side a
540 {
541 const triangle<Point, PointRef> tria(b_, c_, d_);
542
543 if (((p - b_) & tria.unitNormal()) > SMALL)
544 {
545 return false;
546 }
547 }
548
549 // Side b
550 {
551 const triangle<Point, PointRef> tria(a_, d_, c_);
552
553 if (((p - a_) & tria.unitNormal()) > SMALL)
554 {
555 return false;
556 }
557 }
558
559 // Side c
560 {
561 const triangle<Point, PointRef> tria(a_, b_, d_);
562
563 if (((p - a_) & tria.unitNormal()) > SMALL)
564 {
565 return false;
566 }
567 }
568
569 // Side d
570 {
571 const triangle<Point, PointRef> tria(a_, c_, b_);
572
573 if (((p - a_) & tria.unitNormal()) > SMALL)
574 {
575 return false;
576 }
578
579 return true;
580}
581
582
583template<class Point, class PointRef>
585(
586 const tetPoints& tet
588{
589 vol_ += tet.tet().mag();
590}
591
592
593template<class Point, class PointRef>
595(
597 label& nTets
598)
600 tets_(tets),
601 nTets_(nTets)
602{}
603
604
605template<class Point, class PointRef>
607(
608 const tetPoints& tet
609)
610{
611 tets_[nTets_++] = tet;
612}
613
614
615template<class Point, class PointRef>
616inline Foam::point Foam::tetrahedron<Point, PointRef>::planeIntersection
617(
618 const FixedList<scalar, 4>& d,
619 const tetPoints& t,
620 const label negI,
621 const label posI
622)
623{
624 return
625 (d[posI]*t[negI] - d[negI]*t[posI])
626 / (-d[negI]+d[posI]);
627}
628
629
630template<class Point, class PointRef>
631template<class TetOp>
632inline void Foam::tetrahedron<Point, PointRef>::decomposePrism
633(
634 const FixedList<point, 6>& points,
635 TetOp& op
636)
637{
638 op(tetPoints(points[1], points[3], points[2], points[0]));
639 op(tetPoints(points[1], points[2], points[3], points[4]));
640 op(tetPoints(points[4], points[2], points[3], points[5]));
641}
642
643
644template<class Point, class PointRef>
645template<class AboveTetOp, class BelowTetOp>
646inline void Foam::tetrahedron<Point, PointRef>::tetSliceWithPlane
647(
648 const plane& pln,
649 const tetPoints& tet,
650 AboveTetOp& aboveOp,
651 BelowTetOp& belowOp
652)
653{
654 // Distance to plane
656 label nPos = 0;
657 forAll(tet, i)
658 {
659 d[i] = pln.signedDistance(tet[i]);
660 if (d[i] > 0)
661 {
662 ++nPos;
663 }
664 }
665
666 if (nPos == 4)
667 {
668 aboveOp(tet);
669 }
670 else if (nPos == 3)
671 {
672 // Sliced into below tet and above prism.
673 // Prism gets split into two tets.
674
675 // Find the below tet
676 label i0 = -1;
677 forAll(d, i)
678 {
679 if (d[i] <= 0)
680 {
681 i0 = i;
682 break;
683 }
684 }
685
686 label i1 = d.fcIndex(i0);
687 label i2 = d.fcIndex(i1);
688 label i3 = d.fcIndex(i2);
689
690 point p01(planeIntersection(d, tet, i0, i1));
691 point p02(planeIntersection(d, tet, i0, i2));
692 point p03(planeIntersection(d, tet, i0, i3));
693
694 // i0 = tetCell vertex 0: p01,p02,p03 outwards pointing triad
695 // ,, 1 : ,, inwards pointing triad
696 // ,, 2 : ,, outwards pointing triad
697 // ,, 3 : ,, inwards pointing triad
698
699 //Pout<< "Split 3pos tet " << tet << " d:" << d << " into" << nl;
700
701 if (i0 == 0 || i0 == 2)
702 {
703 tetPoints t(tet[i0], p01, p02, p03);
704 //Pout<< " belowtet:" << t << " around i0:" << i0 << endl;
705 //checkTet(t, "nPos 3, belowTet i0==0 or 2");
706 belowOp(t);
707
708 // Prism
710 (
711 {
712 tet[i1],
713 tet[i3],
714 tet[i2],
715 p01,
716 p03,
717 p02
718 }
719 );
720
721 //Pout<< " aboveprism:" << p << endl;
722 decomposePrism(p, aboveOp);
723 }
724 else
725 {
726 tetPoints t(p01, p02, p03, tet[i0]);
727 //Pout<< " belowtet:" << t << " around i0:" << i0 << endl;
728 //checkTet(t, "nPos 3, belowTet i0==1 or 3");
729 belowOp(t);
730
731 // Prism
733 (
734 {
735 tet[i3],
736 tet[i1],
737 tet[i2],
738 p03,
739 p01,
740 p02
741 }
742 );
743 //Pout<< " aboveprism:" << p << endl;
744 decomposePrism(p, aboveOp);
745 }
746 }
747 else if (nPos == 2)
748 {
749 // Tet cut into two prisms. Determine the positive one.
750 label pos0 = -1;
751 label pos1 = -1;
752 forAll(d, i)
753 {
754 if (d[i] > 0)
755 {
756 if (pos0 == -1)
757 {
758 pos0 = i;
759 }
760 else
761 {
762 pos1 = i;
763 }
764 }
765 }
766
767 //Pout<< "Split 2pos tet " << tet << " d:" << d
768 // << " around pos0:" << pos0 << " pos1:" << pos1
769 // << " neg0:" << neg0 << " neg1:" << neg1 << " into" << nl;
770
771 const edge posEdge(pos0, pos1);
772
773 if (posEdge == edge(0, 1))
774 {
775 point p02(planeIntersection(d, tet, 0, 2));
776 point p03(planeIntersection(d, tet, 0, 3));
777 point p12(planeIntersection(d, tet, 1, 2));
778 point p13(planeIntersection(d, tet, 1, 3));
779 // Split the resulting prism
780 {
782 (
783 {
784 tet[0],
785 p02,
786 p03,
787 tet[1],
788 p12,
789 p13
790 }
791 );
792
793 //Pout<< " 01 aboveprism:" << p << endl;
794 decomposePrism(p, aboveOp);
795 }
796 {
798 (
799 {
800 tet[2],
801 p02,
802 p12,
803 tet[3],
804 p03,
805 p13
806 }
807 );
808
809 //Pout<< " 01 belowprism:" << p << endl;
810 decomposePrism(p, belowOp);
811 }
812 }
813 else if (posEdge == edge(1, 2))
814 {
815 point p01(planeIntersection(d, tet, 0, 1));
816 point p13(planeIntersection(d, tet, 1, 3));
817 point p02(planeIntersection(d, tet, 0, 2));
818 point p23(planeIntersection(d, tet, 2, 3));
819 // Split the resulting prism
820 {
822 (
823 {
824 tet[1],
825 p01,
826 p13,
827 tet[2],
828 p02,
829 p23
830 }
831 );
832
833 //Pout<< " 12 aboveprism:" << p << endl;
834 decomposePrism(p, aboveOp);
835 }
836 {
838 (
839 {
840 tet[3],
841 p23,
842 p13,
843 tet[0],
844 p02,
845 p01
846 }
847 );
848
849 //Pout<< " 12 belowprism:" << p << endl;
850 decomposePrism(p, belowOp);
851 }
852 }
853 else if (posEdge == edge(2, 0))
854 {
855 point p01(planeIntersection(d, tet, 0, 1));
856 point p03(planeIntersection(d, tet, 0, 3));
857 point p12(planeIntersection(d, tet, 1, 2));
858 point p23(planeIntersection(d, tet, 2, 3));
859 // Split the resulting prism
860 {
862 (
863 {
864 tet[2],
865 p12,
866 p23,
867 tet[0],
868 p01,
869 p03
870 }
871 );
872
873 //Pout<< " 20 aboveprism:" << p << endl;
874 decomposePrism(p, aboveOp);
875 }
876 {
878 (
879 {
880 tet[1],
881 p12,
882 p01,
883 tet[3],
884 p23,
885 p03
886 }
887 );
888
889 //Pout<< " 20 belowprism:" << p << endl;
890 decomposePrism(p, belowOp);
891 }
892 }
893 else if (posEdge == edge(0, 3))
894 {
895 point p01(planeIntersection(d, tet, 0, 1));
896 point p02(planeIntersection(d, tet, 0, 2));
897 point p13(planeIntersection(d, tet, 1, 3));
898 point p23(planeIntersection(d, tet, 2, 3));
899 // Split the resulting prism
900 {
902 (
903 {
904 tet[3],
905 p23,
906 p13,
907 tet[0],
908 p02,
909 p01
910 }
911 );
912
913 //Pout<< " 03 aboveprism:" << p << endl;
914 decomposePrism(p, aboveOp);
915 }
916 {
918 (
919 {
920 tet[2],
921 p23,
922 p02,
923 tet[1],
924 p13,
925 p01
926 }
927 );
928
929 //Pout<< " 03 belowprism:" << p << endl;
930 decomposePrism(p, belowOp);
931 }
932 }
933 else if (posEdge == edge(1, 3))
934 {
935 point p01(planeIntersection(d, tet, 0, 1));
936 point p12(planeIntersection(d, tet, 1, 2));
937 point p03(planeIntersection(d, tet, 0, 3));
938 point p23(planeIntersection(d, tet, 2, 3));
939 // Split the resulting prism
940 {
942 (
943 {
944 tet[1],
945 p12,
946 p01,
947 tet[3],
948 p23,
949 p03
950 }
951 );
952
953 //Pout<< " 13 aboveprism:" << p << endl;
954 decomposePrism(p, aboveOp);
955 }
956 {
958 (
959 {
960 tet[2],
961 p12,
962 p23,
963 tet[0],
964 p01,
965 p03
966 }
967 );
968
969 //Pout<< " 13 belowprism:" << p << endl;
970 decomposePrism(p, belowOp);
971 }
972 }
973 else if (posEdge == edge(2, 3))
974 {
975 point p02(planeIntersection(d, tet, 0, 2));
976 point p12(planeIntersection(d, tet, 1, 2));
977 point p03(planeIntersection(d, tet, 0, 3));
978 point p13(planeIntersection(d, tet, 1, 3));
979 // Split the resulting prism
980 {
982 (
983 {
984 tet[2],
985 p02,
986 p12,
987 tet[3],
988 p03,
989 p13
990 }
991 );
992
993 //Pout<< " 23 aboveprism:" << p << endl;
994 decomposePrism(p, aboveOp);
995 }
996 {
998 (
999 {
1000 tet[0],
1001 p02,
1002 p03,
1003 tet[1],
1004 p12,
1005 p13
1006 }
1007 );
1008
1009 //Pout<< " 23 belowprism:" << p << endl;
1010 decomposePrism(p, belowOp);
1011 }
1012 }
1013 else
1014 {
1016 << "Missed edge:" << posEdge
1017 << abort(FatalError);
1018 }
1019 }
1020 else if (nPos == 1)
1021 {
1022 // Find the positive tet
1023 label i0 = -1;
1024 forAll(d, i)
1025 {
1026 if (d[i] > 0)
1027 {
1028 i0 = i;
1029 break;
1030 }
1031 }
1032
1033 label i1 = d.fcIndex(i0);
1034 label i2 = d.fcIndex(i1);
1035 label i3 = d.fcIndex(i2);
1036
1037 point p01(planeIntersection(d, tet, i0, i1));
1038 point p02(planeIntersection(d, tet, i0, i2));
1039 point p03(planeIntersection(d, tet, i0, i3));
1040
1041 //Pout<< "Split 1pos tet " << tet << " d:" << d << " into" << nl;
1042
1043 if (i0 == 0 || i0 == 2)
1044 {
1045 tetPoints t(tet[i0], p01, p02, p03);
1046 //Pout<< " abovetet:" << t << " around i0:" << i0 << endl;
1047 //checkTet(t, "nPos 1, aboveTets i0==0 or 2");
1048 aboveOp(t);
1049
1050 // Prism
1052 (
1053 {
1054 tet[i1],
1055 tet[i3],
1056 tet[i2],
1057 p01,
1058 p03,
1059 p02
1060 }
1061 );
1062
1063 //Pout<< " belowprism:" << p << endl;
1064 decomposePrism(p, belowOp);
1065 }
1066 else
1067 {
1068 tetPoints t(p01, p02, p03, tet[i0]);
1069 //Pout<< " abovetet:" << t << " around i0:" << i0 << endl;
1070 //checkTet(t, "nPos 1, aboveTets i0==1 or 3");
1071 aboveOp(t);
1072
1073 // Prism
1075 (
1076 {
1077 tet[i3],
1078 tet[i1],
1079 tet[i2],
1080 p03,
1081 p01,
1082 p02
1083 }
1084 );
1085
1086 //Pout<< " belowprism:" << p << endl;
1087 decomposePrism(p, belowOp);
1088 }
1089 }
1090 else // nPos == 0
1091 {
1092 belowOp(tet);
1093 }
1094}
1095
1096
1097template<class Point, class PointRef>
1098template<class AboveTetOp, class BelowTetOp>
1100(
1101 const plane& pl,
1102 AboveTetOp& aboveOp,
1103 BelowTetOp& belowOp
1104) const
1105{
1106 tetSliceWithPlane(pl, tetPoints(a_, b_, c_, d_), aboveOp, belowOp);
1107}
1108
1109
1110// * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
1111
1112template<class Point, class PointRef>
1113inline Foam::Istream& Foam::operator>>
1114(
1115 Istream& is,
1117)
1118{
1119 is.readBegin("tetrahedron");
1120 is >> t.a_ >> t.b_ >> t.c_ >> t.d_;
1121 is.readEnd("tetrahedron");
1122
1123 is.check(FUNCTION_NAME);
1124
1125 return is;
1126}
1127
1128
1129template<class Point, class PointRef>
1130inline Foam::Ostream& Foam::operator<<
1131(
1132 Ostream& os,
1134)
1135{
1136 // Format like FixedList
1137 os << nl
1139 << t.a_ << token::SPACE
1140 << t.b_ << token::SPACE
1141 << t.c_ << token::SPACE
1142 << t.d_
1143 << token::END_LIST;
1144
1145 return os;
1146}
1147
1148
1149// ************************************************************************* //
CGAL::Point_3< K > Point
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
const Cmpt & b() const noexcept
const Cmpt & d() const noexcept
const Cmpt & a() const noexcept
const Cmpt & c() const noexcept
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
label fcIndex(const label i) const noexcept
Return the forward circular index, i.e. next index which returns to the first at the end of the list.
Definition FixedListI.H:197
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
const point_type & point() const noexcept
Return the point, no checks.
Definition pointHit.H:161
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
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector.H:135
const Cmpt & z() const noexcept
Access to the vector z component.
Definition Vector.H:145
const Cmpt & y() const noexcept
Access to the vector y component.
Definition Vector.H:140
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition plane.H:91
Tet point storage. Default constructable (tetrahedron is not).
Definition tetrahedron.H:85
void flip()
Invert tetrahedron by swapping third and fourth vertices.
const point & c() const noexcept
The third vertex.
const point & d() const noexcept
The fourth vertex.
const point & b() const noexcept
The second vertex.
Pair< point > box() const
The enclosing (bounding) box for the tetrahedron.
tetPointRef tet() const
Return as tetrahedron reference.
const point & a() const noexcept
The first vertex.
tetPoints()=default
Default construct.
treeBoundBox bounds() const
The bounding box for the tetrahedron.
A tetrahedron primitive.
Point circumCentre() const
Return circum-centre.
barycentric pointToBarycentric(const point &pt) const
Calculate the barycentric coordinates from the given point.
vector Sd() const
Face area normal for side d.
triPointRef tri(const label facei) const
Return i-th face.
scalar circumRadius() const
Return circum-radius.
Point centre() const
Return centre (centroid).
vector Sc() const
Face area normal for side c.
bool inside(const point &pt) const
Return true if point is inside tetrahedron.
vector Sb() const
Face area normal for side b.
tetrahedron(const Point &p0, const Point &p1, const Point &p2, const Point &p3)
Construct from four points.
const point & b() const noexcept
scalar mag() const
Return volume.
pointHit nearestPoint(const point &p) const
Return nearest point to p on tetrahedron. Is p itself.
Point randomPoint(Random &rndGen) const
Return a random point in the tetrahedron from a uniform distribution.
Pair< Point > box() const
The enclosing (bounding) box for the tetrahedron.
FixedList< tetPoints, 200 > tetIntersectionList
Storage type for tets originating from intersecting tets.
void sliceWithPlane(const plane &pl, AboveTetOp &aboveOp, BelowTetOp &belowOp) const
Decompose tet into tets above and below plane.
Point barycentricToPoint(const barycentric &bary) const
Calculate the point from the given barycentric coordinates.
treeBoundBox bounds() const
The bounding box for the tetrahedron.
scalar quality() const
Return quality: Ratio of tetrahedron and circum-sphere volume, scaled so that a regular tetrahedron h...
const point & c() const noexcept
const point & a() const noexcept
vector Sa() const
Face area normal for side a.
const Point & d() const noexcept
Return vertex d.
static Pair< Point > box(const Point &p0, const Point &p1, const Point &p2, const Point &p3)
The enclosing (bounding) box for four points.
@ 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
Standard boundBox with extra functionality for use in octree.
A triangle primitive used to calculate face normals and swept volumes. Uses referred points.
Definition triangle.H:234
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
vector areaNormal() const
The area normal - with magnitude equal to area of triangle.
Definition triangleI.H:180
pointHit nearestPoint(const point &p) const
Return nearest point to p on triangle.
Definition triangleI.H:903
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
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
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
#define FUNCTION_NAME
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
barycentric barycentric01(Random &rndGen)
Generate a random barycentric coordinate within the unit tetrahedron.
Definition barycentric.C:66
dimensionedScalar pos0(const dimensionedScalar &ds)
PointHit< point > pointHit
A PointHit with a 3D point.
Definition pointHit.H:51
dimensionedScalar pow3(const dimensionedScalar &ds)
tetrahedron< point, const point & > tetPointRef
A tetrahedron using referred points.
Definition tetrahedron.H:72
Tensor< scalar > tensor
Definition symmTensor.H:57
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
errorManip< error > abort(error &err)
Definition errorManip.H:139
vector point
Point is a vector.
Definition point.H:37
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition triangleFwd.H:46
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition barycentric.H:45
Vector< scalar > vector
Definition vector.H:57
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
const pointField & pts
volScalarField & b
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
tetIntersectionList & tets_
storeOp(tetIntersectionList &, label &)
Random rndGen