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 -------------------------------------------------------------------------------
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 
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 
86 template<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 
100 template<class Point, class PointRef>
102 (
103  const FixedList<Point, 3>& pts
104 )
105 :
106  a_(pts.template get<0>()),
107  b_(pts.template get<1>()),
108  c_(pts.template get<2>())
109 {}
110 
111 
112 template<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 
125 template<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 
140 template<class Point, class PointRef>
142 {
143  is >> *this;
144 }
145 
146 
147 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
148 
149 template<class Point, class PointRef>
151 (
152  const Point& p0,
153  const Point& p1,
154  const Point& p2
155 )
156 {
157  return (1.0/3.0)*(p0 + p1 + p2);
158 }
159 
160 
161 template<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 
174 template<class Point, class PointRef>
176 (
177  const Point& p0,
178  const Point& p1,
179  const Point& p2
180 )
181 {
182  return 0.5*((p1 - p0)^(p2 - p0));
183 }
184 
185 
186 template<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 
199 template<class Point, class PointRef>
201 (
202  const Point& p0,
203  const Point& p1,
204  const Point& p2
205 )
206 {
207  vector n(areaNormal(p0, p1, p2));
208  (void) n.normalise(ROOTVSMALL);
209  return n;
210 }
211 
212 
213 template<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 
226 template<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 
242 template<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 
255 template<class Point, class PointRef>
257 {
258  return Point(c_ - b_);
259 }
260 
261 template<class Point, class PointRef>
263 {
264  return Point(a_ - c_);
265 }
266 
267 template<class Point, class PointRef>
269 {
270  return Point(b_ - a_);
271 }
272 
274 inline Foam::vector Foam::triPoints::vecA() const
275 {
276  return vector(c() - b());
277 }
278 
280 inline Foam::vector Foam::triPoints::vecB() const
281 {
282  return vector(a() - c());
283 }
284 
285 
286 inline Foam::vector Foam::triPoints::vecC() const
287 {
288  return vector(b() - a());
289 }
290 
291 
292 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
295 {
296  return triPointRef(a(), b(), c());
297 }
298 
299 
300 inline void Foam::triPoints::flip()
301 {
302  // swap pt1 <-> pt2
303  point t(b());
304  b() = c();
305  c() = t;
306 }
307 
308 
309 inline Foam::scalar Foam::triPoints::mag() const
310 {
311  return ::Foam::mag(areaNormal());
312 }
313 
314 
315 template<class Point, class PointRef>
316 inline Foam::scalar Foam::triangle<Point, PointRef>::mag() const
317 {
318  return ::Foam::mag(areaNormal());
319 }
320 
321 
322 template<class Point, class PointRef>
324 {
325  scalar d1 = (c_ - a_) & (b_ - a_);
326  scalar d2 = -(c_ - b_) & (b_ - a_);
327  scalar d3 = (c_ - a_) & (c_ - b_);
328 
329  scalar c1 = d2*d3;
330  scalar c2 = d3*d1;
331  scalar c3 = d1*d2;
332 
333  scalar c = c1 + c2 + c3;
334 
335  if (Foam::mag(c) < ROOTVSMALL)
336  {
337  // Degenerate triangle, returning centre instead of circumCentre.
338 
339  return centre();
340  }
341 
342  return
343  (
344  ((c2 + c3)*a_ + (c3 + c1)*b_ + (c1 + c2)*c_)/(2*c)
345  );
346 }
347 
348 
349 template<class Point, class PointRef>
350 inline Foam::scalar Foam::triangle<Point, PointRef>::circumRadius() const
351 {
352  const scalar d1 = (c_ - a_) & (b_ - a_);
353  const scalar d2 = -(c_ - b_) & (b_ - a_);
354  const scalar d3 = (c_ - a_) & (c_ - b_);
355 
356  const scalar denom = d2*d3 + d3*d1 + d1*d2;
357 
358  if (Foam::mag(denom) < VSMALL)
359  {
360  // Degenerate triangle, returning GREAT for circumRadius.
361 
362  return GREAT;
363  }
365  const scalar a = (d1 + d2)*(d2 + d3)*(d3 + d1) / denom;
366  return 0.5*Foam::sqrt(clamp(a, scalar(0), scalar(GREAT)));
367 }
368 
369 
370 template<class Point, class PointRef>
371 inline Foam::scalar Foam::triangle<Point, PointRef>::quality() const
372 {
373  scalar c = circumRadius();
374 
375  if (c < ROOTVSMALL)
376  {
377  // zero circumRadius, something has gone wrong.
378  return SMALL;
379  }
381  return mag()/(Foam::sqr(c)*3.0*sqrt(3.0)/4.0);
382 }
383 
384 
385 template<class Point, class PointRef>
387 (
388  const triangle& t
389 ) const
390 {
391  return (1.0/12.0)*
392  (
393  ((t.a_ - a_) & ((b_ - a_)^(c_ - a_)))
394  + ((t.b_ - b_) & ((c_ - b_)^(t.a_ - b_)))
395  + ((c_ - t.c_) & ((t.b_ - t.c_)^(t.a_ - t.c_)))
396 
397  + ((t.a_ - a_) & ((b_ - a_)^(c_ - a_)))
398  + ((b_ - t.b_) & ((t.a_ - t.b_)^(t.c_ - t.b_)))
399  + ((c_ - t.c_) & ((b_ - t.c_)^(t.a_ - t.c_)))
400  );
401 }
402 
403 
404 template<class Point, class PointRef>
406 (
407  PointRef refPt,
408  scalar density
409 ) const
410 {
411  Point aRel = a_ - refPt;
412  Point bRel = b_ - refPt;
413  Point cRel = c_ - refPt;
414 
415  tensor V
416  (
417  aRel.x(), aRel.y(), aRel.z(),
418  bRel.x(), bRel.y(), bRel.z(),
419  cRel.x(), cRel.y(), cRel.z()
420  );
421 
422  scalar a = Foam::mag((b_ - a_)^(c_ - a_));
423 
424  tensor S = 1/24.0*(tensor::one + I);
425 
426  return
427  (
428  a*I/24.0*
429  (
430  (aRel & aRel)
431  + (bRel & bRel)
432  + (cRel & cRel)
433  + ((aRel + bRel + cRel) & (aRel + bRel + cRel))
434  )
435  - a*(V.T() & S & V)
436  )
437  *density;
438 }
439 
440 
441 template<class Point, class PointRef>
443 {
444  return barycentricToPoint(barycentric2D01(rndGen));
445 }
446 
447 
448 template<class Point, class PointRef>
450 (
451  const barycentric2D& bary
452 ) const
453 {
454  return bary[0]*a_ + bary[1]*b_ + bary[2]*c_;
455 }
456 
457 
458 template<class Point, class PointRef>
460 (
461  const point& pt
462 ) const
463 {
464  barycentric2D bary;
465  pointToBarycentric(pt, bary);
466  return bary;
467 }
468 
469 
470 template<class Point, class PointRef>
472 (
473  const point& pt,
474  barycentric2D& bary
475 ) const
476 {
477  // Reference:
478  // Real-time collision detection, Christer Ericson, 2005, p47-48
479 
480  vector v0 = b_ - a_;
481  vector v1 = c_ - a_;
482  vector v2 = pt - a_;
483 
484  scalar d00 = v0 & v0;
485  scalar d01 = v0 & v1;
486  scalar d11 = v1 & v1;
487  scalar d20 = v2 & v0;
488  scalar d21 = v2 & v1;
489 
490  scalar denom = d00*d11 - d01*d01;
491 
492  if (Foam::mag(denom) < SMALL)
493  {
494  // Degenerate triangle, returning 1/3 barycentric coordinates.
495 
496  bary = barycentric2D(1.0/3.0, 1.0/3.0, 1.0/3.0);
497 
498  return denom;
499  }
500 
501  bary[1] = (d11*d20 - d01*d21)/denom;
502  bary[2] = (d00*d21 - d01*d20)/denom;
503  bary[0] = 1.0 - bary[1] - bary[2];
505  return denom;
506 }
507 
508 
509 template<class Point, class PointRef>
511 (
512  const point& origin,
513  const vector& normal
514 ) const
515 {
516  // Check points cut below(1) or above(2). Mix of above/below == 3
517  // Cf. plane::whichSide()
518  return
519  (
520  (
521  (((a_ - origin) & normal) < 0 ? 1 : 2)
522  | (((b_ - origin) & normal) < 0 ? 1 : 2)
523  | (((c_ - origin) & normal) < 0 ? 1 : 2)
524  ) == 3
525  );
526 }
527 
528 
529 template<class Point, class PointRef>
531 (
532  const point& origin,
533  const vector::components axis
534 ) const
535 {
536  // Direct check of points cut below(1) or above(2)
537  // Cf. plane::whichSide()
538  return
539  (
540  (
541  (a_[axis] < origin[axis] ? 1 : 2)
542  | (b_[axis] < origin[axis] ? 1 : 2)
543  | (c_[axis] < origin[axis] ? 1 : 2)
544  ) == 3
545  );
546 }
547 
548 
549 template<class Point, class PointRef>
551 (
552  const point& p,
553  const vector& q,
554  const intersection::algorithm alg,
555  const intersection::direction dir
556 ) const
557 {
558  // Express triangle in terms of baseVertex (point a_) and
559  // two edge vectors
560  vector E0 = b_ - a_;
561  vector E1 = c_ - a_;
562 
563  // Initialize intersection to miss.
564  pointHit inter(p);
565 
566  vector n(0.5*(E0 ^ E1));
567  scalar area = Foam::mag(n);
568 
569  if (area < VSMALL)
570  {
571  // Ineligible miss.
572  inter.setMiss(false);
573 
574  // The miss point is the nearest point on the triangle. Take any one.
575  inter.setPoint(a_);
576 
577  // The distance to the miss is the distance between the
578  // original point and plane of intersection. No normal so use
579  // distance from p to a_
580  inter.setDistance(Foam::mag(a_ - p));
581 
582  return inter;
583  }
584 
585  vector q1 = q/Foam::mag(q);
586 
587  if (dir == intersection::CONTACT_SPHERE)
588  {
589  n /= area;
590 
591  return ray(p, q1 - n, alg, intersection::VECTOR);
592  }
593 
594  // Intersection point with triangle plane
595  point pInter;
596  // Is intersection point inside triangle
597  bool hit;
598  {
599  // Reuse the fast ray intersection routine below in FULL_RAY
600  // mode since the original intersection routine has rounding problems.
601  pointHit fastInter = intersection(p, q1, intersection::FULL_RAY);
602  hit = fastInter.hit();
603 
604  if (hit)
605  {
606  pInter = fastInter.point();
607  }
608  else
609  {
610  // Calculate intersection of ray with triangle plane
611  vector v = a_ - p;
612  pInter = p + (q1&v)*q1;
613  }
614  }
615 
616  // Distance to intersection point
617  scalar dist = q1 & (pInter - p);
618 
619  const scalar planarPointTol =
620  Foam::min
621  (
622  Foam::min
623  (
624  Foam::mag(E0),
625  Foam::mag(E1)
626  ),
627  Foam::mag(c_ - b_)
629 
630  bool eligible =
632  || (alg == intersection::HALF_RAY && dist > -planarPointTol)
633  || (
634  alg == intersection::VISIBLE
635  && ((q1 & areaNormal()) < -VSMALL)
636  );
637 
638  if (hit && eligible)
639  {
640  // Hit. Set distance to intersection.
641  inter.hitPoint(pInter);
642  inter.setDistance(dist);
643  }
644  else
645  {
646  // Miss or ineligible hit.
647  inter.setMiss(eligible);
648 
649  // The miss point is the nearest point on the triangle
650  inter.setPoint(nearestPoint(p).point());
651 
652  // The distance to the miss is the distance between the
653  // original point and plane of intersection
654  inter.setDistance(Foam::mag(pInter - p));
655  }
656 
657 
658  return inter;
659 }
660 
661 
662 // From "Fast, Minimum Storage Ray/Triangle Intersection"
663 // Moeller/Trumbore.
664 template<class Point, class PointRef>
666 (
667  const point& orig,
668  const vector& dir,
669  const intersection::algorithm alg,
670  const scalar tol
671 ) const
672 {
673  const vector edge1 = b_ - a_;
674  const vector edge2 = c_ - a_;
675 
676  // begin calculating determinant - also used to calculate U parameter
677  const vector pVec = dir ^ edge2;
678 
679  // if determinant is near zero, ray lies in plane of triangle
680  const scalar det = edge1 & pVec;
681 
682  // Initialise to miss
683  pointHit intersection(false, Zero, GREAT, false);
684 
685  if (alg == intersection::VISIBLE)
686  {
687  // Culling branch
688  if (det < ROOTVSMALL)
689  {
690  // Ray on wrong side of triangle. Return miss
691  return intersection;
692  }
693  }
694  else if (alg == intersection::HALF_RAY || alg == intersection::FULL_RAY)
695  {
696  // Non-culling branch
697  if (det > -ROOTVSMALL && det < ROOTVSMALL)
698  {
699  // Ray parallel to triangle. Return miss
700  return intersection;
701  }
702  }
703 
704  const scalar inv_det = 1.0 / det;
705 
706  /* calculate distance from a_ to ray origin */
707  const vector tVec = orig-a_;
708 
709  /* calculate U parameter and test bounds */
710  const scalar u = (tVec & pVec)*inv_det;
711 
712  if (u < -tol || u > 1.0+tol)
713  {
714  // return miss
715  return intersection;
716  }
717 
718  /* prepare to test V parameter */
719  const vector qVec = tVec ^ edge1;
720 
721  /* calculate V parameter and test bounds */
722  const scalar v = (dir & qVec) * inv_det;
723 
724  if (v < -tol || u + v > 1.0+tol)
725  {
726  // return miss
727  return intersection;
728  }
729 
730  /* calculate t, scale parameters, ray intersects triangle */
731  const scalar t = (edge2 & qVec) * inv_det;
732 
733  if (alg == intersection::HALF_RAY && t < -tol)
734  {
735  // Wrong side of orig. Return miss
736  return intersection;
737  }
738 
739  intersection.hitPoint(a_ + u*edge1 + v*edge2);
740  intersection.setDistance(t);
742  return intersection;
743 }
744 
745 
746 template<class Point, class PointRef>
748 (
749  const point& p,
750  label& nearType,
751  label& nearLabel
752 ) const
753 {
754  // Adapted from:
755  // Real-time collision detection, Christer Ericson, 2005, p136-142
756 
757  // Check if P in vertex region outside A
758  vector ab = b_ - a_;
759  vector ac = c_ - a_;
760  vector ap = p - a_;
761 
762  scalar d1 = ab & ap;
763  scalar d2 = ac & ap;
764 
765  if (d1 <= 0.0 && d2 <= 0.0)
766  {
767  // barycentric coordinates (1, 0, 0)
768 
769  nearType = POINT;
770  nearLabel = 0;
771  return pointHit(false, a_, Foam::mag(a_ - p), true);
772  }
773 
774  // Check if P in vertex region outside B
775  vector bp = p - b_;
776  scalar d3 = ab & bp;
777  scalar d4 = ac & bp;
778 
779  if (d3 >= 0.0 && d4 <= d3)
780  {
781  // barycentric coordinates (0, 1, 0)
782 
783  nearType = POINT;
784  nearLabel = 1;
785  return pointHit(false, b_, Foam::mag(b_ - p), true);
786  }
787 
788  // Check if P in edge region of AB, if so return projection of P onto AB
789  scalar vc = d1*d4 - d3*d2;
790 
791  if (vc <= 0.0 && d1 >= 0.0 && d3 <= 0.0)
792  {
793  if ((d1 - d3) < ROOTVSMALL)
794  {
795  // Degenerate triangle, for d1 = d3, a_ and b_ are likely coincident
796  nearType = POINT;
797  nearLabel = 0;
798  return pointHit(false, a_, Foam::mag(a_ - p), true);
799  }
800 
801  // barycentric coordinates (1-v, v, 0)
802  scalar v = d1/(d1 - d3);
803 
804  point nearPt = a_ + v*ab;
805  nearType = EDGE;
806  nearLabel = 0;
807  return pointHit(false, nearPt, Foam::mag(nearPt - p), true);
808  }
809 
810  // Check if P in vertex region outside C
811  vector cp = p - c_;
812  scalar d5 = ab & cp;
813  scalar d6 = ac & cp;
814 
815  if (d6 >= 0.0 && d5 <= d6)
816  {
817  // barycentric coordinates (0, 0, 1)
818 
819  nearType = POINT;
820  nearLabel = 2;
821  return pointHit(false, c_, Foam::mag(c_ - p), true);
822  }
823 
824  // Check if P in edge region of AC, if so return projection of P onto AC
825  scalar vb = d5*d2 - d1*d6;
826 
827  if (vb <= 0.0 && d2 >= 0.0 && d6 <= 0.0)
828  {
829  if ((d2 - d6) < ROOTVSMALL)
830  {
831  // Degenerate triangle, for d2 = d6, a_ and c_ are likely coincident
832  nearType = POINT;
833  nearLabel = 0;
834  return pointHit(false, a_, Foam::mag(a_ - p), true);
835  }
836 
837  // barycentric coordinates (1-w, 0, w)
838  scalar w = d2/(d2 - d6);
839 
840  point nearPt = a_ + w*ac;
841  nearType = EDGE;
842  nearLabel = 2;
843  return pointHit(false, nearPt, Foam::mag(nearPt - p), true);
844  }
845 
846  // Check if P in edge region of BC, if so return projection of P onto BC
847  scalar va = d3*d6 - d5*d4;
848 
849  if (va <= 0.0 && (d4 - d3) >= 0.0 && (d5 - d6) >= 0.0)
850  {
851  if (((d4 - d3) + (d5 - d6)) < ROOTVSMALL)
852  {
853  // Degenerate triangle, for (d4 - d3) = (d6 - d5), b_ and c_ are
854  // likely coincident
855  nearType = POINT;
856  nearLabel = 1;
857  return pointHit(false, b_, Foam::mag(b_ - p), true);
858  }
859 
860  // barycentric coordinates (0, 1-w, w)
861  scalar w = (d4 - d3)/((d4 - d3) + (d5 - d6));
862 
863  point nearPt = b_ + w*(c_ - b_);
864  nearType = EDGE;
865  nearLabel = 1;
866  return pointHit(false, nearPt, Foam::mag(nearPt - p), true);
867  }
868 
869  // P inside face region. Compute Q through its barycentric
870  // coordinates (u, v, w)
871 
872  if ((va + vb + vc) < ROOTVSMALL)
873  {
874  // Degenerate triangle, return the centre because no edge or points are
875  // closest
876  point nearPt = centre();
877  nearType = NONE,
878  nearLabel = -1;
879  return pointHit(true, nearPt, Foam::mag(nearPt - p), false);
880  }
881 
882  scalar denom = 1.0/(va + vb + vc);
883  scalar v = vb * denom;
884  scalar w = vc * denom;
885 
886  // = u*a + v*b + w*c, u = va*denom = 1.0 - v - w
887 
888  point nearPt = a_ + ab*v + ac*w;
889  nearType = NONE,
890  nearLabel = -1;
891  return pointHit(true, nearPt, Foam::mag(nearPt - p), false);
892 }
893 
894 
895 template<class Point, class PointRef>
897 (
898  const point& p
899 ) const
900 {
901  // Dummy labels
902  label nearType = -1;
903  label nearLabel = -1;
905  return nearestPointClassify(p, nearType, nearLabel);
906 }
907 
908 
909 template<class Point, class PointRef>
911 (
912  const point& p,
913  label& nearType,
914  label& nearLabel
915 ) const
916 {
917  return nearestPointClassify(p, nearType, nearLabel).hit();
918 }
919 
920 
921 template<class Point, class PointRef>
923 (
924  const linePointRef& ln,
925  pointHit& lnInfo
926 ) const
927 {
928  vector q = ln.vec();
929  pointHit triInfo
930  (
932  (
933  ln.start(),
934  q,
936  )
937  );
938 
939  if (triInfo.hit())
940  {
941  // Line hits triangle. Find point on line.
942  if (triInfo.distance() > 1)
943  {
944  // Hit beyond endpoint
945  lnInfo.setMiss(true);
946  lnInfo.setPoint(ln.end());
947  scalar dist = triInfo.point().dist(lnInfo.point());
948  lnInfo.setDistance(dist);
949  triInfo.setMiss(true);
950  triInfo.setDistance(dist);
951  }
952  else if (triInfo.distance() < 0)
953  {
954  // Hit beyond startpoint
955  lnInfo.setMiss(true);
956  lnInfo.setPoint(ln.start());
957  scalar dist = triInfo.point().dist(lnInfo.point());
958  lnInfo.setDistance(dist);
959  triInfo.setMiss(true);
960  triInfo.setDistance(dist);
961  }
962  else
963  {
964  // Hit on line
965  lnInfo.hitPoint(triInfo.point());
966  lnInfo.setDistance(0);
967  triInfo.setDistance(0);
968  }
969  }
970  else
971  {
972  // Line skips triangle. See which triangle edge it gets closest to
973 
974  point nearestEdgePoint;
975  point nearestLinePoint;
976  //label minEdgeIndex = 0;
977  scalar minDist = ln.nearestDist
978  (
979  linePointRef(a_, b_),
980  nearestLinePoint,
981  nearestEdgePoint
982  );
983 
984  {
985  point linePoint;
986  point triEdgePoint;
987  scalar dist = ln.nearestDist
988  (
989  linePointRef(b_, c_),
990  linePoint,
991  triEdgePoint
992  );
993  if (dist < minDist)
994  {
995  minDist = dist;
996  nearestEdgePoint = triEdgePoint;
997  nearestLinePoint = linePoint;
998  //minEdgeIndex = 1;
999  }
1000  }
1001 
1002  {
1003  point linePoint;
1004  point triEdgePoint;
1005  scalar dist = ln.nearestDist
1006  (
1007  linePointRef(c_, a_),
1008  linePoint,
1009  triEdgePoint
1010  );
1011  if (dist < minDist)
1012  {
1013  minDist = dist;
1014  nearestEdgePoint = triEdgePoint;
1015  nearestLinePoint = linePoint;
1016  //minEdgeIndex = 2;
1017  }
1018  }
1019 
1020  lnInfo.setDistance(minDist);
1021  triInfo.setDistance(minDist);
1022  triInfo.setMiss(false);
1023  triInfo.setPoint(nearestEdgePoint);
1024 
1025  // Convert point on line to pointHit
1026  if (Foam::mag(nearestLinePoint-ln.start()) < SMALL)
1027  {
1028  lnInfo.setMiss(true);
1029  lnInfo.setPoint(ln.start());
1030  }
1031  else if (Foam::mag(nearestLinePoint-ln.end()) < SMALL)
1032  {
1033  lnInfo.setMiss(true);
1034  lnInfo.setPoint(ln.end());
1035  }
1036  else
1037  {
1038  lnInfo.hitPoint(nearestLinePoint);
1039  }
1040  }
1041  return triInfo;
1042 }
1043 
1044 
1045 template<class Point, class PointRef>
1047 (
1048  const point& p,
1049  const scalar tol
1050 ) const
1051 {
1052  const scalar dist = ((p - a_) & unitNormal());
1054  return ((dist < -tol) ? -1 : (dist > tol) ? +1 : 0);
1055 }
1056 
1057 
1058 template<class Point, class PointRef>
1061  const triPoints&
1062 )
1063 {}
1064 
1065 
1066 template<class Point, class PointRef>
1069  area_(0.0)
1070 {}
1071 
1072 
1073 template<class Point, class PointRef>
1075 (
1076  const triPoints& tri
1077 )
1079  area_ += triangle<Point, const Point&>(tri).mag();
1080 }
1081 
1082 
1083 template<class Point, class PointRef>
1085 (
1086  triIntersectionList& tris,
1087  label& nTris
1088 )
1089 :
1090  tris_(tris),
1091  nTris_(nTris)
1092 {}
1093 
1094 
1095 template<class Point, class PointRef>
1097 (
1098  const triPoints& tri
1099 )
1100 {
1101  tris_[nTris_++] = tri;
1102 }
1103 
1104 
1105 template<class Point, class PointRef>
1107 (
1108  const FixedList<scalar, 3>& d,
1109  const triPoints& t,
1110  const label negI,
1111  const label posI
1112 )
1113 {
1114  return (d[posI]*t[negI] - d[negI]*t[posI])/(-d[negI] + d[posI]);
1115 }
1116 
1117 
1118 // * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
1119 
1120 template<class Point, class PointRef>
1121 inline Foam::Istream& Foam::operator>>
1122 (
1123  Istream& is,
1124  triangle<Point, PointRef>& t
1125 )
1126 {
1127  is.readBegin("triangle");
1128  is >> t.a_ >> t.b_ >> t.c_;
1129  is.readEnd("triangle");
1130 
1131  is.check(FUNCTION_NAME);
1132  return is;
1133 }
1134 
1135 
1136 template<class Point, class PointRef>
1137 inline Foam::Ostream& Foam::operator<<
1138 (
1139  Ostream& os,
1140  const triangle<Point, PointRef>& t
1141 )
1142 {
1143  // Format like FixedList
1144  os << nl
1146  << t.a_ << token::SPACE
1147  << t.b_ << token::SPACE
1148  << t.c_
1149  << token::END_LIST;
1150 
1151  return os;
1152 }
1153 
1154 
1155 // ************************************************************************* //
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:380
const point & c() const noexcept
The third vertex.
Definition: triangle.H:151
A triangle primitive used to calculate face normals and swept volumes. Uses referred points...
Definition: triangle.H:228
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:453
scalar mag() const
The magnitude of the triangle area.
Definition: triangleI.H:309
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:279
scalar quality() const
Return quality: Ratio of triangle and circum-circle.
Definition: triangleI.H:364
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:236
vector vecB() const
Edge vector opposite point b(): from c() to a()
Definition: triangleI.H:273
tensor inertia(PointRef refPt=Zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
Definition: triangleI.H:399
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:741
Begin list [isseparator].
Definition: token.H:161
triPointRef tri() const
Return as triangle reference.
Definition: triangleI.H:287
void flip()
Flip triangle orientation by swapping second and third vertices.
Definition: triangleI.H:293
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:1040
Point vecB() const
Edge vector opposite point b(): from c() to a()
Definition: triangleI.H:255
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:81
bool classify(const point &p, label &nearType, label &nearLabel) const
Classify nearest point to p in triangle plane.
Definition: triangleI.H:904
pointHit nearestPoint(const point &p) const
Return nearest point to p on triangle.
Definition: triangleI.H:890
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:207
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:544
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:249
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:155
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:302
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:343
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:1078
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:186
const point & a() const noexcept
The first vertex.
Definition: triangle.H:141
Point barycentricToPoint(const barycentric2D &bary) const
Calculate the point from the given barycentric coordinates.
Definition: triangleI.H:443
bool intersects(const point &origin, const vector &normal) const
Fast intersection detection with a plane.
Definition: triangleI.H:504
CGAL::Point_3< K > Point
Pair< point > box() const
The enclosing (bounding) box for the triangle.
Definition: triangleI.H:242
Point vecC() const
Edge vector opposite point c(): from a() to b()
Definition: triangleI.H:261
Point circumCentre() const
Return circum-centre.
Definition: triangleI.H:316
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:267
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:435
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:659
Tensor of scalars, i.e. Tensor<scalar>.
const point & b() const noexcept
The second vertex.
Definition: triangle.H:146
vector unitNormal() const
Return unit normal.
Definition: triangleI.H:213
const volScalarField & p0
Definition: EEqn.H:36
point centre() const
Return centre (centroid)
Definition: triangleI.H:161
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:180
const pointField & pts
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127