|
| template<typename DerivedV , typename DerivedE , typename DerivedWV , typename DerivedWF , typename DerivedWE , typename DerivedJ > |
| void | cdt (const Eigen::MatrixBase< DerivedV > &V, const Eigen::MatrixBase< DerivedE > &E, const std::string &flags, Eigen::PlainObjectBase< DerivedWV > &WV, Eigen::PlainObjectBase< DerivedWF > &WF, Eigen::PlainObjectBase< DerivedWE > &WE, Eigen::PlainObjectBase< DerivedJ > &J) |
| | Construct the constrained delaunay triangulation of the convex hull of a given set of points and segments in 2D.
|
| |
| void | scaf_precompute (const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, const Eigen::MatrixXd &V_init, triangle::SCAFData &data, MappingEnergyType slim_energy, Eigen::VectorXi &b, Eigen::MatrixXd &bc, double soft_p) |
| | Compute necessary information to start using SCAF.
|
| |
| Eigen::MatrixXd | scaf_solve (triangle::SCAFData &data, int iter_num) |
| | Run iter_num iterations of SCAF, with precomputed data.
|
| |
| void | scaf_system (triangle::SCAFData &s, Eigen::SparseMatrix< double > &L, Eigen::VectorXd &rhs) |
| | Set up the SCAF system L * uv = rhs, without solving it.
|
| |
| template<typename DerivedV , typename DerivedE , typename DerivedH , typename DerivedVM , typename DerivedEM , typename DerivedV2 , typename DerivedF2 , typename DerivedVM2 , typename DerivedE2 , typename DerivedEM2 > |
| void | triangulate (const Eigen::MatrixBase< DerivedV > &V, const Eigen::MatrixBase< DerivedE > &E, const Eigen::MatrixBase< DerivedH > &H, const Eigen::MatrixBase< DerivedVM > &VM, const Eigen::MatrixBase< DerivedEM > &EM, const std::string flags, Eigen::PlainObjectBase< DerivedV2 > &V2, Eigen::PlainObjectBase< DerivedF2 > &F2, Eigen::PlainObjectBase< DerivedVM2 > &VM2, Eigen::PlainObjectBase< DerivedE2 > &E2, Eigen::PlainObjectBase< DerivedEM2 > &EM2) |
| | Triangulate the interior of a polygon using the triangle library.
|
| |
| template<typename DerivedV , typename DerivedE , typename DerivedH , typename DerivedV2 , typename DerivedF2 > |
| void | triangulate (const Eigen::MatrixBase< DerivedV > &V, const Eigen::MatrixBase< DerivedE > &E, const Eigen::MatrixBase< DerivedH > &H, const std::string flags, Eigen::PlainObjectBase< DerivedV2 > &V2, Eigen::PlainObjectBase< DerivedF2 > &F2) |
| | 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 DerivedE , typename DerivedWV , typename DerivedWF , typename DerivedWE , typename DerivedJ >
| void igl::triangle::cdt |
( |
const Eigen::MatrixBase< DerivedV > & |
V, |
|
|
const Eigen::MatrixBase< DerivedE > & |
E, |
|
|
const std::string & |
flags, |
|
|
Eigen::PlainObjectBase< DerivedWV > & |
WV, |
|
|
Eigen::PlainObjectBase< DerivedWF > & |
WF, |
|
|
Eigen::PlainObjectBase< DerivedWE > & |
WE, |
|
|
Eigen::PlainObjectBase< DerivedJ > & |
J |
|
) |
| |
Construct the constrained delaunay triangulation of the convex hull of a given set of points and segments in 2D.
This differs from a direct call to triangulate because it will preprocess the input to remove duplicates and return an adjusted segment list on the output.
- Parameters
-
| [in] | V | #V by 2 list of texture mesh vertices |
| [in] | E | #E by 2 list of constraint edge indices into V |
| [in] | flags | string of triangle flags should contain "-c" unless the some subset of segments are known to enclose all other points/segments. |
| [out] | WV | #WV by 2 list of background mesh vertices |
| [out] | WF | #WF by 2 list of background mesh triangle indices into WV |
| [out] | WE | #WE by 2 list of constraint edge indices into WV (might be smaller than E because degenerate constraints have been removed) |
| [out] | J | #V list of indices into WF/WE for each vertex in V |
template<typename DerivedV , typename DerivedE , typename DerivedH , typename DerivedVM , typename DerivedEM , typename DerivedV2 , typename DerivedF2 , typename DerivedVM2 , typename DerivedE2 , typename DerivedEM2 >
| void igl::triangle::triangulate |
( |
const Eigen::MatrixBase< DerivedV > & |
V, |
|
|
const Eigen::MatrixBase< DerivedE > & |
E, |
|
|
const Eigen::MatrixBase< DerivedH > & |
H, |
|
|
const Eigen::MatrixBase< DerivedVM > & |
VM, |
|
|
const Eigen::MatrixBase< DerivedEM > & |
EM, |
|
|
const std::string |
flags, |
|
|
Eigen::PlainObjectBase< DerivedV2 > & |
V2, |
|
|
Eigen::PlainObjectBase< DerivedF2 > & |
F2, |
|
|
Eigen::PlainObjectBase< DerivedVM2 > & |
VM2, |
|
|
Eigen::PlainObjectBase< DerivedE2 > & |
E2, |
|
|
Eigen::PlainObjectBase< DerivedEM2 > & |
EM2 |
|
) |
| |
Triangulate the interior of a polygon using the triangle library.
- Parameters
-
| [in] | V | #V by 2 list of 2D vertex positions |
| [in] | E | #E by 2 list of vertex ids forming unoriented edges of the boundary of the polygon |
| [in] | H | #H by 2 coordinates of points contained inside holes of the polygon |
| [in] | VM | #V list of markers for input vertices |
| [in] | EM | #E list of markers for input edges |
| [in] | flags | string of options pass to triangle (see triangle documentation) |
| [out] | V2 | #V2 by 2 coordinates of the vertives of the generated triangulation |
| [out] | F2 | #F2 by 3 list of indices forming the faces of the generated triangulation |
| [out] | VM2 | #V2 list of markers for output vertices |
| [out] | E2 | #E2 by 2 list of output edges |
| [out] | EM2 | #E2 list of markers for output edges |
template<typename DerivedV , typename DerivedE , typename DerivedH , typename DerivedV2 , typename DerivedF2 >
| void igl::triangle::triangulate |
( |
const Eigen::MatrixBase< DerivedV > & |
V, |
|
|
const Eigen::MatrixBase< DerivedE > & |
E, |
|
|
const Eigen::MatrixBase< DerivedH > & |
H, |
|
|
const std::string |
flags, |
|
|
Eigen::PlainObjectBase< DerivedV2 > & |
V2, |
|
|
Eigen::PlainObjectBase< DerivedF2 > & |
F2 |
|
) |
| |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.