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