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