libigl v2.5.0
Loading...
Searching...
No Matches
igl::predicates Namespace Reference

Enumerations

enum class  Orientation {
  POSITIVE =1 , INSIDE =1 , NEGATIVE =-1 , OUTSIDE =-1 ,
  COLLINEAR =0 , COPLANAR =0 , COCIRCULAR =0 , COSPHERICAL =0 ,
  DEGENERATE =0
}
 Types of orientations and other predicate results. More...
 

Functions

template<typename DerivedV , typename DerivedF >
void delaunay_triangulation (const Eigen::MatrixBase< DerivedV > &V, Eigen::PlainObjectBase< DerivedF > &F)
 Given a set of points in 2D, return a Delaunay triangulation of these points using predicates.
 
template<typename DerivedP , typename DerivedRT , typename DerivedF , typename DerivedI >
void ear_clipping (const Eigen::MatrixBase< DerivedP > &P, const Eigen::MatrixBase< DerivedRT > &RT, Eigen::PlainObjectBase< DerivedF > &eF, Eigen::PlainObjectBase< DerivedI > &I)
 Implementation of ear clipping triangulation algorithm for a 2D polygon.
 
template<typename DerivedP , typename DerivedF >
bool ear_clipping (const Eigen::MatrixBase< DerivedP > &P, Eigen::PlainObjectBase< DerivedF > &eF)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
template<typename DerivedV , typename DerivedF >
void lexicographic_triangulation (const Eigen::MatrixBase< DerivedV > &V, Eigen::PlainObjectBase< DerivedF > &F)
 Given a set of points in 2D, return a lexicographic triangulation of these points using predicates.
 
template<typename DerivedP , typename DerivedQ >
bool point_inside_convex_polygon (const Eigen::MatrixBase< DerivedP > &P, const Eigen::MatrixBase< DerivedQ > &q)
 check whether 2d point lies inside 2d convex polygon
 
template<typename DerivedV , typename DerivedI , typename DerivedC , typename DerivedF , typename DerivedJ >
void polygons_to_triangles (const Eigen::MatrixBase< DerivedV > &V, const Eigen::MatrixBase< DerivedI > &I, const Eigen::MatrixBase< DerivedC > &C, Eigen::PlainObjectBase< DerivedF > &F, Eigen::PlainObjectBase< DerivedJ > &J)
 Given a polygon mesh, trivially triangulate each polygon with a fan.
 
void exactinit ()
 Initialize internal variable used by predciates.
 
template<typename Vector2D >
Orientation orient2d (const Eigen::MatrixBase< Vector2D > &pa, const Eigen::MatrixBase< Vector2D > &pb, const Eigen::MatrixBase< Vector2D > &pc)
 Compute the orientation of the triangle formed by pa, pb, pc.
 
template<typename Vector3D >
Orientation orient3d (const Eigen::MatrixBase< Vector3D > &pa, const Eigen::MatrixBase< Vector3D > &pb, const Eigen::MatrixBase< Vector3D > &pc, const Eigen::MatrixBase< Vector3D > &pd)
 Compute the orientation of the tetrahedron formed by pa, pb, pc, pd.
 
template<typename Vector2D >
Orientation incircle (const Eigen::MatrixBase< Vector2D > &pa, const Eigen::MatrixBase< Vector2D > &pb, const Eigen::MatrixBase< Vector2D > &pc, const Eigen::MatrixBase< Vector2D > &pd)
 Decide whether a point is inside/outside/on a circle.
 
template<typename Vector3D >
Orientation insphere (const Eigen::MatrixBase< Vector3D > &pa, const Eigen::MatrixBase< Vector3D > &pb, const Eigen::MatrixBase< Vector3D > &pc, const Eigen::MatrixBase< Vector3D > &pd, const Eigen::MatrixBase< Vector3D > &pe)
 Decide whether a point is inside/outside/on a sphere.
 
template<typename DerivedP >
bool segment_segment_intersect (const Eigen::MatrixBase< DerivedP > &A, const Eigen::MatrixBase< DerivedP > &B, const Eigen::MatrixBase< DerivedP > &C, const Eigen::MatrixBase< DerivedP > &D)
 Given two segments in 2d test whether they intersect each other using predicates orient2d.
 

Enumeration Type Documentation

◆ Orientation

enum class igl::predicates::Orientation
strong

Types of orientations and other predicate results.

include/igl/predicates/predicates.h

Enumerator
POSITIVE 
INSIDE 
NEGATIVE 
OUTSIDE 
COLLINEAR 
COPLANAR 
COCIRCULAR 
COSPHERICAL 
DEGENERATE 

Function Documentation

◆ delaunay_triangulation()

template<typename DerivedV , typename DerivedF >
void igl::predicates::delaunay_triangulation ( const Eigen::MatrixBase< DerivedV > &  V,
Eigen::PlainObjectBase< DerivedF > &  F 
)

Given a set of points in 2D, return a Delaunay triangulation of these points using predicates.

Parameters
[in]V#V by 2 list of vertex positions
[out]F#F by 3 of faces in Delaunay triangulation.

◆ ear_clipping() [1/2]

template<typename DerivedP , typename DerivedRT , typename DerivedF , typename DerivedI >
void igl::predicates::ear_clipping ( const Eigen::MatrixBase< DerivedP > &  P,
const Eigen::MatrixBase< DerivedRT > &  RT,
Eigen::PlainObjectBase< DerivedF > &  eF,
Eigen::PlainObjectBase< DerivedI > &  I 
)

Implementation of ear clipping triangulation algorithm for a 2D polygon.

https://www.geometrictools.com/Documentation/TriangulationByEarClipping.pdf If the polygon is simple and oriented counter-clockwise, all vertices will be clipped and the result mesh is (P,eF) Otherwise, the function will try to clip as many ears as possible.

Parameters
[in]P: n*2, size n 2D polygon
[in]RTn*1, preserved vertices (do not clip) marked as 1, otherwise 0
[out]eFclipped ears, in original index of P
[out]I: size #nP vector, maps index from nP to P, e.g. nP's ith vertex is origianlly I(i) in P
Precondition
To result in a proper mesh, P should be oriented counter-clockwise with no self-intersections.
Note
This implementation does not handle polygons with holes.

◆ ear_clipping() [2/2]

template<typename DerivedP , typename DerivedF >
bool igl::predicates::ear_clipping ( const Eigen::MatrixBase< DerivedP > &  P,
Eigen::PlainObjectBase< DerivedF > &  eF 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Reverses P if necessary. Orientation of output will match input

Returns
true if mesh is proper (should correspond to input being a simple polygon in either orientation), false otherwise.

◆ lexicographic_triangulation()

template<typename DerivedV , typename DerivedF >
void igl::predicates::lexicographic_triangulation ( const Eigen::MatrixBase< DerivedV > &  V,
Eigen::PlainObjectBase< DerivedF > &  F 
)

Given a set of points in 2D, return a lexicographic triangulation of these points using predicates.

Parameters
[in]V#V by 2 list of vertex positions
[out]F#F by 3 of faces in Delaunay triangulation.

◆ point_inside_convex_polygon()

template<typename DerivedP , typename DerivedQ >
bool igl::predicates::point_inside_convex_polygon ( const Eigen::MatrixBase< DerivedP > &  P,
const Eigen::MatrixBase< DerivedQ > &  q 
)

check whether 2d point lies inside 2d convex polygon

Parameters
[in]Pn*2 polygon, n >= 3
[in]q2d query point
Returns
true if point is inside polygon

◆ polygons_to_triangles()

template<typename DerivedV , typename DerivedI , typename DerivedC , typename DerivedF , typename DerivedJ >
void igl::predicates::polygons_to_triangles ( const Eigen::MatrixBase< DerivedV > &  V,
const Eigen::MatrixBase< DerivedI > &  I,
const Eigen::MatrixBase< DerivedC > &  C,
Eigen::PlainObjectBase< DerivedF > &  F,
Eigen::PlainObjectBase< DerivedJ > &  J 
)

Given a polygon mesh, trivially triangulate each polygon with a fan.

This purely combinatorial triangulation will work well for convex/flat polygons and degrade otherwise.

Parameters
[in]V#V by dim list of vertex positions
[in]I#I vectorized list of polygon corner indices into rows of some matrix V
[in]C#polygons+1 list of cumulative polygon sizes so that C(i+1)-C(i) = size of the ith polygon, and so I(C(i)) through I(C(i+1)-1) are the indices of the ith polygon
[out]F#F by 3 list of triangle indices into rows of V
[out]J#F list of indices into 0:#P-1 of corresponding polygon

◆ exactinit()

void igl::predicates::exactinit ( )

Initialize internal variable used by predciates.

Must be called before using exact predicates. It is safe to call this function from multiple threads.

include/igl/predicates/predicates.h

◆ orient2d()

template<typename Vector2D >
Orientation igl::predicates::orient2d ( const Eigen::MatrixBase< Vector2D > &  pa,
const Eigen::MatrixBase< Vector2D > &  pb,
const Eigen::MatrixBase< Vector2D > &  pc 
)

Compute the orientation of the triangle formed by pa, pb, pc.

Parameters
[in]pa2D point on line
[in]pb2D point on line
[in]pc2D query point.
Returns
POSITIVE if pa, pb, pc are counterclockwise oriented. NEGATIVE if they are clockwise oriented. COLLINEAR if they are collinear.

include/igl/predicates/predicates.h

◆ orient3d()

template<typename Vector3D >
Orientation igl::predicates::orient3d ( const Eigen::MatrixBase< Vector3D > &  pa,
const Eigen::MatrixBase< Vector3D > &  pb,
const Eigen::MatrixBase< Vector3D > &  pc,
const Eigen::MatrixBase< Vector3D > &  pd 
)

Compute the orientation of the tetrahedron formed by pa, pb, pc, pd.

Parameters
[in]pa2D point on plane
[in]pb2D point on plane
[in]pc2D point on plane
[in]pd2D query point
Returns
POSITIVE if pd is "below" the oriented plane formed by pa, pb and pc. NEGATIVE if pd is "above" the plane. COPLANAR if pd is on the plane.

include/igl/predicates/predicates.h

◆ incircle()

template<typename Vector2D >
Orientation igl::predicates::incircle ( const Eigen::MatrixBase< Vector2D > &  pa,
const Eigen::MatrixBase< Vector2D > &  pb,
const Eigen::MatrixBase< Vector2D > &  pc,
const Eigen::MatrixBase< Vector2D > &  pd 
)

Decide whether a point is inside/outside/on a circle.

Parameters
[in]pa2D point on circle
[in]pb2D point on circle
[in]pc2D point on circle
[in]pd2D point query
Returns
INSIDE if pd is inside of the circle defined by pa, pb and pc. OUSIDE if pd is outside of the circle. COCIRCULAR pd is exactly on the circle.

include/igl/predicates/predicates.h

◆ insphere()

template<typename Vector3D >
Orientation igl::predicates::insphere ( const Eigen::MatrixBase< Vector3D > &  pa,
const Eigen::MatrixBase< Vector3D > &  pb,
const Eigen::MatrixBase< Vector3D > &  pc,
const Eigen::MatrixBase< Vector3D > &  pd,
const Eigen::MatrixBase< Vector3D > &  pe 
)

Decide whether a point is inside/outside/on a sphere.

Parameters
[in]pa2D point on sphere
[in]pb2D point on sphere
[in]pc2D point on sphere
[in]pd2D point on sphere
[in]pe2D point query
Returns
INSIDE if pe is inside of the sphere defined by pa, pb, pc and pd. OUSIDE if pe is outside of the sphere. COSPHERICAL pd is exactly on the sphere.

include/igl/predicates/predicates.h

◆ segment_segment_intersect()

template<typename DerivedP >
bool igl::predicates::segment_segment_intersect ( const Eigen::MatrixBase< DerivedP > &  A,
const Eigen::MatrixBase< DerivedP > &  B,
const Eigen::MatrixBase< DerivedP > &  C,
const Eigen::MatrixBase< DerivedP > &  D 
)

Given two segments in 2d test whether they intersect each other using predicates orient2d.

Parameters
[in]A1st endpoint of segment 1
[in]B2st endpoint of segment 1
[in]C1st endpoint of segment 2
[in]D2st endpoint of segment 2
Returns
true if they intersect