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

Namespaces

namespace  cgal
 
namespace  comiso
 
namespace  tetgen
 

Functions

template<typename DerivedValues , typename DerivedPoints , typename DerivedVertices , typename DerivedFaces >
void marching_cubes (const Eigen::MatrixBase< DerivedValues > &values, const Eigen::MatrixBase< DerivedPoints > &points, const unsigned x_res, const unsigned y_res, const unsigned z_res, const double isovalue, Eigen::PlainObjectBase< DerivedVertices > &vertices, Eigen::PlainObjectBase< DerivedFaces > &faces)
 Performs marching cubes reconstruction on a grid defined by values, and points, and generates a mesh defined by vertices and faces.
 
template<typename DerivedValues , typename DerivedPoints , typename DerivedVertices , typename DerivedFaces >
void marching_cubes (const Eigen::MatrixBase< DerivedValues > &values, const Eigen::MatrixBase< DerivedPoints > &points, const unsigned x_res, const unsigned y_res, const unsigned z_res, Eigen::PlainObjectBase< DerivedVertices > &vertices, Eigen::PlainObjectBase< DerivedFaces > &faces)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts. Overload of the above function where the isovalue defaults to 0.0.
 
template<typename DerivedValue , typename DerivedPoint , typename DerivedPoints , typename DerivedVertices , typename DerivedFaces >
void marching_cubes (const std::function< DerivedValue(const DerivedPoint &) > &value_fun, const Eigen::MatrixBase< DerivedPoints > &points, const unsigned x_res, const unsigned y_res, const unsigned z_res, const double isovalue, Eigen::PlainObjectBase< DerivedVertices > &vertices, Eigen::PlainObjectBase< DerivedFaces > &faces)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
template<typename DerivedValues , typename DerivedPoints , typename DerivedVertices , typename DerivedIndices , typename DerivedFaces >
void marching_cubes (const Eigen::MatrixBase< DerivedValues > &values, const Eigen::MatrixBase< DerivedPoints > &points, const Eigen::MatrixBase< DerivedIndices > &indices, const double isovalue, Eigen::PlainObjectBase< DerivedVertices > &vertices, Eigen::PlainObjectBase< DerivedFaces > &faces)
 Perform marching cubes reconstruction on the sparse grid cells defined by (indices, points).
 
template<typename DerivedValues , typename DerivedPoints , typename DerivedVertices , typename DerivedIndices , typename DerivedFaces >
void marching_cubes (const Eigen::MatrixBase< DerivedValues > &values, const Eigen::MatrixBase< DerivedPoints > &points, const Eigen::MatrixBase< DerivedIndices > &indices, Eigen::PlainObjectBase< DerivedVertices > &vertices, Eigen::PlainObjectBase< DerivedFaces > &faces)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
bool progressive_hulls (const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, const size_t max_m, Eigen::MatrixXd &U, Eigen::MatrixXi &G, Eigen::VectorXi &J)
 Collapses edges until desired number of faces is achieved but ensures that new vertices are placed outside all previous meshes as per "progressive hulls" in "Silhouette clipping" [Sander et al.
 
void progressive_hulls_cost_and_placement (const int e, const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, const Eigen::MatrixXi &E, const Eigen::VectorXi &EMAP, const Eigen::MatrixXi &EF, const Eigen::MatrixXi &EI, double &cost, Eigen::RowVectorXd &p)
 A "cost and placement" compatible with igl::decimate implementing the "progressive hulls" algorithm in "Silhouette clipping" [Sander et al.
 
bool quadprog (const Eigen::MatrixXd &G, const Eigen::VectorXd &g0, const Eigen::MatrixXd &CE, const Eigen::VectorXd &ce0, const Eigen::MatrixXd &CI, const Eigen::VectorXd &ci0, Eigen::VectorXd &x)
 Solve a (dense) convex quadratric program.
 

Function Documentation

◆ marching_cubes() [1/5]

template<typename DerivedValues , typename DerivedPoints , typename DerivedVertices , typename DerivedFaces >
void igl::copyleft::marching_cubes ( const Eigen::MatrixBase< DerivedValues > &  values,
const Eigen::MatrixBase< DerivedPoints > &  points,
const unsigned  x_res,
const unsigned  y_res,
const unsigned  z_res,
const double  isovalue,
Eigen::PlainObjectBase< DerivedVertices > &  vertices,
Eigen::PlainObjectBase< DerivedFaces > &  faces 
)

Performs marching cubes reconstruction on a grid defined by values, and points, and generates a mesh defined by vertices and faces.

Parameters
[in]values#number_of_grid_points x 1 array – the scalar values of an implicit function defined on the grid points (<0 in the inside of the surface, 0 on the border, >0 outside)
[in]points#number_of_grid_points x 3 array – 3-D positions of the grid points, ordered in x,y,z order: points[index] = the point at (x,y,z) where : x = (index % (xres -1), y = (index / (xres-1)) %(yres-1), z = index / (xres -1) / (yres -1) ). where x,y,z index x, y, z dimensions i.e. index = x + y*xres + z*xres*yres
[in]xresresolutions of the grid in x dimension
[in]yresresolutions of the grid in y dimension
[in]zresresolutions of the grid in z dimension
[in]isovaluethe isovalue of the surface to reconstruct
[out]vertices#V by 3 list of mesh vertex positions
[out]faces#F by 3 list of mesh triangle indices
See also
igl::marching_cubes

◆ marching_cubes() [2/5]

template<typename DerivedValues , typename DerivedPoints , typename DerivedVertices , typename DerivedFaces >
void igl::copyleft::marching_cubes ( const Eigen::MatrixBase< DerivedValues > &  values,
const Eigen::MatrixBase< DerivedPoints > &  points,
const unsigned  x_res,
const unsigned  y_res,
const unsigned  z_res,
Eigen::PlainObjectBase< DerivedVertices > &  vertices,
Eigen::PlainObjectBase< DerivedFaces > &  faces 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts. Overload of the above function where the isovalue defaults to 0.0.

◆ marching_cubes() [3/5]

template<typename DerivedValue , typename DerivedPoint , typename DerivedPoints , typename DerivedVertices , typename DerivedFaces >
void igl::copyleft::marching_cubes ( const std::function< DerivedValue(const DerivedPoint &) > &  value_fun,
const Eigen::MatrixBase< DerivedPoints > &  points,
const unsigned  x_res,
const unsigned  y_res,
const unsigned  z_res,
const double  isovalue,
Eigen::PlainObjectBase< DerivedVertices > &  vertices,
Eigen::PlainObjectBase< DerivedFaces > &  faces 
)

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

Parameters
[in]value_funa function that takes a 3D point and returns a scalar value

◆ marching_cubes() [4/5]

template<typename DerivedValues , typename DerivedPoints , typename DerivedVertices , typename DerivedIndices , typename DerivedFaces >
void igl::copyleft::marching_cubes ( const Eigen::MatrixBase< DerivedValues > &  values,
const Eigen::MatrixBase< DerivedPoints > &  points,
const Eigen::MatrixBase< DerivedIndices > &  indices,
const double  isovalue,
Eigen::PlainObjectBase< DerivedVertices > &  vertices,
Eigen::PlainObjectBase< DerivedFaces > &  faces 
)

Perform marching cubes reconstruction on the sparse grid cells defined by (indices, points).

The indices parameter is an nx8 dense array of index values into the points and values arrays. Each row of indices represents a cube for which to generate vertices and faces over.

Parameters
[in]values#number_of_grid_points x 1 array – the scalar values of an implicit function defined on the grid points (<0 in the inside of the surface, 0 on the border, >0 outside)
[in]points#number_of_grid_points x 3 array – 3-D positions of the grid points, ordered in x,y,z order:
[in]indices#cubes x 8 array – one row for each cube where each value is the index of a vertex in points and a scalar in values. i.e. points[indices[i, j]] = the position of the j'th vertex of the i'th cube
[out]vertices#V by 3 list of mesh vertex positions
[out]faces#F by 3 list of mesh triangle indices
Note
The winding direction of the cube indices will affect the output winding of the faces

◆ marching_cubes() [5/5]

template<typename DerivedValues , typename DerivedPoints , typename DerivedVertices , typename DerivedIndices , typename DerivedFaces >
void igl::copyleft::marching_cubes ( const Eigen::MatrixBase< DerivedValues > &  values,
const Eigen::MatrixBase< DerivedPoints > &  points,
const Eigen::MatrixBase< DerivedIndices > &  indices,
Eigen::PlainObjectBase< DerivedVertices > &  vertices,
Eigen::PlainObjectBase< DerivedFaces > &  faces 
)

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

isovalue defaults to 0.0

◆ progressive_hulls()

bool igl::copyleft::progressive_hulls ( const Eigen::MatrixXd &  V,
const Eigen::MatrixXi &  F,
const size_t  max_m,
Eigen::MatrixXd &  U,
Eigen::MatrixXi &  G,
Eigen::VectorXi &  J 
)

Collapses edges until desired number of faces is achieved but ensures that new vertices are placed outside all previous meshes as per "progressive hulls" in "Silhouette clipping" [Sander et al.

2000].

Precondition
Assumes (V,F) is a closed manifold mesh
Parameters
[in]V#V by dim list of vertex positions
[in]F#F by 3 list of face indices into V.
[in]max_mdesired number of output faces
[out]U#U by dim list of output vertex posistions (can be same ref as V)
[out]G#G by 3 list of output face indices into U (can be same ref as G)
[out]J#G list of indices into F of birth faces
Returns
true if m was reached (otherwise #G > m)

◆ progressive_hulls_cost_and_placement()

void igl::copyleft::progressive_hulls_cost_and_placement ( const int  e,
const Eigen::MatrixXd &  V,
const Eigen::MatrixXi &  F,
const Eigen::MatrixXi &  E,
const Eigen::VectorXi &  EMAP,
const Eigen::MatrixXi &  EF,
const Eigen::MatrixXi &  EI,
double &  cost,
Eigen::RowVectorXd &  p 
)

A "cost and placement" compatible with igl::decimate implementing the "progressive hulls" algorithm in "Silhouette clipping" [Sander et al.

2000]. This implementation fixes an issue that the original linear program becomes unstable for flat patches by introducing a small quadratic energy term pulling the collapsed edge toward its midpoint. This function is not really meant to be called directly but rather passed to igl::decimate as a handle.

Parameters
[in]eindex of edge to be collapsed
[in]V#V by 3 list of vertex positions
[in]F#F by 3 list of faces indices into V
[in]E#E by 3 list of edges indices into V
[in]EMAP#F*3 list of indices into E, mapping each directed edge to unique unique edge in E
[in]EF#E by 2 list of edge flaps, EF(e,0)=f means e=(i-->j) is the edge of F(f,:) opposite the vth corner, where EI(e,0)=v. Similarly EF(e,1) " e=(j->i)
[in]EI#E by 2 list of edge flap corners (see above).
[out]costcost of collapsing edge e
[out]pposition to place collapsed vertex
See also
igl::decimate, igl::collapse_edge

◆ quadprog()

bool igl::copyleft::quadprog ( const Eigen::MatrixXd &  G,
const Eigen::VectorXd &  g0,
const Eigen::MatrixXd &  CE,
const Eigen::VectorXd &  ce0,
const Eigen::MatrixXd &  CI,
const Eigen::VectorXd &  ci0,
Eigen::VectorXd &  x 
)

Solve a (dense) convex quadratric program.

Given in the form

 min  0.5 x G x + g0 x
 s.t. CE' x + ce0  = 0
 and  CI' x + ci0 >= 0
Parameters
[in]G#x by #x matrix of quadratic coefficients
[in]g0#x vector of linear coefficients
[in]CE#x by #CE list of linear equality coefficients
[in]ce0#CE list of linear equality right-hand sides
[in]CI#x by #CI list of linear equality coefficients
[in]ci0#CI list of linear equality right-hand sides
[out]x#x vector of solution values
Returns
true iff success