A triangle primitive used to calculate face normals and swept volumes. More...
Data Structures | |
class | dummyOp |
Dummy. More... | |
class | storeOp |
Store resulting tris. More... | |
class | sumAreaOp |
Sum resulting areas. More... | |
Public Types | |
enum | proxType { NONE, POINT, EDGE } |
Return types for classify. More... | |
typedef FixedList< triPoints, 27 > | triIntersectionList |
Storage type for triangles originating from intersecting triangle. More... | |
Public Member Functions | |
triangle (const Point &a, const Point &b, const Point &c) | |
Construct from three points. More... | |
triangle (const FixedList< Point, 3 > &) | |
Construct from three points. More... | |
triangle (const UList< Point > &, const FixedList< label, 3 > &indices) | |
Construct from three points in the list of points. More... | |
triangle (Istream &) | |
Construct from Istream. More... | |
const Point & | a () const |
Return first vertex. More... | |
const Point & | b () const |
Return second vertex. More... | |
const Point & | c () const |
Return third vertex. More... | |
Point | centre () const |
Return centre (centroid) More... | |
scalar | mag () const |
Return scalar magnitude. More... | |
vector | normal () const |
Return vector normal. More... | |
Point | circumCentre () const |
Return circum-centre. More... | |
scalar | circumRadius () const |
Return circum-radius. More... | |
scalar | quality () const |
Return quality: Ratio of triangle and circum-circle. More... | |
scalar | sweptVol (const triangle &t) const |
Return swept-volume. More... | |
tensor | inertia (PointRef refPt=vector::zero, scalar density=1.0) const |
Return the inertia tensor, with optional reference. More... | |
Point | randomPoint (Random &rndGen) const |
Return a random point on the triangle from a uniform. More... | |
Point | randomPoint (cachedRandom &rndGen) const |
Return a random point on the triangle from a uniform. More... | |
scalar | barycentric (const point &pt, List< scalar > &bary) const |
Calculate the barycentric coordinates of the given. More... | |
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. More... | |
pointHit | intersection (const point &p, const vector &q, const intersection::algorithm alg, const scalar tol=0.0) const |
Fast intersection with a ray. More... | |
pointHit | nearestPointClassify (const point &p, label &nearType, label &nearLabel) const |
Find the nearest point to p on the triangle and classify it: More... | |
pointHit | nearestPoint (const point &p) const |
Return nearest point to p on triangle. More... | |
bool | classify (const point &p, label &nearType, label &nearLabel) const |
Classify nearest point to p in triangle plane. More... | |
pointHit | nearestPoint (const linePointRef &edge, pointHit &edgePoint) const |
Return nearest point to line on triangle. Returns hit if. More... | |
template<class AboveOp , class BelowOp > | |
void | sliceWithPlane (const plane &pl, AboveOp &aboveOp, BelowOp &belowOp) const |
Decompose triangle into triangles above and below plane. More... | |
template<class InsideOp , class OutsideOp > | |
void | triangleOverlap (const vector &n, const triangle< Point, PointRef > &tri, InsideOp &insideOp, OutsideOp &outsideOp) const |
Decompose triangle into triangles inside and outside. More... | |
Static Private Member Functions | |
static point | planeIntersection (const FixedList< scalar, 3 > &d, const triPoints &t, const label negI, const label posI) |
Helper: calculate intersection point. More... | |
template<class AboveOp , class BelowOp > | |
static void | triSliceWithPlane (const plane &pl, const triPoints &tri, AboveOp &aboveOp, BelowOp &belowOp) |
Helper: slice triangle with plane. More... | |
Private Attributes | |
PointRef | a_ |
PointRef | b_ |
PointRef | c_ |
Friends | |
Istream & | operator>> (Istream &, triangle &) |
Ostream & | operator (Ostream &, const triangle &) |
A triangle primitive used to calculate face normals and swept volumes.
Definition at line 59 of file triangle.H.
typedef FixedList<triPoints, 27> triIntersectionList |
Storage type for triangles originating from intersecting triangle.
with another triangle
Definition at line 91 of file triangle.H.
enum proxType |
Construct from three points.
Definition at line 35 of file triangleI.H.
Construct from three points in the list of points.
The indices could be from triFace etc.
|
inlinestaticprivate |
Helper: calculate intersection point.
Definition at line 869 of file triangleI.H.
|
inlinestaticprivate |
Helper: slice triangle with plane.
Definition at line 35 of file triangle.C.
|
inline |
Return first vertex.
Definition at line 83 of file triangleI.H.
Referenced by triSurfaceTools::calcInterpolationWeights(), meshSurfaceCheckInvertedVertices::checkVertices(), Foam::help::doTrianglesIntersect(), Foam::help::doTrianglesOverlap(), particle< Type >::hitWallFaces(), Foam::help::nearestPointOnTheTriangle(), offsetSurface::operator()(), triangle< Foam::Vector, Foam::Vector >::triangleOverlap(), Foam::help::triLineIntersection(), and isoSurface::trimToPlanes().
|
inline |
Return second vertex.
Definition at line 89 of file triangleI.H.
Referenced by triSurfaceTools::calcInterpolationWeights(), meshSurfaceCheckInvertedVertices::checkVertices(), Foam::help::doTrianglesIntersect(), Foam::help::doTrianglesOverlap(), particle< Type >::hitWallFaces(), Foam::help::nearestPointOnTheTriangle(), offsetSurface::operator()(), triangle< Foam::Vector, Foam::Vector >::triangleOverlap(), Foam::help::triLineIntersection(), and isoSurface::trimToPlanes().
|
inline |
Return third vertex.
Definition at line 95 of file triangleI.H.
Referenced by triSurfaceTools::calcInterpolationWeights(), meshSurfaceCheckInvertedVertices::checkVertices(), Foam::help::doTrianglesIntersect(), Foam::help::doTrianglesOverlap(), particle< Type >::hitWallFaces(), Foam::help::nearestPointOnTheTriangle(), offsetSurface::operator()(), triangle< Foam::Vector, Foam::Vector >::triangleOverlap(), Foam::help::triLineIntersection(), and isoSurface::trimToPlanes().
|
inline |
Return centre (centroid)
Definition at line 102 of file triangleI.H.
Referenced by Foam::help::doTrianglesOverlap(), wallBoundedStreamLine::findNearestTet(), triangulateNonPlanarBaseFaces::findNonPlanarBoundaryFaces(), knuppMetric::knuppMetric(), momentOfInertia::massPropertiesShell(), tetMeshOptimisation::optimiseBoundarySurfaceLaplace(), polyMesh::pointInCell(), quadricMetric::quadricMetric(), and wallBoundedStreamLine::track().
|
inline |
Return scalar magnitude.
Definition at line 109 of file triangleI.H.
Referenced by pointLinear< Type >::correction(), FreeStream< CloudType >::inflow(), momentOfInertia::massPropertiesShell(), and triangle::sumAreaOp::operator()().
|
inline |
Return vector normal.
Definition at line 116 of file triangleI.H.
Referenced by primitiveMeshGeometry::checkFaceTwist(), polyMeshGeometry::checkFaceTwist(), polyMeshGeometry::checkTriangleTwist(), meshSurfaceCheckInvertedVertices::checkVertices(), Foam::help::doTrianglesIntersect(), Foam::help::doTrianglesOverlap(), triangulateNonPlanarBaseFaces::findNonPlanarBoundaryFaces(), particle< Type >::hitWallFaces(), solidParticle::hitWallPatch(), knuppMetric::knuppMetric(), tetMeshOptimisation::optimiseBoundaryVolumeOptimizer(), volumeOptimizer::optimiseSteepestDescent(), KinematicCloud< CloudType >::patchData(), polyMesh::pointInCell(), quadricMetric::quadricMetric(), partTet::Sa(), tetrahedron::Sa(), partTet::Sb(), tetrahedron::Sb(), partTet::Sc(), tetrahedron::Sc(), partTet::Sd(), tetrahedron::Sd(), wallBoundedParticle::trackToEdge(), and particle< Type >::trackToFace().
|
inline |
Return circum-centre.
Definition at line 123 of file triangleI.H.
Referenced by main().
|
inline |
Return circum-radius.
Definition at line 150 of file triangleI.H.
Referenced by main().
|
inline |
Return quality: Ratio of triangle and circum-circle.
area, scaled so that an equilateral triangle has a quality of 1
Definition at line 174 of file triangleI.H.
|
inline |
Return swept-volume.
Definition at line 190 of file triangleI.H.
|
inline |
Return the inertia tensor, with optional reference.
point and density specification
Definition at line 209 of file triangleI.H.
Return a random point on the triangle from a uniform.
distribution
Definition at line 245 of file triangleI.H.
Referenced by FreeStream< CloudType >::inflow(), and patchInjectionBase::setPositionAndCell().
|
inline |
Return a random point on the triangle from a uniform.
distribution
Definition at line 262 of file triangleI.H.
Calculate the barycentric coordinates of the given.
point, in the same order as a, b, c. Returns the determinant of the solution.
Definition at line 281 of file triangleI.H.
Referenced by cellPointWeight::findTriangle(), offsetSurface::operator()(), and isoSurface::trimToBox().
|
inline |
Return point intersection with a ray.
For a hit, the distance is signed. Positive number represents the point in front of triangle. In case of miss pointHit is set to nearest point on triangle and its distance to the distance between the original point and the plane intersection point
Definition at line 322 of file triangleI.H.
|
inline |
Fast intersection with a ray.
For a hit, the pointHit.distance() is the line parameter t : intersection=p+t*q. Only defined for VISIBLE, FULL_RAY or HALF_RAY. tol increases the virtual size of the triangle by a relative factor.
Definition at line 438 of file triangleI.H.
Referenced by mappedPatchBase::facePoint(), and particle< Type >::hitWallFaces().
Foam::pointHit nearestPointClassify | ( | const point & | p, |
label & | nearType, | ||
label & | nearLabel | ||
) | const |
Find the nearest point to p on the triangle and classify it:
+ near point (nearType=POINT, nearLabel=0, 1, 2) + near edge (nearType=EDGE, nearLabel=0, 1, 2) Note: edges are counted from starting vertex so e.g. edge 2 is from f[2] to f[0]
Definition at line 521 of file triangleI.H.
Referenced by triSurfaceTools::calcInterpolationWeights(), and face::nearestPointClassify().
|
inline |
Return nearest point to p on triangle.
Definition at line 670 of file triangleI.H.
Referenced by cellPointWeight::findTriangle(), and tetrahedron::nearestPoint().
Classify nearest point to p in triangle plane.
w.r.t. triangle edges and points. Returns inside (true)/outside (false).
Definition at line 684 of file triangleI.H.
Referenced by ParticleCollector< CloudType >::collectParcelPolygon().
|
inline |
Return nearest point to line on triangle. Returns hit if.
point is inside triangle. Sets edgePoint to point on edge (hit if nearest is inside line)
|
inline |
Decompose triangle into triangles above and below plane.
Definition at line 113 of file triangle.C.
|
inline |
Decompose triangle into triangles inside and outside.
(with respect to user provided normal) other triangle.
Definition at line 126 of file triangle.C.
|
private |
Definition at line 142 of file triangle.H.
Referenced by triangle< Foam::Vector, Foam::Vector >::sweptVol().
|
private |
Definition at line 142 of file triangle.H.
Referenced by triangle< Foam::Vector, Foam::Vector >::sweptVol().
|
private |
Definition at line 142 of file triangle.H.
Referenced by triangle< Foam::Vector, Foam::Vector >::sweptVol().
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.