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 | Copyright (C) 2011-2017 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>
96 inline Foam::scalar Foam::triangle<Point, PointRef>::mag() const
97 {
98  return Foam::mag(normal());
99 }
100 
101 
102 template<class Point, class PointRef>
104 {
105  return 0.5*((b_ - a_)^(c_ - a_));
106 }
107 
108 
109 template<class Point, class PointRef>
111 {
112  scalar d1 = (c_ - a_) & (b_ - a_);
113  scalar d2 = -(c_ - b_) & (b_ - a_);
114  scalar d3 = (c_ - a_) & (c_ - b_);
115 
116  scalar c1 = d2*d3;
117  scalar c2 = d3*d1;
118  scalar c3 = d1*d2;
119 
120  scalar c = c1 + c2 + c3;
121 
122  if (Foam::mag(c) < ROOTVSMALL)
123  {
124  // Degenerate triangle, returning centre instead of circumCentre.
125 
126  return centre();
127  }
128 
129  return
130  (
131  ((c2 + c3)*a_ + (c3 + c1)*b_ + (c1 + c2)*c_)/(2*c)
132  );
133 }
134 
135 
136 template<class Point, class PointRef>
138 {
139  scalar d1 = (c_ - a_) & (b_ - a_);
140  scalar d2 = -(c_ - b_) & (b_ - a_);
141  scalar d3 = (c_ - a_) & (c_ - b_);
142 
143  scalar denom = d2*d3 + d3*d1 + d1*d2;
144 
145  if (Foam::mag(denom) < VSMALL)
146  {
147  // Degenerate triangle, returning GREAT for circumRadius.
148 
149  return GREAT;
150  }
151  else
152  {
153  scalar a = (d1 + d2)*(d2 + d3)*(d3 + d1) / denom;
154 
155  return 0.5*Foam::sqrt(min(GREAT, max(0, a)));
156  }
157 }
158 
159 
160 template<class Point, class PointRef>
161 inline Foam::scalar Foam::triangle<Point, PointRef>::quality() const
162 {
163  scalar c = circumRadius();
164 
165  if (c < ROOTVSMALL)
166  {
167  // zero circumRadius, something has gone wrong.
168  return SMALL;
169  }
170 
171  return mag()/(Foam::sqr(c)*3.0*sqrt(3.0)/4.0);
172 }
173 
174 
175 template<class Point, class PointRef>
177 (
178  const triangle& t
179 ) const
180 {
181  return (1.0/12.0)*
182  (
183  ((t.a_ - a_) & ((b_ - a_)^(c_ - a_)))
184  + ((t.b_ - b_) & ((c_ - b_)^(t.a_ - b_)))
185  + ((c_ - t.c_) & ((t.b_ - t.c_)^(t.a_ - t.c_)))
186 
187  + ((t.a_ - a_) & ((b_ - a_)^(c_ - a_)))
188  + ((b_ - t.b_) & ((t.a_ - t.b_)^(t.c_ - t.b_)))
189  + ((c_ - t.c_) & ((b_ - t.c_)^(t.a_ - t.c_)))
190  );
191 }
192 
193 
194 template<class Point, class PointRef>
196 (
197  PointRef refPt,
198  scalar density
199 ) const
200 {
201  Point aRel = a_ - refPt;
202  Point bRel = b_ - refPt;
203  Point cRel = c_ - refPt;
204 
205  tensor V
206  (
207  aRel.x(), aRel.y(), aRel.z(),
208  bRel.x(), bRel.y(), bRel.z(),
209  cRel.x(), cRel.y(), cRel.z()
210  );
211 
212  scalar a = Foam::mag((b_ - a_)^(c_ - a_));
213 
214  tensor S = 1/24.0*(tensor::one + I);
215 
216  return
217  (
218  a*I/24.0*
219  (
220  (aRel & aRel)
221  + (bRel & bRel)
222  + (cRel & cRel)
223  + ((aRel + bRel + cRel) & (aRel + bRel + cRel))
224  )
225  - a*(V.T() & S & V)
226  )
227  *density;
228 }
229 
230 
231 template<class Point, class PointRef>
233 {
234  return barycentricToPoint(barycentric2D01(rndGen));
235 }
236 
237 
238 template<class Point, class PointRef>
240 (
241  cachedRandom& rndGen
242 ) const
243 {
244  return barycentricToPoint(barycentric2D01(rndGen));
245 }
246 
247 
248 template<class Point, class PointRef>
250 (
251  const barycentric2D& bary
252 ) const
253 {
254  return bary[0]*a_ + bary[1]*b_ + bary[2]*c_;
255 }
256 
257 
258 template<class Point, class PointRef>
260 (
261  const point& pt
262 ) const
263 {
264  barycentric2D bary;
265  pointToBarycentric(pt, bary);
266  return bary;
267 }
268 
269 
270 template<class Point, class PointRef>
272 (
273  const point& pt,
274  barycentric2D& bary
275 ) const
276 {
277  // Reference:
278  // Real-time collision detection, Christer Ericson, 2005, p47-48
279 
280  vector v0 = b_ - a_;
281  vector v1 = c_ - a_;
282  vector v2 = pt - a_;
283 
284  scalar d00 = v0 & v0;
285  scalar d01 = v0 & v1;
286  scalar d11 = v1 & v1;
287  scalar d20 = v2 & v0;
288  scalar d21 = v2 & v1;
289 
290  scalar denom = d00*d11 - d01*d01;
291 
292  if (Foam::mag(denom) < SMALL)
293  {
294  // Degenerate triangle, returning 1/3 barycentric coordinates.
295 
296  bary = barycentric2D(1.0/3.0, 1.0/3.0, 1.0/3.0);
297 
298  return denom;
299  }
300 
301  bary[1] = (d11*d20 - d01*d21)/denom;
302  bary[2] = (d00*d21 - d01*d20)/denom;
303  bary[0] = 1.0 - bary[1] - bary[2];
304 
305  return denom;
306 }
307 
308 
309 template<class Point, class PointRef>
311 (
312  const point& p,
313  const vector& q,
314  const intersection::algorithm alg,
315  const intersection::direction dir
316 ) const
317 {
318  // Express triangle in terms of baseVertex (point a_) and
319  // two edge vectors
320  vector E0 = b_ - a_;
321  vector E1 = c_ - a_;
322 
323  // Initialize intersection to miss.
324  pointHit inter(p);
325 
326  vector n(0.5*(E0 ^ E1));
327  scalar area = Foam::mag(n);
328 
329  if (area < VSMALL)
330  {
331  // Ineligible miss.
332  inter.setMiss(false);
333 
334  // The miss point is the nearest point on the triangle. Take any one.
335  inter.setPoint(a_);
336 
337  // The distance to the miss is the distance between the
338  // original point and plane of intersection. No normal so use
339  // distance from p to a_
340  inter.setDistance(Foam::mag(a_ - p));
341 
342  return inter;
343  }
344 
345  vector q1 = q/Foam::mag(q);
346 
347  if (dir == intersection::CONTACT_SPHERE)
348  {
349  n /= area;
350 
351  return ray(p, q1 - n, alg, intersection::VECTOR);
352  }
353 
354  // Intersection point with triangle plane
355  point pInter;
356  // Is intersection point inside triangle
357  bool hit;
358  {
359  // Reuse the fast ray intersection routine below in FULL_RAY
360  // mode since the original intersection routine has rounding problems.
361  pointHit fastInter = intersection(p, q1, intersection::FULL_RAY);
362  hit = fastInter.hit();
363 
364  if (hit)
365  {
366  pInter = fastInter.rawPoint();
367  }
368  else
369  {
370  // Calculate intersection of ray with triangle plane
371  vector v = a_ - p;
372  pInter = p + (q1&v)*q1;
373  }
374  }
375 
376  // Distance to intersection point
377  scalar dist = q1 & (pInter - p);
378 
379  const scalar planarPointTol =
380  Foam::min
381  (
382  Foam::min
383  (
384  Foam::mag(E0),
385  Foam::mag(E1)
386  ),
387  Foam::mag(c_ - b_)
388  )*intersection::planarTol();
389 
390  bool eligible =
391  alg == intersection::FULL_RAY
392  || (alg == intersection::HALF_RAY && dist > -planarPointTol)
393  || (
394  alg == intersection::VISIBLE
395  && ((q1 & normal()) < -VSMALL)
396  );
397 
398  if (hit && eligible)
399  {
400  // Hit. Set distance to intersection.
401  inter.setHit();
402  inter.setPoint(pInter);
403  inter.setDistance(dist);
404  }
405  else
406  {
407  // Miss or ineligible hit.
408  inter.setMiss(eligible);
409 
410  // The miss point is the nearest point on the triangle
411  inter.setPoint(nearestPoint(p).rawPoint());
412 
413  // The distance to the miss is the distance between the
414  // original point and plane of intersection
415  inter.setDistance(Foam::mag(pInter - p));
416  }
417 
418 
419  return inter;
420 }
421 
422 
423 // From "Fast, Minimum Storage Ray/Triangle Intersection"
424 // Moeller/Trumbore.
425 template<class Point, class PointRef>
427 (
428  const point& orig,
429  const vector& dir,
430  const intersection::algorithm alg,
431  const scalar tol
432 ) const
433 {
434  const vector edge1 = b_ - a_;
435  const vector edge2 = c_ - a_;
436 
437  // begin calculating determinant - also used to calculate U parameter
438  const vector pVec = dir ^ edge2;
439 
440  // if determinant is near zero, ray lies in plane of triangle
441  const scalar det = edge1 & pVec;
442 
443  // Initialise to miss
444  pointHit intersection(false, Zero, GREAT, false);
445 
446  if (alg == intersection::VISIBLE)
447  {
448  // Culling branch
449  if (det < ROOTVSMALL)
450  {
451  // Ray on wrong side of triangle. Return miss
452  return intersection;
453  }
454  }
455  else if (alg == intersection::HALF_RAY || alg == intersection::FULL_RAY)
456  {
457  // Non-culling branch
458  if (det > -ROOTVSMALL && det < ROOTVSMALL)
459  {
460  // Ray parallel to triangle. Return miss
461  return intersection;
462  }
463  }
464 
465  const scalar inv_det = 1.0 / det;
466 
467  /* calculate distance from a_ to ray origin */
468  const vector tVec = orig-a_;
469 
470  /* calculate U parameter and test bounds */
471  const scalar u = (tVec & pVec)*inv_det;
472 
473  if (u < -tol || u > 1.0+tol)
474  {
475  // return miss
476  return intersection;
477  }
478 
479  /* prepare to test V parameter */
480  const vector qVec = tVec ^ edge1;
481 
482  /* calculate V parameter and test bounds */
483  const scalar v = (dir & qVec) * inv_det;
484 
485  if (v < -tol || u + v > 1.0+tol)
486  {
487  // return miss
488  return intersection;
489  }
490 
491  /* calculate t, scale parameters, ray intersects triangle */
492  const scalar t = (edge2 & qVec) * inv_det;
493 
494  if (alg == intersection::HALF_RAY && t < -tol)
495  {
496  // Wrong side of orig. Return miss
497  return intersection;
498  }
499 
500  intersection.setHit();
501  intersection.setPoint(a_ + u*edge1 + v*edge2);
502  intersection.setDistance(t);
503 
504  return intersection;
505 }
506 
507 
508 template<class Point, class PointRef>
510 (
511  const point& p,
512  label& nearType,
513  label& nearLabel
514 ) const
515 {
516  // Adapted from:
517  // Real-time collision detection, Christer Ericson, 2005, p136-142
518 
519  // Check if P in vertex region outside A
520  vector ab = b_ - a_;
521  vector ac = c_ - a_;
522  vector ap = p - a_;
523 
524  scalar d1 = ab & ap;
525  scalar d2 = ac & ap;
526 
527  if (d1 <= 0.0 && d2 <= 0.0)
528  {
529  // barycentric coordinates (1, 0, 0)
530 
531  nearType = POINT;
532  nearLabel = 0;
533  return pointHit(false, a_, Foam::mag(a_ - p), true);
534  }
535 
536  // Check if P in vertex region outside B
537  vector bp = p - b_;
538  scalar d3 = ab & bp;
539  scalar d4 = ac & bp;
540 
541  if (d3 >= 0.0 && d4 <= d3)
542  {
543  // barycentric coordinates (0, 1, 0)
544 
545  nearType = POINT;
546  nearLabel = 1;
547  return pointHit(false, b_, Foam::mag(b_ - p), true);
548  }
549 
550  // Check if P in edge region of AB, if so return projection of P onto AB
551  scalar vc = d1*d4 - d3*d2;
552 
553  if (vc <= 0.0 && d1 >= 0.0 && d3 <= 0.0)
554  {
555  if ((d1 - d3) < ROOTVSMALL)
556  {
557  // Degenerate triangle, for d1 = d3, a_ and b_ are likely coincident
558  nearType = POINT;
559  nearLabel = 0;
560  return pointHit(false, a_, Foam::mag(a_ - p), true);
561  }
562 
563  // barycentric coordinates (1-v, v, 0)
564  scalar v = d1/(d1 - d3);
565 
566  point nearPt = a_ + v*ab;
567  nearType = EDGE;
568  nearLabel = 0;
569  return pointHit(false, nearPt, Foam::mag(nearPt - p), true);
570  }
571 
572  // Check if P in vertex region outside C
573  vector cp = p - c_;
574  scalar d5 = ab & cp;
575  scalar d6 = ac & cp;
576 
577  if (d6 >= 0.0 && d5 <= d6)
578  {
579  // barycentric coordinates (0, 0, 1)
580 
581  nearType = POINT;
582  nearLabel = 2;
583  return pointHit(false, c_, Foam::mag(c_ - p), true);
584  }
585 
586  // Check if P in edge region of AC, if so return projection of P onto AC
587  scalar vb = d5*d2 - d1*d6;
588 
589  if (vb <= 0.0 && d2 >= 0.0 && d6 <= 0.0)
590  {
591  if ((d2 - d6) < ROOTVSMALL)
592  {
593  // Degenerate triangle, for d2 = d6, a_ and c_ are likely coincident
594  nearType = POINT;
595  nearLabel = 0;
596  return pointHit(false, a_, Foam::mag(a_ - p), true);
597  }
598 
599  // barycentric coordinates (1-w, 0, w)
600  scalar w = d2/(d2 - d6);
601 
602  point nearPt = a_ + w*ac;
603  nearType = EDGE;
604  nearLabel = 2;
605  return pointHit(false, nearPt, Foam::mag(nearPt - p), true);
606  }
607 
608  // Check if P in edge region of BC, if so return projection of P onto BC
609  scalar va = d3*d6 - d5*d4;
610 
611  if (va <= 0.0 && (d4 - d3) >= 0.0 && (d5 - d6) >= 0.0)
612  {
613  if (((d4 - d3) + (d5 - d6)) < ROOTVSMALL)
614  {
615  // Degenerate triangle, for (d4 - d3) = (d6 - d5), b_ and c_ are
616  // likely coincident
617  nearType = POINT;
618  nearLabel = 1;
619  return pointHit(false, b_, Foam::mag(b_ - p), true);
620  }
621 
622  // barycentric coordinates (0, 1-w, w)
623  scalar w = (d4 - d3)/((d4 - d3) + (d5 - d6));
624 
625  point nearPt = b_ + w*(c_ - b_);
626  nearType = EDGE;
627  nearLabel = 1;
628  return pointHit(false, nearPt, Foam::mag(nearPt - p), true);
629  }
630 
631  // P inside face region. Compute Q through its barycentric
632  // coordinates (u, v, w)
633 
634  if ((va + vb + vc) < ROOTVSMALL)
635  {
636  // Degenerate triangle, return the centre because no edge or points are
637  // closest
638  point nearPt = centre();
639  nearType = NONE,
640  nearLabel = -1;
641  return pointHit(true, nearPt, Foam::mag(nearPt - p), false);
642  }
643 
644  scalar denom = 1.0/(va + vb + vc);
645  scalar v = vb * denom;
646  scalar w = vc * denom;
647 
648  // = u*a + v*b + w*c, u = va*denom = 1.0 - v - w
649 
650  point nearPt = a_ + ab*v + ac*w;
651  nearType = NONE,
652  nearLabel = -1;
653  return pointHit(true, nearPt, Foam::mag(nearPt - p), false);
654 }
655 
656 
657 template<class Point, class PointRef>
659 (
660  const point& p
661 ) const
662 {
663  // Dummy labels
664  label nearType = -1;
665  label nearLabel = -1;
666 
667  return nearestPointClassify(p, nearType, nearLabel);
668 }
669 
670 
671 template<class Point, class PointRef>
673 (
674  const point& p,
675  label& nearType,
676  label& nearLabel
677 ) const
678 {
679  return nearestPointClassify(p, nearType, nearLabel).hit();
680 }
681 
682 
683 template<class Point, class PointRef>
685 (
686  const linePointRef& ln,
687  pointHit& lnInfo
688 ) const
689 {
690  vector q = ln.vec();
691  pointHit triInfo
692  (
694  (
695  ln.start(),
696  q,
697  intersection::FULL_RAY
698  )
699  );
700 
701  if (triInfo.hit())
702  {
703  // Line hits triangle. Find point on line.
704  if (triInfo.distance() > 1)
705  {
706  // Hit beyond endpoint
707  lnInfo.setMiss(true);
708  lnInfo.setPoint(ln.end());
709  scalar dist = Foam::mag(triInfo.hitPoint()-lnInfo.missPoint());
710  lnInfo.setDistance(dist);
711  triInfo.setMiss(true);
712  triInfo.setDistance(dist);
713  }
714  else if (triInfo.distance() < 0)
715  {
716  // Hit beyond startpoint
717  lnInfo.setMiss(true);
718  lnInfo.setPoint(ln.start());
719  scalar dist = Foam::mag(triInfo.hitPoint()-lnInfo.missPoint());
720  lnInfo.setDistance(dist);
721  triInfo.setMiss(true);
722  triInfo.setDistance(dist);
723  }
724  else
725  {
726  // Hit on line
727  lnInfo.setHit();
728  lnInfo.setPoint(triInfo.hitPoint());
729  lnInfo.setDistance(0.0);
730  triInfo.setDistance(0.0);
731  }
732  }
733  else
734  {
735  // Line skips triangle. See which triangle edge it gets closest to
736 
737  point nearestEdgePoint;
738  point nearestLinePoint;
739  //label minEdgeIndex = 0;
740  scalar minDist = ln.nearestDist
741  (
742  linePointRef(a_, b_),
743  nearestLinePoint,
744  nearestEdgePoint
745  );
746 
747  {
748  point linePoint;
749  point triEdgePoint;
750  scalar dist = ln.nearestDist
751  (
752  linePointRef(b_, c_),
753  linePoint,
754  triEdgePoint
755  );
756  if (dist < minDist)
757  {
758  minDist = dist;
759  nearestEdgePoint = triEdgePoint;
760  nearestLinePoint = linePoint;
761  //minEdgeIndex = 1;
762  }
763  }
764 
765  {
766  point linePoint;
767  point triEdgePoint;
768  scalar dist = ln.nearestDist
769  (
770  linePointRef(c_, a_),
771  linePoint,
772  triEdgePoint
773  );
774  if (dist < minDist)
775  {
776  minDist = dist;
777  nearestEdgePoint = triEdgePoint;
778  nearestLinePoint = linePoint;
779  //minEdgeIndex = 2;
780  }
781  }
782 
783  lnInfo.setDistance(minDist);
784  triInfo.setDistance(minDist);
785  triInfo.setMiss(false);
786  triInfo.setPoint(nearestEdgePoint);
787 
788  // Convert point on line to pointHit
789  if (Foam::mag(nearestLinePoint-ln.start()) < SMALL)
790  {
791  lnInfo.setMiss(true);
792  lnInfo.setPoint(ln.start());
793  }
794  else if (Foam::mag(nearestLinePoint-ln.end()) < SMALL)
795  {
796  lnInfo.setMiss(true);
797  lnInfo.setPoint(ln.end());
798  }
799  else
800  {
801  lnInfo.setHit();
802  lnInfo.setPoint(nearestLinePoint);
803  }
804  }
805  return triInfo;
806 }
807 
808 
809 // * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
810 
811 template<class Point, class PointRef>
812 inline Foam::Istream& Foam::operator>>
813 (
814  Istream& is,
816 )
817 {
818  is.readBegin("triangle");
819  is >> t.a_ >> t.b_ >> t.c_;
820  is.readEnd("triangle");
821 
822  is.check("Istream& operator>>(Istream&, triangle&)");
823  return is;
824 }
825 
826 
827 template<class Point, class PointRef>
828 inline Foam::Ostream& Foam::operator<<
829 (
830  Ostream& os,
832 )
833 {
834  os << nl
835  << token::BEGIN_LIST
836  << t.a_ << token::SPACE
837  << t.b_ << token::SPACE
838  << t.c_
839  << token::END_LIST;
840 
841  return os;
842 }
843 
844 
845 // ************************************************************************* //
scalar sweptVol(const triangle &t) const
Return swept-volume.
Definition: triangleI.H:177
Istream & readBegin(const char *funcName)
Definition: Istream.C:86
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:58
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
A line primitive.
Definition: line.H:56
const Point & c() const
Return third vertex.
Definition: triangleI.H:82
Barycentric2D< scalar > barycentric2D
A scalar version of the templated Barycentric2D.
Definition: barycentric2D.H:46
barycentric2D pointToBarycentric(const point &pt) const
Calculate the barycentric coordinates from the given point.
Definition: triangleI.H:260
scalar mag() const
Return scalar magnitude.
Definition: triangleI.H:96
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
scalar quality() const
Return quality: Ratio of triangle and circum-circle.
Definition: triangleI.H:161
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
const dimensionedScalar c2
Second radiation constant: default SI units: [m.K].
dimensionedScalar sqrt(const dimensionedScalar &ds)
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:510
Random number generator.
Definition: cachedRandom.H:63
dimensionedScalar det(const dimensionedSphericalTensor &dt)
const Point & missPoint() const
Return miss point.
Definition: PointHit.H:145
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
bool classify(const point &p, label &nearType, label &nearLabel) const
Classify nearest point to p in triangle plane.
Definition: triangleI.H:673
pointHit nearestPoint(const point &p) const
Return nearest point to p on triangle.
Definition: triangleI.H:659
void setHit()
Definition: PointHit.H:169
Foam::intersection.
Definition: intersection.H:49
void setDistance(const scalar d)
Definition: PointHit.H:186
line< point, const point & > linePointRef
Line using referred points.
Definition: linePointRef.H:45
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:311
static const Identity< scalar > I
Definition: Identity.H:93
const Point & hitPoint() const
Return hit point.
Definition: PointHit.H:126
PointRef start() const
Return first vertex.
Definition: lineI.H:60
Point centre() const
Return centre (centroid)
Definition: triangleI.H:89
bool hit() const
Is there a hit.
Definition: PointHit.H:120
static const zero Zero
Definition: zero.H:91
scalar circumRadius() const
Return circum-radius.
Definition: triangleI.H:137
tensor inertia(PointRef refPt=vector::zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
Definition: triangleI.H:196
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:61
Simple random number generator.
Definition: Random.H:49
vector normal() const
Return vector normal.
Definition: triangleI.H:103
void setMiss(const bool eligible)
Definition: PointHit.H:175
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
triangle(const Point &a, const Point &b, const Point &c)
Construct from three points.
Definition: triangleI.H:34
const Point & rawPoint() const
Return point with no checking.
Definition: PointHit.H:158
static const char nl
Definition: Ostream.H:262
const volScalarField & cp
Point barycentricToPoint(const barycentric2D &bary) const
Calculate the point from the given barycentric coordinates.
Definition: triangleI.H:250
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::barycentric2D barycentric2D01(Foam::scalar s, Foam::scalar t)
Definition: barycentric2D.C:33
const Point & a() const
Return first vertex.
Definition: triangleI.H:70
A normal distribution model.
Point vec() const
Return start-end vector.
Definition: lineI.H:87
Point circumCentre() const
Return circum-centre.
Definition: triangleI.H:110
const Point & b() const
Return second vertex.
Definition: triangleI.H:76
PointRef end() const
Return second vertex.
Definition: lineI.H:66
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
Point randomPoint(Random &rndGen) const
Return a random point on the triangle from a uniform.
Definition: triangleI.H:232
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
PointHit< point > pointHit
Definition: pointHit.H:41
scalar distance() const
Return distance to hit.
Definition: PointHit.H:139
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:427
void setPoint(const Point &p)
Definition: PointHit.H:181
PointHit< Point > nearestDist(const Point &p) const
Return nearest distance to line from a given point.
Definition: lineI.H:95