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

Namespaces

namespace  scaf
 

Classes

struct  SCAFData
 Use a similar interface to igl::slim Implement ready-to-use 2D version of the algorithm described in SCAF: Simplicial Complex Augmentation Framework for Bijective Maps Zhongshi Jiang, Scott Schaefer, Daniele Panozzo, ACM Trancaction on Graphics (Proc. More...
 

Functions

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.
 

Function Documentation

◆ cdt()

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]flagsstring 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

◆ scaf_precompute()

void igl::triangle::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.

Parameters
[in]V#V by 3 list of mesh vertex positions
[in]F#F by 3/3 list of mesh faces (triangles/tets)
[in]V_init#V by 3 list of initial mesh vertex positions
[in,out]dataresulting precomputed data
[in]slim_energyEnergy type to minimize
[in]blist of boundary indices into V (soft constraint)
[in]bc#b by dim list of boundary conditions (soft constraint)
[in]soft_pSoft penalty factor (can be zero)
Note
Why aren't slim_energy, b, bc, soft_p const?

include/igl/triangle/scaf.h

◆ scaf_solve()

Eigen::MatrixXd igl::triangle::scaf_solve ( triangle::SCAFData data,
int  iter_num 
)

Run iter_num iterations of SCAF, with precomputed data.

Parameters
[in]dataprecomputed data
[in]iter_numnumber of iterations to run
Returns
resulting V_o (in SLIMData): #V by dim list of mesh vertex positions

◆ scaf_system()

void igl::triangle::scaf_system ( triangle::SCAFData s,
Eigen::SparseMatrix< double > &  L,
Eigen::VectorXd &  rhs 
)

Set up the SCAF system L * uv = rhs, without solving it.

Parameters
[in]sigl::SCAFData. Will be modified by energy and Jacobian computation.
[out]Lm by m matrix
[out]rhsm by 1 vector with m = dim * (#V_mesh + #V_scaf - #V_frame)

◆ triangulate() [1/2]

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]flagsstring 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

◆ triangulate() [2/2]

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.