Go to the documentation of this file.
45 return i ? i-1 : size-1;
65 vector vec(nextPt - thisPt);
66 vec /=
mag(vec) + VSMALL;
90 if (
sin < -ROOTVSMALL)
110 const point& rayOrigin,
126 if ((posOnEdge < 0) || (posOnEdge > 1))
133 point intersectPt = p1 + posOnEdge * (p2 - p1);
135 if (((intersectPt - rayOrigin) & rayDir) < 0)
166 if ((area01Pt > 0) && (area12Pt > 0) && (area20Pt > 0))
170 else if ((area01Pt < 0) && (area12Pt < 0) && (area20Pt < 0))
190 const label startIndex,
198 const vector& rightE = edges[right(
f.
size(), startIndex)];
199 const vector leftE = -edges[left(
f.
size(), startIndex)];
202 scalar cosHalfAngle = GREAT;
203 scalar sinHalfAngle = GREAT;
204 calcHalfAngle(
normal, rightE, leftE, cosHalfAngle, sinHalfAngle);
208 + sinHalfAngle*(
normal ^ rightE)
212 rayDir /=
mag(rayDir) + VSMALL;
219 label faceVertI =
f.fcIndex(startIndex);
223 scalar minPosOnEdge = GREAT;
242 minIndex = faceVertI;
243 minPosOnEdge = posOnEdge;
246 faceVertI =
f.fcIndex(faceVertI);
261 const label leftIndex = minIndex;
262 const label rightIndex =
f.fcIndex(minIndex);
268 if (
mag(minPosOnEdge) < edgeRelTol &&
f.fcIndex(startIndex) != leftIndex)
277 mag(minPosOnEdge - 1) < edgeRelTol
278 &&
f.fcIndex(rightIndex) != startIndex
295 scalar maxCos = -GREAT;
298 faceVertI =
f.fcIndex(
f.fcIndex(startIndex));
305 (faceVertI == leftIndex)
306 || (faceVertI == rightIndex)
307 || (triangleContainsPoint(
normal, startPt, leftPt, rightPt, pt))
312 vector edgePt0 = pt - startPt;
313 edgePt0 /=
mag(edgePt0);
315 scalar
cos = rayDir & edgePt0;
319 minIndex = faceVertI;
322 faceVertI =
f.fcIndex(faceVertI);
331 if (
f.fcIndex(startIndex) != leftIndex)
360 scalar minCos = GREAT;
365 const vector& rightEdge = edges[right(size, fp)];
366 const vector leftEdge = -edges[left(size, fp)];
368 if (((rightEdge ^ leftEdge) &
normal) < ROOTVSMALL)
370 scalar
cos = rightEdge & leftEdge;
386 const vector& rightEdge = edges[right(size, fp)];
387 const vector leftEdge = -edges[left(size, fp)];
389 scalar
cos = rightEdge & leftEdge;
420 <<
"Illegal face:" <<
f
429 triFace& tri = operator[](triI++);
463 if (index1 != -1 && index2 != -1)
470 startIndex =
f.fcIndex(startIndex);
473 if (index1 == -1 || index2 == -1)
480 scalar maxCos = -GREAT;
484 const vector& rightEdge = edges[right(size, fp)];
485 const vector leftEdge = -edges[left(size, fp)];
487 scalar
cos = rightEdge & leftEdge;
496 <<
"Cannot find valid diagonal on face " <<
f
499 <<
"Returning naive triangulation starting from "
500 <<
f[maxIndex] <<
" which might not be correct for a"
501 <<
" concave or warped face" <<
endl;
504 label fp =
f.fcIndex(maxIndex);
506 for (
label i = 0; i < size-2; i++)
508 label nextFp =
f.fcIndex(fp);
510 triFace& tri = operator[](triI++);
511 tri[0] =
f[maxIndex];
523 <<
"Cannot find valid diagonal on face " <<
f
526 <<
"Returning empty triFaceList" <<
endl;
541 diff = index2 - index1;
546 diff = index2 + size - index1;
552 if (nPoints1 == size || nPoints2 == size)
555 <<
"Illegal split of face:" <<
f
557 <<
" at indices " << index1 <<
" and " << index2
563 face face1(nPoints1);
565 label faceVertI = index1;
566 for (
int i = 0; i < nPoints1; i++)
568 face1[i] =
f[faceVertI];
569 faceVertI =
f.fcIndex(faceVertI);
573 face face2(nPoints2);
576 for (
int i = 0; i < nPoints2; i++)
578 face2[i] =
f[faceVertI];
579 faceVertI =
f.fcIndex(faceVertI);
619 avgNormal /=
mag(avgNormal) + VSMALL;
623 bool valid = split(fallBack,
points,
f, avgNormal, triI);
645 bool valid = split(fallBack,
points,
f,
n, triI);
points setSize(newPointi)
bool hit() const
Is there a hit.
bool split(const bool fallBack, const pointField &points, const face &f, const vector &normal, label &triI)
Split face f into triangles. Handles all simple (convex & concave)
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
This class describes the interaction of a face and a point. It carries the info of a successful hit a...
dimensionedScalar sin(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< scalar > mag(const dimensioned< Type > &)
static label left(const label size, label i)
Edge to the left of face vertex i.
Field< vector > vectorField
Specialisation of Field<T> for vector.
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
static void findDiagonal(const pointField &points, const face &f, const vectorField &edges, const vector &normal, const label startIndex, label &index1, label &index2)
Starting from startIndex find diagonal. Return in index1, index2.
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
static pointHit rayEdgeIntersect(const vector &normal, const point &rayOrigin, const vector &rayDir, const point &p1, const point &p2, scalar &posOnEdge)
Calculate intersection point between edge p1-p2 and ray (in 2D).
faceTriangulation()
Construct null.
scalar distance() const
Return distance to hit.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
void setPoint(const Point &p)
static void calcHalfAngle(const vector &normal, const vector &e0, const vector &e1, scalar &cosHalfAngle, scalar &sinHalfAngle)
Calculates half angle components of angle from e0 to e1.
static bool triangleContainsPoint(const vector &n, const point &p0, const point &p1, const point &p2, const point &pt)
static const scalar edgeRelTol
Relative tolerance on edge.
static tmp< vectorField > calcEdges(const face &, const pointField &)
Calculate normalized edge vectors.
static label findStart(const face &f, const vectorField &edges, const vector &normal)
Find label of vertex to start splitting from. This will be the.
List< triFace > triFaceList
list of triFaces
errorManip< error > abort(error &err)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
A triangular face using a FixedList of labels corresponding to mesh vertices.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionedScalar e
Elementary charge.
static label right(const label size, label i)
Edge to the right of face vertex i.
A List with indirect addressing.
A face is a list of labels corresponding to mesh vertices.
void size(const label)
Override size to be inconsistent with allocated storage.
scalar normalIntersect(const point &pnt0, const vector &dir) const
Return cut coefficient for plane and line defined by.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
#define WarningInFunction
Report a warning using Foam::Warning.
A normal distribution model.
triangle< point, const point & > triPointRef
dimensionedScalar cos(const dimensionedScalar &ds)
void setDistance(const scalar d)