| 
    libigl v2.5.0
    
   | 
 
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.   | |
      
  | 
  strong | 
Types of orientations and other predicate results.
include/igl/predicates/predicates.h
| Enumerator | |
|---|---|
| POSITIVE | |
| INSIDE | |
| NEGATIVE | |
| OUTSIDE | |
| COLLINEAR | |
| COPLANAR | |
| COCIRCULAR | |
| COSPHERICAL | |
| DEGENERATE | |
| 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.
| [in] | V | #V by 2 list of vertex positions | 
| [out] | F | #F by 3 of faces in Delaunay triangulation. | 
| 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.
| [in] | P | : n*2, size n 2D polygon | 
| [in] | RT | n*1, preserved vertices (do not clip) marked as 1, otherwise 0 | 
| [out] | eF | clipped 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 | 
| 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
| 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.
| [in] | V | #V by 2 list of vertex positions | 
| [out] | F | #F by 3 of faces in Delaunay triangulation. | 
| 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
| [in] | P | n*2 polygon, n >= 3 | 
| [in] | q | 2d query point | 
| 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.
| [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 | 
| 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.
| 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.
| [in] | pa | 2D point on line | 
| [in] | pb | 2D point on line | 
| [in] | pc | 2D query point. | 
| 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.
| [in] | pa | 2D point on plane | 
| [in] | pb | 2D point on plane | 
| [in] | pc | 2D point on plane | 
| [in] | pd | 2D query point | 
| 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.
| [in] | pa | 2D point on circle | 
| [in] | pb | 2D point on circle | 
| [in] | pc | 2D point on circle | 
| [in] | pd | 2D point query | 
| 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.
| [in] | pa | 2D point on sphere | 
| [in] | pb | 2D point on sphere | 
| [in] | pc | 2D point on sphere | 
| [in] | pd | 2D point on sphere | 
| [in] | pe | 2D point query | 
| 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.
| [in] | A | 1st endpoint of segment 1 | 
| [in] | B | 2st endpoint of segment 1 | 
| [in] | C | 1st endpoint of segment 2 | 
| [in] | D | 2st endpoint of segment 2 |