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