34const Foam::scalar Foam::faceTriangulation::edgeRelTol = 1
e-6;
38Foam::label Foam::faceTriangulation::right(
const label, label i)
45Foam::label Foam::faceTriangulation::left(
const label size, label i)
47 return i ? i-1 : size-1;
60 auto& edges = tedges.ref();
67 vector vec(nextPt - thisPt);
69 edges[i] = vec.normalise();
77void Foam::faceTriangulation::calcHalfAngle
89 scalar
sin = (e0 ^ e1) & normal;
91 if (
sin < -ROOTVSMALL)
111 const point& rayOrigin,
122 const vector y = normal ^ rayDir;
124 posOnEdge = plane(rayOrigin,
y).normalIntersect(p1, (p2-p1));
127 if ((posOnEdge < 0) || (posOnEdge > 1))
134 point intersectPt = p1 + posOnEdge * (p2 - p1);
136 if (((intersectPt - rayOrigin) & rayDir) < 0)
144 result.setPoint(intersectPt);
145 result.setDistance(
mag(intersectPt - rayOrigin));
154bool Foam::faceTriangulation::triangleContainsPoint
164 const scalar area12Pt =
triPointRef(p1, p2, pt).areaNormal() &
n;
167 if ((area01Pt > 0) && (area12Pt > 0) && (area20Pt > 0))
171 else if ((area01Pt < 0) && (area12Pt < 0) && (area20Pt < 0))
183void Foam::faceTriangulation::findDiagonal
189 const label startIndex,
197 const vector& rightE = edges[right(
f.
size(), startIndex)];
198 const vector leftE = -edges[left(
f.
size(), startIndex)];
201 scalar cosHalfAngle = GREAT;
202 scalar sinHalfAngle = GREAT;
203 calcHalfAngle(normal, rightE, leftE, cosHalfAngle, sinHalfAngle);
207 + sinHalfAngle*(normal ^ rightE)
218 label faceVertI =
f.
fcIndex(startIndex);
222 scalar minPosOnEdge = GREAT;
224 for (label i = 0; i <
f.
size() - 2; i++)
238 if (inter.hit() && inter.distance() < minInter.distance())
241 minIndex = faceVertI;
242 minPosOnEdge = posOnEdge;
260 const label leftIndex = minIndex;
261 const label rightIndex =
f.
fcIndex(minIndex);
267 if (
mag(minPosOnEdge) < edgeRelTol &&
f.
fcIndex(startIndex) != leftIndex)
276 mag(minPosOnEdge - 1) < edgeRelTol
277 &&
f.
fcIndex(rightIndex) != startIndex
294 scalar maxCos = -GREAT;
298 for (label i = 0; i <
f.
size() - 3; i++)
304 (faceVertI == leftIndex)
305 || (faceVertI == rightIndex)
306 || (triangleContainsPoint(normal, startPt, leftPt, rightPt, pt))
311 vector edgePt0 = (pt - startPt);
314 scalar
cos = rayDir & edgePt0;
318 minIndex = faceVertI;
330 if (
f.
fcIndex(startIndex) != leftIndex)
350Foam::label Foam::faceTriangulation::findStart
357 const label size =
f.
size();
359 scalar minCos = GREAT;
364 const vector& rightEdge = edges[right(size, fp)];
365 const vector leftEdge = -edges[left(size, fp)];
367 if (((rightEdge ^ leftEdge) & normal) < ROOTVSMALL)
369 scalar
cos = rightEdge & leftEdge;
385 const vector& rightEdge = edges[right(size, fp)];
386 const vector leftEdge = -edges[left(size, fp)];
388 scalar
cos = rightEdge & leftEdge;
405bool Foam::faceTriangulation::split
414 const label size =
f.
size();
419 <<
"Illegal face:" <<
f
420 <<
" with points " << UIndirectList<point>(
points,
f)
428 triFace& tri = operator[](triI++);
440 tmp<vectorField> tedges(calcEdges(
f,
points));
443 label startIndex = findStart(
f, edges, normal);
462 if (index1 != -1 && index2 != -1)
472 if (index1 == -1 || index2 == -1)
479 scalar maxCos = -GREAT;
483 const vector& rightEdge = edges[right(size, fp)];
484 const vector leftEdge = -edges[left(size, fp)];
486 scalar
cos = rightEdge & leftEdge;
495 <<
"Cannot find valid diagonal on face " <<
f
496 <<
" with points " << UIndirectList<point>(
points,
f)
498 <<
"Returning naive triangulation starting from "
499 <<
f[maxIndex] <<
" which might not be correct for a"
500 <<
" concave or warped face" <<
endl;
505 for (label i = 0; i < size-2; i++)
509 triFace& tri = operator[](triI++);
510 tri[0] =
f[maxIndex];
522 <<
"Cannot find valid diagonal on face " <<
f
523 <<
" with points " << UIndirectList<point>(
points,
f)
525 <<
"Returning empty triFaceList" <<
endl;
540 diff = index2 - index1;
545 diff = index2 + size - index1;
548 label nPoints1 =
diff + 1;
549 label nPoints2 = size -
diff + 1;
551 if (nPoints1 == size || nPoints2 == size)
554 <<
"Illegal split of face:" <<
f
555 <<
" with points " << UIndirectList<point>(
points,
f)
556 <<
" at indices " << index1 <<
" and " << index2
562 face face1(nPoints1);
564 label faceVertI = index1;
565 for (
int i = 0; i < nPoints1; i++)
567 face1[i] =
f[faceVertI];
572 face face2(nPoints2);
575 for (
int i = 0; i < nPoints2; i++)
577 face2[i] =
f[faceVertI];
619 bool valid = split(fallBack,
points,
f, avgNormal, triI);
640 bool valid = split(fallBack,
points,
f,
n, triI);
static bool split(const std::string &line, std::string &key, std::string &val)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void setSize(const label n)
Alias for resize()
Describes the interaction of a face and a point. It carries the info of a successful hit and (if succ...
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
void size(const label n)
Older name for setAddressableSize.
label fcIndex(const label i) const noexcept
faceTriangulation()
Construct null.
A face is a list of labels corresponding to mesh vertices.
A class for managing temporary objects.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const volScalarField & p0
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define WarningInFunction
Report a warning using Foam::Warning.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
dimensionedScalar sin(const dimensionedScalar &ds)
vectorField pointField
pointField is a vectorField.
vector point
Point is a vector.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
triangle< point, const point & > triPointRef
A triangle using referred points.
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
errorManip< error > abort(error &err)
PointHit< point > pointHit
A PointIndexHit for 3D points.
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0)
dimensionedScalar cos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a)
#define forAll(list, i)
Loop across all elements in list.