triangleI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "IOstreams.H"
27 #include "pointHit.H"
28 #include "mathematicalConstants.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class Point, class PointRef>
34 (
35  const Point& a,
36  const Point& b,
37  const Point& c
38 )
39 :
40  a_(a),
41  b_(b),
42  c_(c)
43 {}
44 
45 
46 template<class Point, class PointRef>
48 (
49  const UList<Point>& points,
50  const FixedList<label, 3>& indices
51 )
52 :
53  a_(points[indices[0]]),
54  b_(points[indices[1]]),
55  c_(points[indices[2]])
56 {}
57 
58 
59 
60 template<class Point, class PointRef>
62 {
63  is >> *this;
64 }
65 
66 
67 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68 
69 template<class Point, class PointRef>
70 inline const Point& Foam::triangle<Point, PointRef>::a() const
71 {
72  return a_;
73 }
74 
75 template<class Point, class PointRef>
76 inline const Point& Foam::triangle<Point, PointRef>::b() const
77 {
78  return b_;
79 }
80 
81 template<class Point, class PointRef>
82 inline const Point& Foam::triangle<Point, PointRef>::c() const
83 {
84  return c_;
85 }
86 
87 
88 template<class Point, class PointRef>
90 {
91  return (1.0/3.0)*(a_ + b_ + c_);
92 }
93 
94 
95 template<class Point, class PointRef>
97 {
98  return 0.5*((b_ - a_)^(c_ - a_));
99 }
100 
101 
102 template<class Point, class PointRef>
103 inline Foam::scalar Foam::triangle<Point, PointRef>::mag() const
104 {
105  return Foam::mag(area());
106 }
107 
108 
109 template<class Point, class PointRef>
111 {
112  const vector a = area();
113  const scalar maga = Foam::mag(a);
114  return maga > 0 ? a/maga : Zero;
115 }
116 
117 
118 template<class Point, class PointRef>
121 {
122  const scalar d1 = (c_ - a_) & (b_ - a_);
123  const scalar d2 = -(c_ - b_) & (b_ - a_);
124  const scalar d3 = (c_ - a_) & (c_ - b_);
125 
126  const scalar c1 = d2*d3;
127  const scalar c2 = d3*d1;
128  const scalar c3 = d1*d2;
129 
130  const scalar denom = c1 + c2 + c3;
131 
132  if (Foam::mag(denom) < rootVSmall)
133  {
134  // Degenerate. Return a point far away.
135  static const scalar sqrt3 = sqrt(scalar(3));
136  return Tuple2<Point, scalar>(point::uniform(great), sqrt3*great);
137  }
138  else
139  {
140  return
142  (
143  ((c2 + c3)*a_ + (c3 + c1)*b_ + (c1 + c2)*c_)/(2*denom),
144  sqrt(max(0, (d1 + d2)*(d2 + d3)*(d3 + d1)/(4*denom)))
145  );
146  }
147 }
148 
149 
150 template<class Point, class PointRef>
151 inline Foam::scalar Foam::triangle<Point, PointRef>::quality() const
152 {
153  const scalar r = circumCircle().second();
154  static const scalar sqrt3 = sqrt(scalar(3));
155  return mag()/(0.75*sqrt3*sqr(min(r, great)) + rootVSmall);
156 }
157 
158 
159 template<class Point, class PointRef>
161 (
162  const triangle& t
163 ) const
164 {
165  return (1.0/12.0)*
166  (
167  ((t.a_ - a_) & ((b_ - a_)^(c_ - a_)))
168  + ((t.b_ - b_) & ((c_ - b_)^(t.a_ - b_)))
169  + ((c_ - t.c_) & ((t.b_ - t.c_)^(t.a_ - t.c_)))
170 
171  + ((t.a_ - a_) & ((b_ - a_)^(c_ - a_)))
172  + ((b_ - t.b_) & ((t.a_ - t.b_)^(t.c_ - t.b_)))
173  + ((c_ - t.c_) & ((b_ - t.c_)^(t.a_ - t.c_)))
174  );
175 }
176 
177 
178 template<class Point, class PointRef>
180 (
181  PointRef refPt,
182  scalar density
183 ) const
184 {
185  const Point aRel = a_ - refPt;
186  const Point bRel = b_ - refPt;
187  const Point cRel = c_ - refPt;
188 
189  const tensor V
190  (
191  aRel.x(), aRel.y(), aRel.z(),
192  bRel.x(), bRel.y(), bRel.z(),
193  cRel.x(), cRel.y(), cRel.z()
194  );
195 
196  const scalar a = Foam::mag((b_ - a_)^(c_ - a_));
197  const symmTensor S = 1/24.0*(symmTensor::one + I);
198 
199  return
200  (
201  a*I/24.0*
202  (
203  (aRel & aRel)
204  + (bRel & bRel)
205  + (cRel & cRel)
206  + ((aRel + bRel + cRel) & (aRel + bRel + cRel))
207  )
208  - a*symm(V.T() & S & V)
209  )*density;
210 }
211 
212 
213 template<class Point, class PointRef>
215 (
217 ) const
218 {
219  return barycentricToPoint(barycentric2D01(rndGen));
220 }
221 
222 
223 template<class Point, class PointRef>
225 (
226  const barycentric2D& bary
227 ) const
228 {
229  return bary[0]*a_ + bary[1]*b_ + bary[2]*c_;
230 }
231 
232 
233 template<class Point, class PointRef>
235 (
236  const point& pt
237 ) const
238 {
239  barycentric2D bary;
240  pointToBarycentric(pt, bary);
241  return bary;
242 }
243 
244 
245 template<class Point, class PointRef>
247 (
248  const point& pt,
249  barycentric2D& bary
250 ) const
251 {
252  // Reference:
253  // Real-time collision detection, Christer Ericson, 2005, p47-48
254 
255  vector v0 = b_ - a_;
256  vector v1 = c_ - a_;
257  vector v2 = pt - a_;
258 
259  scalar d00 = v0 & v0;
260  scalar d01 = v0 & v1;
261  scalar d11 = v1 & v1;
262  scalar d20 = v2 & v0;
263  scalar d21 = v2 & v1;
264 
265  scalar denom = d00*d11 - d01*d01;
266 
267  if (Foam::mag(denom) < small)
268  {
269  // Degenerate triangle, returning 1/3 barycentric coordinates.
270 
271  bary = barycentric2D(1.0/3.0, 1.0/3.0, 1.0/3.0);
272 
273  return denom;
274  }
275 
276  bary[1] = (d11*d20 - d01*d21)/denom;
277  bary[2] = (d00*d21 - d01*d20)/denom;
278  bary[0] = 1.0 - bary[1] - bary[2];
279 
280  return denom;
281 }
282 
283 
284 template<class Point, class PointRef>
286 (
287  const point& p,
288  const vector& q,
289  const intersection::algorithm alg,
290  const intersection::direction dir
291 ) const
292 {
293  // Express triangle in terms of baseVertex (point a_) and
294  // two edge vectors
295  const vector E0 = b_ - a_;
296  const vector E1 = c_ - a_;
297 
298  // Initialise intersection to miss.
299  pointHit inter(p);
300 
301  vector n(0.5*(E0 ^ E1));
302  const scalar magArea = Foam::mag(n);
303 
304  if (magArea < vSmall)
305  {
306  // Ineligible miss.
307  inter.setMiss(false);
308 
309  // The miss point is the nearest point on the triangle. Take any one.
310  inter.setPoint(a_);
311 
312  // The distance to the miss is the distance between the
313  // original point and plane of intersection. No normal so use
314  // distance from p to a_
315  inter.setDistance(Foam::mag(a_ - p));
316 
317  return inter;
318  }
319 
320  const vector q1 = q/Foam::mag(q);
321 
323  {
324  n /= magArea;
325 
326  return ray(p, q1 - n, alg, intersection::direction::vector);
327  }
328 
329  // Intersection point with triangle plane
330  point pInter;
331 
332  // Is intersection point inside triangle
333  bool hit;
334  {
335  // Reuse the fast ray intersection routine below in algorithm::fullRay
336  // mode since the original intersection routine has rounding problems.
337  pointHit fastInter = intersection
338  (
339  p,
340  q1,
342  );
343  hit = fastInter.hit();
344 
345  if (hit)
346  {
347  pInter = fastInter.rawPoint();
348  }
349  else
350  {
351  // Calculate intersection of ray with triangle plane
352  vector v = a_ - p;
353  pInter = p + (q1&v)*q1;
354  }
355  }
356 
357  // Distance to intersection point
358  const scalar dist = q1 & (pInter - p);
359 
360  const scalar planarPointTol =
361  Foam::min
362  (
363  Foam::min
364  (
365  Foam::mag(E0),
366  Foam::mag(E1)
367  ),
368  Foam::mag(c_ - b_)
370 
371  bool eligible =
373  || (alg == intersection::algorithm::halfRay && dist > -planarPointTol)
374  || (
376  && ((q1 & area()) < -vSmall)
377  );
378 
379  if (hit && eligible)
380  {
381  // Hit. Set distance to intersection.
382  inter.setHit();
383  inter.setPoint(pInter);
384  inter.setDistance(dist);
385  }
386  else
387  {
388  // Miss or ineligible hit.
389  inter.setMiss(eligible);
390 
391  // The miss point is the nearest point on the triangle
392  inter.setPoint(nearestPoint(p).rawPoint());
393 
394  // The distance to the miss is the distance between the
395  // original point and plane of intersection
396  inter.setDistance(Foam::mag(pInter - p));
397  }
398 
399  return inter;
400 }
401 
402 
403 // From "Fast, Minimum Storage Ray/Triangle Intersection"
404 // Moeller/Trumbore.
405 template<class Point, class PointRef>
407 (
408  const point& orig,
409  const vector& dir,
410  const intersection::algorithm alg,
411  const scalar tol
412 ) const
413 {
414  const vector edge1 = b_ - a_;
415  const vector edge2 = c_ - a_;
416 
417  // Begin calculating determinant - also used to calculate U parameter
418  const vector pVec = dir ^ edge2;
419 
420  // Ff determinant is near zero, ray lies in plane of triangle
421  const scalar det = edge1 & pVec;
422 
423  // Initialise to miss
424  pointHit intersection(false, Zero, great, false);
425 
427  {
428  // Culling branch
429  if (det < rootVSmall)
430  {
431  // Ray on wrong side of triangle. Return miss
432  return intersection;
433  }
434  }
435  else if
436  (
439  )
440  {
441  // Non-culling branch
442  if (det > -rootVSmall && det < rootVSmall)
443  {
444  // Ray parallel to triangle. Return miss
445  return intersection;
446  }
447  }
448 
449  const scalar inv_det = 1.0 / det;
450 
451  // Calculate distance from a_ to ray origin
452  const vector tVec = orig-a_;
453 
454  // Calculate U parameter and test bounds
455  const scalar u = (tVec & pVec)*inv_det;
456 
457  if (u < -tol || u > 1.0+tol)
458  {
459  // return miss
460  return intersection;
461  }
462 
463  // Prepare to test V parameter
464  const vector qVec = tVec ^ edge1;
465 
466  // Calculate V parameter and test bounds
467  const scalar v = (dir & qVec) * inv_det;
468 
469  if (v < -tol || u + v > 1.0+tol)
470  {
471  // return miss
472  return intersection;
473  }
474 
475  // Calculate t, scale parameters, ray intersects triangle
476  const scalar t = (edge2 & qVec) * inv_det;
477 
478  if (alg == intersection::algorithm::halfRay && t < -tol)
479  {
480  // Wrong side of orig. Return miss
481  return intersection;
482  }
483 
484  intersection.setHit();
485  intersection.setPoint(a_ + u*edge1 + v*edge2);
486  intersection.setDistance(t);
487 
488  return intersection;
489 }
490 
491 
492 template<class Point, class PointRef>
494 (
495  const point& p,
496  label& nearType,
497  label& nearLabel
498 ) const
499 {
500  // Adapted from:
501  // Real-time collision detection, Christer Ericson, 2005, p136-142
502 
503  // Check if P in vertex region outside A
504  vector ab = b_ - a_;
505  vector ac = c_ - a_;
506  vector ap = p - a_;
507 
508  scalar d1 = ab & ap;
509  scalar d2 = ac & ap;
510 
511  if (d1 <= 0.0 && d2 <= 0.0)
512  {
513  // barycentric coordinates (1, 0, 0)
514 
515  nearType = POINT;
516  nearLabel = 0;
517  return pointHit(false, a_, Foam::mag(a_ - p), true);
518  }
519 
520  // Check if P in vertex region outside B
521  vector bp = p - b_;
522  scalar d3 = ab & bp;
523  scalar d4 = ac & bp;
524 
525  if (d3 >= 0.0 && d4 <= d3)
526  {
527  // barycentric coordinates (0, 1, 0)
528 
529  nearType = POINT;
530  nearLabel = 1;
531  return pointHit(false, b_, Foam::mag(b_ - p), true);
532  }
533 
534  // Check if P in edge region of AB, if so return projection of P onto AB
535  scalar vc = d1*d4 - d3*d2;
536 
537  if (vc <= 0.0 && d1 >= 0.0 && d3 <= 0.0)
538  {
539  if ((d1 - d3) < rootVSmall)
540  {
541  // Degenerate triangle, for d1 = d3, a_ and b_ are likely coincident
542  nearType = POINT;
543  nearLabel = 0;
544  return pointHit(false, a_, Foam::mag(a_ - p), true);
545  }
546 
547  // barycentric coordinates (1-v, v, 0)
548  scalar v = d1/(d1 - d3);
549 
550  point nearPt = a_ + v*ab;
551  nearType = EDGE;
552  nearLabel = 0;
553  return pointHit(false, nearPt, Foam::mag(nearPt - p), true);
554  }
555 
556  // Check if P in vertex region outside C
557  vector cp = p - c_;
558  scalar d5 = ab & cp;
559  scalar d6 = ac & cp;
560 
561  if (d6 >= 0.0 && d5 <= d6)
562  {
563  // barycentric coordinates (0, 0, 1)
564 
565  nearType = POINT;
566  nearLabel = 2;
567  return pointHit(false, c_, Foam::mag(c_ - p), true);
568  }
569 
570  // Check if P in edge region of AC, if so return projection of P onto AC
571  scalar vb = d5*d2 - d1*d6;
572 
573  if (vb <= 0.0 && d2 >= 0.0 && d6 <= 0.0)
574  {
575  if ((d2 - d6) < rootVSmall)
576  {
577  // Degenerate triangle, for d2 = d6, a_ and c_ are likely coincident
578  nearType = POINT;
579  nearLabel = 0;
580  return pointHit(false, a_, Foam::mag(a_ - p), true);
581  }
582 
583  // barycentric coordinates (1-w, 0, w)
584  scalar w = d2/(d2 - d6);
585 
586  point nearPt = a_ + w*ac;
587  nearType = EDGE;
588  nearLabel = 2;
589  return pointHit(false, nearPt, Foam::mag(nearPt - p), true);
590  }
591 
592  // Check if P in edge region of BC, if so return projection of P onto BC
593  scalar va = d3*d6 - d5*d4;
594 
595  if (va <= 0.0 && (d4 - d3) >= 0.0 && (d5 - d6) >= 0.0)
596  {
597  if (((d4 - d3) + (d5 - d6)) < rootVSmall)
598  {
599  // Degenerate triangle, for (d4 - d3) = (d6 - d5), b_ and c_ are
600  // likely coincident
601  nearType = POINT;
602  nearLabel = 1;
603  return pointHit(false, b_, Foam::mag(b_ - p), true);
604  }
605 
606  // barycentric coordinates (0, 1-w, w)
607  scalar w = (d4 - d3)/((d4 - d3) + (d5 - d6));
608 
609  point nearPt = b_ + w*(c_ - b_);
610  nearType = EDGE;
611  nearLabel = 1;
612  return pointHit(false, nearPt, Foam::mag(nearPt - p), true);
613  }
614 
615  // P inside face region. Compute Q through its barycentric
616  // coordinates (u, v, w)
617 
618  if ((va + vb + vc) < rootVSmall)
619  {
620  // Degenerate triangle, return the centre because no edge or points are
621  // closest
622  point nearPt = centre();
623  nearType = NONE,
624  nearLabel = -1;
625  return pointHit(true, nearPt, Foam::mag(nearPt - p), false);
626  }
627 
628  scalar denom = 1.0/(va + vb + vc);
629  scalar v = vb * denom;
630  scalar w = vc * denom;
631 
632  // = u*a + v*b + w*c, u = va*denom = 1.0 - v - w
633 
634  point nearPt = a_ + ab*v + ac*w;
635  nearType = NONE,
636  nearLabel = -1;
637  return pointHit(true, nearPt, Foam::mag(nearPt - p), false);
638 }
639 
640 
641 template<class Point, class PointRef>
643 (
644  const point& p
645 ) const
646 {
647  // Dummy labels
648  label nearType = -1;
649  label nearLabel = -1;
650 
651  return nearestPointClassify(p, nearType, nearLabel);
652 }
653 
654 
655 template<class Point, class PointRef>
657 (
658  const point& p,
659  label& nearType,
660  label& nearLabel
661 ) const
662 {
663  return nearestPointClassify(p, nearType, nearLabel).hit();
664 }
665 
666 
667 template<class Point, class PointRef>
669 (
670  const linePointRef& ln,
671  pointHit& lnInfo
672 ) const
673 {
674  vector q = ln.vec();
675  pointHit triInfo
676  (
678  (
679  ln.start(),
680  q,
682  )
683  );
684 
685  if (triInfo.hit())
686  {
687  // Line hits triangle. Find point on line.
688  if (triInfo.distance() > 1)
689  {
690  // Hit beyond endpoint
691  lnInfo.setMiss(true);
692  lnInfo.setPoint(ln.end());
693  scalar dist = Foam::mag(triInfo.hitPoint()-lnInfo.missPoint());
694  lnInfo.setDistance(dist);
695  triInfo.setMiss(true);
696  triInfo.setDistance(dist);
697  }
698  else if (triInfo.distance() < 0)
699  {
700  // Hit beyond startpoint
701  lnInfo.setMiss(true);
702  lnInfo.setPoint(ln.start());
703  scalar dist = Foam::mag(triInfo.hitPoint()-lnInfo.missPoint());
704  lnInfo.setDistance(dist);
705  triInfo.setMiss(true);
706  triInfo.setDistance(dist);
707  }
708  else
709  {
710  // Hit on line
711  lnInfo.setHit();
712  lnInfo.setPoint(triInfo.hitPoint());
713  lnInfo.setDistance(0.0);
714  triInfo.setDistance(0.0);
715  }
716  }
717  else
718  {
719  // Line skips triangle. See which triangle edge it gets closest to
720 
721  point nearestEdgePoint;
722  point nearestLinePoint;
723  // label minEdgeIndex = 0;
724  scalar minDist = ln.nearestDist
725  (
726  linePointRef(a_, b_),
727  nearestLinePoint,
728  nearestEdgePoint
729  );
730 
731  {
732  point linePoint;
733  point triEdgePoint;
734  scalar dist = ln.nearestDist
735  (
736  linePointRef(b_, c_),
737  linePoint,
738  triEdgePoint
739  );
740  if (dist < minDist)
741  {
742  minDist = dist;
743  nearestEdgePoint = triEdgePoint;
744  nearestLinePoint = linePoint;
745  // minEdgeIndex = 1;
746  }
747  }
748 
749  {
750  point linePoint;
751  point triEdgePoint;
752  scalar dist = ln.nearestDist
753  (
754  linePointRef(c_, a_),
755  linePoint,
756  triEdgePoint
757  );
758  if (dist < minDist)
759  {
760  minDist = dist;
761  nearestEdgePoint = triEdgePoint;
762  nearestLinePoint = linePoint;
763  // minEdgeIndex = 2;
764  }
765  }
766 
767  lnInfo.setDistance(minDist);
768  triInfo.setDistance(minDist);
769  triInfo.setMiss(false);
770  triInfo.setPoint(nearestEdgePoint);
771 
772  // Convert point on line to pointHit
773  if (Foam::mag(nearestLinePoint-ln.start()) < small)
774  {
775  lnInfo.setMiss(true);
776  lnInfo.setPoint(ln.start());
777  }
778  else if (Foam::mag(nearestLinePoint-ln.end()) < small)
779  {
780  lnInfo.setMiss(true);
781  lnInfo.setPoint(ln.end());
782  }
783  else
784  {
785  lnInfo.setHit();
786  lnInfo.setPoint(nearestLinePoint);
787  }
788  }
789  return triInfo;
790 }
791 
792 
793 // * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
794 
795 template<class Point, class PointRef>
796 inline Foam::Istream& Foam::operator>>
797 (
798  Istream& is,
800 )
801 {
802  is.readBegin("triangle");
803  is >> t.a_ >> t.b_ >> t.c_;
804  is.readEnd("triangle");
805 
806  is.check("Istream& operator>>(Istream&, triangle&)");
807  return is;
808 }
809 
810 
811 template<class Point, class PointRef>
812 inline Foam::Ostream& Foam::operator<<
813 (
814  Ostream& os,
815  const triangle<Point, PointRef>& t
816 )
817 {
818  os << nl
820  << t.a_ << token::SPACE
821  << t.b_ << token::SPACE
822  << t.c_
823  << token::END_LIST;
824 
825  return os;
826 }
827 
828 
829 // ************************************************************************* //
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
label n
Templated 2D Barycentric derived from VectorSpace. Has 3 components, one of which is redundant.
Definition: Barycentric2D.H:53
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
Istream & readEnd(const char *funcName)
Definition: Istream.C:103
Istream & readBegin(const char *funcName)
Definition: Istream.C:86
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
const Point & hitPoint() const
Return hit point.
Definition: PointHit.H:126
void setPoint(const Point &p)
Definition: PointHit.H:181
void setDistance(const scalar d)
Definition: PointHit.H:186
scalar distance() const
Return distance to hit.
Definition: PointHit.H:139
const Point & rawPoint() const
Return point with no checking.
Definition: PointHit.H:158
const Point & missPoint() const
Return miss point.
Definition: PointHit.H:145
void setMiss(const bool eligible)
Definition: PointHit.H:175
bool hit() const
Is there a hit.
Definition: PointHit.H:120
void setHit()
Definition: PointHit.H:169
Tensor< Cmpt > T() const
Return transpose.
Definition: TensorI.H:331
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:66
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:74
static const Form one
Definition: VectorSpace.H:119
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:168
Foam::intersection.
Definition: intersection.H:50
static scalar planarTol()
Return planar tolerance.
Definition: intersection.H:83
A line primitive.
Definition: line.H:71
Random number generator.
@ BEGIN_LIST
Definition: token.H:109
@ END_LIST
Definition: token.H:110
A triangle primitive used to calculate face areas and swept volumes.
Definition: triangle.H:80
barycentric2D pointToBarycentric(const point &pt) const
Calculate the barycentric coordinates from the given point.
Definition: triangleI.H:235
triangle(const Point &a, const Point &b, const Point &c)
Construct from three points.
Definition: triangleI.H:34
vector area() const
Return vector area.
Definition: triangleI.H:96
Point barycentricToPoint(const barycentric2D &bary) const
Calculate the point from the given barycentric coordinates.
Definition: triangleI.H:225
scalar sweptVol(const triangle &t) const
Return swept-volume.
Definition: triangleI.H:161
const Point & a() const
Return first vertex.
Definition: triangleI.H:70
Point centre() const
Return centre (centroid)
Definition: triangleI.H:89
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:494
pointHit ray(const point &p, const vector &q, const intersection::algorithm=intersection::algorithm::fullRay, const intersection::direction dir=intersection::direction::vector) const
Return point intersection with a ray.
Definition: triangleI.H:286
const Point & c() const
Return third vertex.
Definition: triangleI.H:82
scalar mag() const
Return scalar magnitude.
Definition: triangleI.H:103
symmTensor inertia(PointRef refPt=vector::zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
Definition: triangleI.H:180
pointHit nearestPoint(const point &p) const
Return nearest point to p on triangle.
Definition: triangleI.H:643
Tuple2< Point, scalar > circumCircle() const
Return the circum centre and radius.
Definition: triangleI.H:120
Point randomPoint(randomGenerator &rndGen) const
Return a random point on the triangle from a uniform.
Definition: triangleI.H:215
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:407
bool classify(const point &p, label &nearType, label &nearLabel) const
Classify nearest point to p in triangle plane.
Definition: triangleI.H:657
vector normal() const
Return unit normal.
Definition: triangleI.H:110
const Point & b() const
Return second vertex.
Definition: triangleI.H:76
scalar quality() const
Return quality: Ratio of triangle and circum-circle.
Definition: triangleI.H:151
const pointField & points
volScalarField & b
Definition: createFields.H:27
const dimensionedScalar c1
First radiation constant: default SI units: [W/m^2].
const dimensionedScalar c2
Second radiation constant: default SI units: [m K].
const dimensionedScalar c
Speed of light in a vacuum.
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
static const zero Zero
Definition: zero.H:97
scalar minDist(const List< pointIndexHit > &hitList)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
line< point, const point & > linePointRef
Line using referred points.
Definition: linePointRef.H:45
static const Identity< scalar > I
Definition: Identity.H:93
Barycentric2D< scalar > barycentric2D
A scalar version of the templated Barycentric2D.
Definition: barycentric2D.H:45
void det(LagrangianPatchField< scalar > &f, const LagrangianPatchField< tensor > &f1)
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void symm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
PointHit< point > pointHit
Definition: pointHit.H:41
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy, recursively if necessary, the source to the destination.
Definition: POSIX.C:753
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
static const char nl
Definition: Ostream.H:267
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: POSIX.C:908
barycentric2D barycentric2D01(randomGenerator &rndGen)
Generate a random barycentric coordinate within the unit triangle.
Definition: barycentric2D.C:31
randomGenerator rndGen(653213)
volScalarField & p