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-2024 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  Point aRel = a_ - refPt;
186  Point bRel = b_ - refPt;
187  Point cRel = c_ - refPt;
188 
189  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  scalar a = Foam::mag((b_ - a_)^(c_ - a_));
197 
198  tensor S = 1/24.0*(tensor::one + I);
199 
200  return
201  (
202  a*I/24.0*
203  (
204  (aRel & aRel)
205  + (bRel & bRel)
206  + (cRel & cRel)
207  + ((aRel + bRel + cRel) & (aRel + bRel + cRel))
208  )
209  - a*(V.T() & S & V)
210  )
211  *density;
212 }
213 
214 
215 template<class Point, class PointRef>
217 (
219 ) const
220 {
221  return barycentricToPoint(barycentric2D01(rndGen));
222 }
223 
224 
225 template<class Point, class PointRef>
227 (
228  const barycentric2D& bary
229 ) const
230 {
231  return bary[0]*a_ + bary[1]*b_ + bary[2]*c_;
232 }
233 
234 
235 template<class Point, class PointRef>
237 (
238  const point& pt
239 ) const
240 {
241  barycentric2D bary;
242  pointToBarycentric(pt, bary);
243  return bary;
244 }
245 
246 
247 template<class Point, class PointRef>
249 (
250  const point& pt,
251  barycentric2D& bary
252 ) const
253 {
254  // Reference:
255  // Real-time collision detection, Christer Ericson, 2005, p47-48
256 
257  vector v0 = b_ - a_;
258  vector v1 = c_ - a_;
259  vector v2 = pt - a_;
260 
261  scalar d00 = v0 & v0;
262  scalar d01 = v0 & v1;
263  scalar d11 = v1 & v1;
264  scalar d20 = v2 & v0;
265  scalar d21 = v2 & v1;
266 
267  scalar denom = d00*d11 - d01*d01;
268 
269  if (Foam::mag(denom) < small)
270  {
271  // Degenerate triangle, returning 1/3 barycentric coordinates.
272 
273  bary = barycentric2D(1.0/3.0, 1.0/3.0, 1.0/3.0);
274 
275  return denom;
276  }
277 
278  bary[1] = (d11*d20 - d01*d21)/denom;
279  bary[2] = (d00*d21 - d01*d20)/denom;
280  bary[0] = 1.0 - bary[1] - bary[2];
281 
282  return denom;
283 }
284 
285 
286 template<class Point, class PointRef>
288 (
289  const point& p,
290  const vector& q,
291  const intersection::algorithm alg,
292  const intersection::direction dir
293 ) const
294 {
295  // Express triangle in terms of baseVertex (point a_) and
296  // two edge vectors
297  const vector E0 = b_ - a_;
298  const vector E1 = c_ - a_;
299 
300  // Initialise intersection to miss.
301  pointHit inter(p);
302 
303  vector n(0.5*(E0 ^ E1));
304  const scalar magArea = Foam::mag(n);
305 
306  if (magArea < vSmall)
307  {
308  // Ineligible miss.
309  inter.setMiss(false);
310 
311  // The miss point is the nearest point on the triangle. Take any one.
312  inter.setPoint(a_);
313 
314  // The distance to the miss is the distance between the
315  // original point and plane of intersection. No normal so use
316  // distance from p to a_
317  inter.setDistance(Foam::mag(a_ - p));
318 
319  return inter;
320  }
321 
322  const vector q1 = q/Foam::mag(q);
323 
325  {
326  n /= magArea;
327 
328  return ray(p, q1 - n, alg, intersection::direction::vector);
329  }
330 
331  // Intersection point with triangle plane
332  point pInter;
333 
334  // Is intersection point inside triangle
335  bool hit;
336  {
337  // Reuse the fast ray intersection routine below in algorithm::fullRay
338  // mode since the original intersection routine has rounding problems.
339  pointHit fastInter = intersection
340  (
341  p,
342  q1,
344  );
345  hit = fastInter.hit();
346 
347  if (hit)
348  {
349  pInter = fastInter.rawPoint();
350  }
351  else
352  {
353  // Calculate intersection of ray with triangle plane
354  vector v = a_ - p;
355  pInter = p + (q1&v)*q1;
356  }
357  }
358 
359  // Distance to intersection point
360  const scalar dist = q1 & (pInter - p);
361 
362  const scalar planarPointTol =
363  Foam::min
364  (
365  Foam::min
366  (
367  Foam::mag(E0),
368  Foam::mag(E1)
369  ),
370  Foam::mag(c_ - b_)
372 
373  bool eligible =
375  || (alg == intersection::algorithm::halfRay && dist > -planarPointTol)
376  || (
378  && ((q1 & area()) < -vSmall)
379  );
380 
381  if (hit && eligible)
382  {
383  // Hit. Set distance to intersection.
384  inter.setHit();
385  inter.setPoint(pInter);
386  inter.setDistance(dist);
387  }
388  else
389  {
390  // Miss or ineligible hit.
391  inter.setMiss(eligible);
392 
393  // The miss point is the nearest point on the triangle
394  inter.setPoint(nearestPoint(p).rawPoint());
395 
396  // The distance to the miss is the distance between the
397  // original point and plane of intersection
398  inter.setDistance(Foam::mag(pInter - p));
399  }
400 
401  return inter;
402 }
403 
404 
405 // From "Fast, Minimum Storage Ray/Triangle Intersection"
406 // Moeller/Trumbore.
407 template<class Point, class PointRef>
409 (
410  const point& orig,
411  const vector& dir,
412  const intersection::algorithm alg,
413  const scalar tol
414 ) const
415 {
416  const vector edge1 = b_ - a_;
417  const vector edge2 = c_ - a_;
418 
419  // Begin calculating determinant - also used to calculate U parameter
420  const vector pVec = dir ^ edge2;
421 
422  // Ff determinant is near zero, ray lies in plane of triangle
423  const scalar det = edge1 & pVec;
424 
425  // Initialise to miss
426  pointHit intersection(false, Zero, great, false);
427 
429  {
430  // Culling branch
431  if (det < rootVSmall)
432  {
433  // Ray on wrong side of triangle. Return miss
434  return intersection;
435  }
436  }
437  else if
438  (
441  )
442  {
443  // Non-culling branch
444  if (det > -rootVSmall && det < rootVSmall)
445  {
446  // Ray parallel to triangle. Return miss
447  return intersection;
448  }
449  }
450 
451  const scalar inv_det = 1.0 / det;
452 
453  // Calculate distance from a_ to ray origin
454  const vector tVec = orig-a_;
455 
456  // Calculate U parameter and test bounds
457  const scalar u = (tVec & pVec)*inv_det;
458 
459  if (u < -tol || u > 1.0+tol)
460  {
461  // return miss
462  return intersection;
463  }
464 
465  // Prepare to test V parameter
466  const vector qVec = tVec ^ edge1;
467 
468  // Calculate V parameter and test bounds
469  const scalar v = (dir & qVec) * inv_det;
470 
471  if (v < -tol || u + v > 1.0+tol)
472  {
473  // return miss
474  return intersection;
475  }
476 
477  // Calculate t, scale parameters, ray intersects triangle
478  const scalar t = (edge2 & qVec) * inv_det;
479 
480  if (alg == intersection::algorithm::halfRay && t < -tol)
481  {
482  // Wrong side of orig. Return miss
483  return intersection;
484  }
485 
486  intersection.setHit();
487  intersection.setPoint(a_ + u*edge1 + v*edge2);
488  intersection.setDistance(t);
489 
490  return intersection;
491 }
492 
493 
494 template<class Point, class PointRef>
496 (
497  const point& p,
498  label& nearType,
499  label& nearLabel
500 ) const
501 {
502  // Adapted from:
503  // Real-time collision detection, Christer Ericson, 2005, p136-142
504 
505  // Check if P in vertex region outside A
506  vector ab = b_ - a_;
507  vector ac = c_ - a_;
508  vector ap = p - a_;
509 
510  scalar d1 = ab & ap;
511  scalar d2 = ac & ap;
512 
513  if (d1 <= 0.0 && d2 <= 0.0)
514  {
515  // barycentric coordinates (1, 0, 0)
516 
517  nearType = POINT;
518  nearLabel = 0;
519  return pointHit(false, a_, Foam::mag(a_ - p), true);
520  }
521 
522  // Check if P in vertex region outside B
523  vector bp = p - b_;
524  scalar d3 = ab & bp;
525  scalar d4 = ac & bp;
526 
527  if (d3 >= 0.0 && d4 <= d3)
528  {
529  // barycentric coordinates (0, 1, 0)
530 
531  nearType = POINT;
532  nearLabel = 1;
533  return pointHit(false, b_, Foam::mag(b_ - p), true);
534  }
535 
536  // Check if P in edge region of AB, if so return projection of P onto AB
537  scalar vc = d1*d4 - d3*d2;
538 
539  if (vc <= 0.0 && d1 >= 0.0 && d3 <= 0.0)
540  {
541  if ((d1 - d3) < rootVSmall)
542  {
543  // Degenerate triangle, for d1 = d3, a_ and b_ are likely coincident
544  nearType = POINT;
545  nearLabel = 0;
546  return pointHit(false, a_, Foam::mag(a_ - p), true);
547  }
548 
549  // barycentric coordinates (1-v, v, 0)
550  scalar v = d1/(d1 - d3);
551 
552  point nearPt = a_ + v*ab;
553  nearType = EDGE;
554  nearLabel = 0;
555  return pointHit(false, nearPt, Foam::mag(nearPt - p), true);
556  }
557 
558  // Check if P in vertex region outside C
559  vector cp = p - c_;
560  scalar d5 = ab & cp;
561  scalar d6 = ac & cp;
562 
563  if (d6 >= 0.0 && d5 <= d6)
564  {
565  // barycentric coordinates (0, 0, 1)
566 
567  nearType = POINT;
568  nearLabel = 2;
569  return pointHit(false, c_, Foam::mag(c_ - p), true);
570  }
571 
572  // Check if P in edge region of AC, if so return projection of P onto AC
573  scalar vb = d5*d2 - d1*d6;
574 
575  if (vb <= 0.0 && d2 >= 0.0 && d6 <= 0.0)
576  {
577  if ((d2 - d6) < rootVSmall)
578  {
579  // Degenerate triangle, for d2 = d6, a_ and c_ are likely coincident
580  nearType = POINT;
581  nearLabel = 0;
582  return pointHit(false, a_, Foam::mag(a_ - p), true);
583  }
584 
585  // barycentric coordinates (1-w, 0, w)
586  scalar w = d2/(d2 - d6);
587 
588  point nearPt = a_ + w*ac;
589  nearType = EDGE;
590  nearLabel = 2;
591  return pointHit(false, nearPt, Foam::mag(nearPt - p), true);
592  }
593 
594  // Check if P in edge region of BC, if so return projection of P onto BC
595  scalar va = d3*d6 - d5*d4;
596 
597  if (va <= 0.0 && (d4 - d3) >= 0.0 && (d5 - d6) >= 0.0)
598  {
599  if (((d4 - d3) + (d5 - d6)) < rootVSmall)
600  {
601  // Degenerate triangle, for (d4 - d3) = (d6 - d5), b_ and c_ are
602  // likely coincident
603  nearType = POINT;
604  nearLabel = 1;
605  return pointHit(false, b_, Foam::mag(b_ - p), true);
606  }
607 
608  // barycentric coordinates (0, 1-w, w)
609  scalar w = (d4 - d3)/((d4 - d3) + (d5 - d6));
610 
611  point nearPt = b_ + w*(c_ - b_);
612  nearType = EDGE;
613  nearLabel = 1;
614  return pointHit(false, nearPt, Foam::mag(nearPt - p), true);
615  }
616 
617  // P inside face region. Compute Q through its barycentric
618  // coordinates (u, v, w)
619 
620  if ((va + vb + vc) < rootVSmall)
621  {
622  // Degenerate triangle, return the centre because no edge or points are
623  // closest
624  point nearPt = centre();
625  nearType = NONE,
626  nearLabel = -1;
627  return pointHit(true, nearPt, Foam::mag(nearPt - p), false);
628  }
629 
630  scalar denom = 1.0/(va + vb + vc);
631  scalar v = vb * denom;
632  scalar w = vc * denom;
633 
634  // = u*a + v*b + w*c, u = va*denom = 1.0 - v - w
635 
636  point nearPt = a_ + ab*v + ac*w;
637  nearType = NONE,
638  nearLabel = -1;
639  return pointHit(true, nearPt, Foam::mag(nearPt - p), false);
640 }
641 
642 
643 template<class Point, class PointRef>
645 (
646  const point& p
647 ) const
648 {
649  // Dummy labels
650  label nearType = -1;
651  label nearLabel = -1;
652 
653  return nearestPointClassify(p, nearType, nearLabel);
654 }
655 
656 
657 template<class Point, class PointRef>
659 (
660  const point& p,
661  label& nearType,
662  label& nearLabel
663 ) const
664 {
665  return nearestPointClassify(p, nearType, nearLabel).hit();
666 }
667 
668 
669 template<class Point, class PointRef>
671 (
672  const linePointRef& ln,
673  pointHit& lnInfo
674 ) const
675 {
676  vector q = ln.vec();
677  pointHit triInfo
678  (
680  (
681  ln.start(),
682  q,
684  )
685  );
686 
687  if (triInfo.hit())
688  {
689  // Line hits triangle. Find point on line.
690  if (triInfo.distance() > 1)
691  {
692  // Hit beyond endpoint
693  lnInfo.setMiss(true);
694  lnInfo.setPoint(ln.end());
695  scalar dist = Foam::mag(triInfo.hitPoint()-lnInfo.missPoint());
696  lnInfo.setDistance(dist);
697  triInfo.setMiss(true);
698  triInfo.setDistance(dist);
699  }
700  else if (triInfo.distance() < 0)
701  {
702  // Hit beyond startpoint
703  lnInfo.setMiss(true);
704  lnInfo.setPoint(ln.start());
705  scalar dist = Foam::mag(triInfo.hitPoint()-lnInfo.missPoint());
706  lnInfo.setDistance(dist);
707  triInfo.setMiss(true);
708  triInfo.setDistance(dist);
709  }
710  else
711  {
712  // Hit on line
713  lnInfo.setHit();
714  lnInfo.setPoint(triInfo.hitPoint());
715  lnInfo.setDistance(0.0);
716  triInfo.setDistance(0.0);
717  }
718  }
719  else
720  {
721  // Line skips triangle. See which triangle edge it gets closest to
722 
723  point nearestEdgePoint;
724  point nearestLinePoint;
725  // label minEdgeIndex = 0;
726  scalar minDist = ln.nearestDist
727  (
728  linePointRef(a_, b_),
729  nearestLinePoint,
730  nearestEdgePoint
731  );
732 
733  {
734  point linePoint;
735  point triEdgePoint;
736  scalar dist = ln.nearestDist
737  (
738  linePointRef(b_, c_),
739  linePoint,
740  triEdgePoint
741  );
742  if (dist < minDist)
743  {
744  minDist = dist;
745  nearestEdgePoint = triEdgePoint;
746  nearestLinePoint = linePoint;
747  // minEdgeIndex = 1;
748  }
749  }
750 
751  {
752  point linePoint;
753  point triEdgePoint;
754  scalar dist = ln.nearestDist
755  (
756  linePointRef(c_, a_),
757  linePoint,
758  triEdgePoint
759  );
760  if (dist < minDist)
761  {
762  minDist = dist;
763  nearestEdgePoint = triEdgePoint;
764  nearestLinePoint = linePoint;
765  // minEdgeIndex = 2;
766  }
767  }
768 
769  lnInfo.setDistance(minDist);
770  triInfo.setDistance(minDist);
771  triInfo.setMiss(false);
772  triInfo.setPoint(nearestEdgePoint);
773 
774  // Convert point on line to pointHit
775  if (Foam::mag(nearestLinePoint-ln.start()) < small)
776  {
777  lnInfo.setMiss(true);
778  lnInfo.setPoint(ln.start());
779  }
780  else if (Foam::mag(nearestLinePoint-ln.end()) < small)
781  {
782  lnInfo.setMiss(true);
783  lnInfo.setPoint(ln.end());
784  }
785  else
786  {
787  lnInfo.setHit();
788  lnInfo.setPoint(nearestLinePoint);
789  }
790  }
791  return triInfo;
792 }
793 
794 
795 // * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
796 
797 template<class Point, class PointRef>
798 inline Foam::Istream& Foam::operator>>
799 (
800  Istream& is,
802 )
803 {
804  is.readBegin("triangle");
805  is >> t.a_ >> t.b_ >> t.c_;
806  is.readEnd("triangle");
807 
808  is.check("Istream& operator>>(Istream&, triangle&)");
809  return is;
810 }
811 
812 
813 template<class Point, class PointRef>
814 inline Foam::Ostream& Foam::operator<<
815 (
816  Ostream& os,
817  const triangle<Point, PointRef>& t
818 )
819 {
820  os << nl
822  << t.a_ << token::SPACE
823  << t.b_ << token::SPACE
824  << t.c_
825  << token::END_LIST;
826 
827  return os;
828 }
829 
830 
831 // ************************************************************************* //
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 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:237
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:227
tensor inertia(PointRef refPt=vector::zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
Definition: triangleI.H:180
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:496
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:288
const Point & c() const
Return third vertex.
Definition: triangleI.H:82
scalar mag() const
Return scalar magnitude.
Definition: triangleI.H:103
pointHit nearestPoint(const point &p) const
Return nearest point to p on triangle.
Definition: triangleI.H:645
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:217
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:409
bool classify(const point &p, label &nearType, label &nearLabel) const
Classify nearest point to p in triangle plane.
Definition: triangleI.H:659
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:25
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 > &)
dimensionedScalar det(const dimensionedSphericalTensor &dt)
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
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
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
PointHit< point > pointHit
Definition: pointHit.H:41
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
static const char nl
Definition: Ostream.H:266
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