86template<
class Po
int,
class Po
intRef>
100template<
class Po
int,
class Po
intRef>
106 a_(
pts.template get<0>()),
107 b_(
pts.template get<1>()),
108 c_(
pts.template get<2>())
112template<
class Po
int,
class Po
intRef>
120 b_(
points[indices.template get<1>()]),
121 c_(
points[indices.template get<2>()])
125template<
class Po
int,
class Po
intRef>
140template<
class Po
int,
class Po
intRef>
149template<
class Po
int,
class Po
intRef>
157 return (1.0/3.0)*(
p0 + p1 + p2);
161template<
class Po
int,
class Po
intRef>
164 return (1.0/3.0)*(a_ + b_ + c_);
174template<
class Po
int,
class Po
intRef>
182 return 0.5*((p1 -
p0)^(p2 -
p0));
186template<
class Po
int,
class Po
intRef>
189 return 0.5*((b_ - a_)^(c_ - a_));
199template<
class Po
int,
class Po
intRef>
208 (void)
n.normalise(ROOTVSMALL);
213template<
class Po
int,
class Po
intRef>
226template<
class Po
int,
class Po
intRef>
242template<
class Po
int,
class Po
intRef>
255template<
class Po
int,
class Po
intRef>
258 return Point(c_ - b_);
261template<
class Po
int,
class Po
intRef>
264 return Point(a_ - c_);
267template<
class Po
int,
class Po
intRef>
270 return Point(b_ - a_);
315template<
class Po
int,
class Po
intRef>
328template<
class Po
int,
class Po
intRef>
335template<
class Po
int,
class Po
intRef>
338 scalar d1 = (c_ - a_) & (b_ - a_);
339 scalar d2 = -(c_ - b_) & (b_ - a_);
340 scalar d3 = (c_ - a_) & (c_ - b_);
346 scalar c = c1 + c2 + c3;
357 ((c2 + c3)*a_ + (c3 + c1)*b_ + (c1 + c2)*c_)/(2*
c)
362template<
class Po
int,
class Po
intRef>
365 const scalar d1 = (c_ - a_) & (b_ - a_);
366 const scalar d2 = -(c_ - b_) & (b_ - a_);
367 const scalar d3 = (c_ - a_) & (c_ - b_);
369 const scalar denom = d2*d3 + d3*d1 + d1*d2;
378 const scalar
a = (d1 + d2)*(d2 + d3)*(d3 + d1) / denom;
383template<
class Po
int,
class Po
intRef>
386 scalar c = circumRadius();
398template<
class Po
int,
class Po
intRef>
406 ((t.a_ - a_) & ((b_ - a_)^(c_ - a_)))
407 + ((t.b_ - b_) & ((c_ - b_)^(t.a_ - b_)))
408 + ((c_ - t.c_) & ((t.b_ - t.c_)^(t.a_ - t.c_)))
410 + ((t.a_ - a_) & ((b_ - a_)^(c_ - a_)))
411 + ((b_ - t.b_) & ((t.a_ - t.b_)^(t.c_ - t.b_)))
412 + ((c_ - t.c_) & ((b_ - t.c_)^(t.a_ - t.c_)))
417template<
class Po
int,
class Po
intRef>
424 Point aRel = a_ - refPt;
425 Point bRel = b_ - refPt;
426 Point cRel = c_ - refPt;
430 aRel.x(), aRel.y(), aRel.z(),
431 bRel.x(), bRel.y(), bRel.z(),
432 cRel.x(), cRel.y(), cRel.z()
435 scalar a =
Foam::mag((b_ - a_)^(c_ - a_));
446 + ((aRel + bRel + cRel) & (aRel + bRel + cRel))
454template<
class Po
int,
class Po
intRef>
461template<
class Po
int,
class Po
intRef>
467 return bary[0]*a_ + bary[1]*b_ + bary[2]*c_;
471template<
class Po
int,
class Po
intRef>
483template<
class Po
int,
class Po
intRef>
497 scalar d00 = v0 & v0;
498 scalar d01 = v0 & v1;
499 scalar d11 = v1 & v1;
500 scalar d20 = v2 & v0;
501 scalar d21 = v2 & v1;
503 scalar denom = d00*d11 - d01*d01;
514 bary[1] = (d11*d20 - d01*d21)/denom;
515 bary[2] = (d00*d21 - d01*d20)/denom;
516 bary[0] = 1.0 - bary[1] - bary[2];
522template<
class Po
int,
class Po
intRef>
534 (((a_ - origin) & normal) < 0 ? 1 : 2)
535 | (((b_ - origin) & normal) < 0 ? 1 : 2)
536 | (((c_ - origin) & normal) < 0 ? 1 : 2)
542template<
class Po
int,
class Po
intRef>
554 (a_[axis] < origin[axis] ? 1 : 2)
555 | (b_[axis] < origin[axis] ? 1 : 2)
556 | (c_[axis] < origin[axis] ? 1 : 2)
562template<
class Po
int,
class Po
intRef>
615 hit = fastInter.hit();
619 pInter = fastInter.point();
625 pInter =
p + (q1&v)*q1;
630 scalar dist = q1 & (pInter -
p);
632 const scalar planarPointTol =
648 && ((q1 & areaNormal()) < -VSMALL)
677template<
class Po
int,
class Po
intRef>
686 const vector edge1 = b_ - a_;
687 const vector edge2 = c_ - a_;
690 const vector pVec = dir ^ edge2;
693 const scalar
det = edge1 & pVec;
701 if (
det < ROOTVSMALL)
710 if (
det > -ROOTVSMALL &&
det < ROOTVSMALL)
717 const scalar inv_det = 1.0 /
det;
720 const vector tVec = orig-a_;
723 const scalar u = (tVec & pVec)*inv_det;
725 if (u < -tol || u > 1.0+tol)
732 const vector qVec = tVec ^ edge1;
735 const scalar v = (dir & qVec) * inv_det;
737 if (v < -tol || u + v > 1.0+tol)
744 const scalar t = (edge2 & qVec) * inv_det;
759template<
class Po
int,
class Po
intRef>
778 if (d1 <= 0.0 && d2 <= 0.0)
792 if (d3 >= 0.0 && d4 <= d3)
802 scalar vc = d1*d4 - d3*d2;
804 if (vc <= 0.0 && d1 >= 0.0 && d3 <= 0.0)
806 if ((d1 - d3) < ROOTVSMALL)
815 scalar v = d1/(d1 - d3);
817 point nearPt = a_ + v*ab;
828 if (d6 >= 0.0 && d5 <= d6)
838 scalar vb = d5*d2 - d1*d6;
840 if (vb <= 0.0 && d2 >= 0.0 && d6 <= 0.0)
842 if ((d2 - d6) < ROOTVSMALL)
851 scalar w = d2/(d2 - d6);
853 point nearPt = a_ + w*ac;
860 scalar va = d3*d6 - d5*d4;
862 if (va <= 0.0 && (d4 - d3) >= 0.0 && (d5 - d6) >= 0.0)
864 if (((d4 - d3) + (d5 - d6)) < ROOTVSMALL)
874 scalar w = (d4 - d3)/((d4 - d3) + (d5 - d6));
876 point nearPt = b_ + w*(c_ - b_);
885 if ((va + vb + vc) < ROOTVSMALL)
889 point nearPt = centre();
895 scalar denom = 1.0/(va + vb + vc);
896 scalar v = vb * denom;
897 scalar w = vc * denom;
901 point nearPt = a_ + ab*v + ac*w;
908template<
class Po
int,
class Po
intRef>
916 label nearLabel = -1;
922template<
class Po
int,
class Po
intRef>
934template<
class Po
int,
class Po
intRef>
955 if (triInfo.distance() > 1)
960 scalar dist = triInfo.point().dist(lnInfo.
point());
962 triInfo.setMiss(
true);
963 triInfo.setDistance(dist);
965 else if (triInfo.distance() < 0)
970 scalar dist = triInfo.point().dist(lnInfo.
point());
972 triInfo.setMiss(
true);
973 triInfo.setDistance(dist);
980 triInfo.setDistance(0);
987 point nearestEdgePoint;
988 point nearestLinePoint;
990 scalar minDist =
ln.nearestDist
1000 scalar dist =
ln.nearestDist
1009 nearestEdgePoint = triEdgePoint;
1010 nearestLinePoint = linePoint;
1018 scalar dist =
ln.nearestDist
1027 nearestEdgePoint = triEdgePoint;
1028 nearestLinePoint = linePoint;
1034 triInfo.setDistance(minDist);
1035 triInfo.setMiss(
false);
1036 triInfo.setPoint(nearestEdgePoint);
1039 if (
Foam::mag(nearestLinePoint-
ln.start()) < SMALL)
1044 else if (
Foam::mag(nearestLinePoint-
ln.end()) < SMALL)
1058template<
class Po
int,
class Po
intRef>
1067 return ((dist < -tol) ? -1 : (dist > tol) ? +1 : 0);
1071template<
class Po
int,
class Po
intRef>
1081template<
class Po
int,
class Po
intRef>
1093template<
class Po
int,
class Po
intRef>
1099 tris_[nTris_++] = tri;
1103template<
class Po
int,
class Po
intRef>
1104inline Foam::point Foam::triangle<Point, PointRef>::planeIntersection
1112 return (d[posI]*t[negI] - d[negI]*t[posI])/(-d[negI] + d[posI]);
1118template<
class Po
int,
class Po
intRef>
1122 triangle<Point, PointRef>& t
1125 is.readBegin(
"triangle");
1126 is >> t.a_ >> t.b_ >> t.c_;
1127 is.readEnd(
"triangle");
1134template<
class Po
int,
class Po
intRef>
1135inline Foam::Ostream& Foam::operator<<
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
A 1D vector of objects of type <T> with a fixed length <N>.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
An ordered pair of two objects of type <T> with first() and second() elements.
scalar distance() const noexcept
Return distance to hit.
void setPoint(const point_type &p)
Set the point.
void setDistance(const scalar d) noexcept
Set the distance.
bool hit() const noexcept
Is there a hit.
void setMiss(const bool eligible) noexcept
Set the hit status off and set the eligible miss status.
const point_type & point() const noexcept
Return the point, no checks.
const point_type & hitPoint() const
Return the hit point. Fatal if not hit.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
scalar dist(const Vector< Cmpt > &v2) const
The L2-norm distance from another vector. The mag() of the difference.
static scalar planarTol()
Return planar tolerance.
@ BEGIN_LIST
Begin list [isseparator].
@ END_LIST
End list [isseparator].
Triangle point storage. Default constructable (triangle is not).
void flip()
Flip triangle orientation by swapping second and third vertices.
triPointRef tri() const
Return as triangle reference.
const point & c() const noexcept
The third vertex.
vector vecA() const
Edge vector opposite point a(): from b() to c().
vector areaNormal() const
The area normal - with magnitude equal to area of triangle.
const point & b() const noexcept
The second vertex.
Pair< point > box() const
The enclosing (bounding) box for the triangle.
vector unitNormal() const
Return unit normal.
scalar mag() const
The magnitude of the triangle area.
vector vecC() const
Edge vector opposite point c(): from a() to b().
triPoints()=default
Default construct.
const point & a() const noexcept
The first vertex.
scalar magSqr() const
The magnitude squared of the triangle area.
point centre() const
Return centre (centroid).
vector vecB() const
Edge vector opposite point b(): from c() to a().
A triangle primitive used to calculate face normals and swept volumes. Uses referred points.
Point circumCentre() const
Return circum-centre.
barycentric2D pointToBarycentric(const point &pt) const
Calculate the barycentric coordinates from the given point.
static vector areaNormal(const Point &p0, const Point &p1, const Point &p2)
The area normal for a triangle defined by three points (right-hand rule). Magnitude equal to area of ...
scalar circumRadius() const
Return circum-radius.
Point barycentricToPoint(const barycentric2D &bary) const
Calculate the point from the given barycentric coordinates.
scalar sweptVol(const triangle &t) const
Return swept-volume.
vector areaNormal() const
The area normal - with magnitude equal to area of triangle.
Point centre() const
Return centre (centroid).
Point vecB() const
Edge vector opposite point b(): from c() to a().
FixedList< triPoints, 27 > triIntersectionList
Storage type for triangles originating from intersecting triangle with another triangle.
pointHit nearestPointClassify(const point &p, label &nearType, label &nearLabel) const
Find the nearest point to p on the triangle and classify it:
triangle(const Point &p0, const Point &p1, const Point &p2)
Construct from three points.
static Point centre(const Point &p0, const Point &p1, const Point &p2)
The centre (centroid) of three points.
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.
vector unitNormal() const
Return unit normal.
tensor inertia(PointRef refPt=Zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
scalar mag() const
The magnitude of the triangle area.
bool intersects(const point &origin, const vector &normal) const
Fast intersection detection with a plane.
pointHit nearestPoint(const point &p) const
Return nearest point to p on triangle.
Point randomPoint(Random &rndGen) const
Return a random point on the triangle from a uniform distribution.
Pair< Point > box() const
The enclosing (bounding) box for the triangle.
pointHit intersection(const point &p, const vector &q, const intersection::algorithm alg, const scalar tol=0.0) const
Fast intersection with a ray.
bool classify(const point &p, label &nearType, label &nearLabel) const
Classify nearest point to p in triangle plane.
scalar magSqr() const
The magnitude squared of the triangle area.
int sign(const point &p, const scalar tol=SMALL) const
The sign for which side of the face plane the point is on.
static vector unitNormal(const Point &p0, const Point &p1, const Point &p2)
The unit normal for a triangle defined by three points (right-hand rule).
scalar quality() const
Return quality: Ratio of triangle and circum-circle.
const point & c() const noexcept
Point vecA() const
Edge vector opposite point a(): from b() to c().
const Point & a() const noexcept
The first vertex.
Point vecC() const
Edge vector opposite point c(): from a() to b().
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const volScalarField & p0
OBJstream os(runTime.globalPath()/outputName)
@ NONE
No type, or default initialized type.
dimensionedScalar det(const dimensionedSphericalTensor &dt)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
barycentric2D barycentric2D01(Random &rndGen)
Generate a random barycentric coordinate within the unit triangle.
PointHit< point > pointHit
A PointHit with a 3D point.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Barycentric2D< scalar > barycentric2D
A scalar version of the templated Barycentric2D.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
static const Identity< scalar > I
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
line< point, const point & > linePointRef
A line using referred points.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
vector point
Point is a vector.
static constexpr const zero Zero
Global zero (0).
triangle< point, const point & > triPointRef
A triangle using referred points.
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy the source to the destination (recursively if necessary).
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
constexpr char nl
The newline '\n' character (0x0a).
storeOp(triIntersectionList &, label &)
triIntersectionList & tris_