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