Mesh processing algorithms. More...
Enumerations | |
enum class | Curvature { min , max , mean , gauss , max_abs } |
Type of curvature to be computed. More... | |
Functions | |
void | curvature (SurfaceMesh &mesh, Curvature c=Curvature::mean, int smoothing_steps=0, bool use_tensor=false, bool use_two_ring=false) |
Compute per-vertex curvature (min,max,mean,Gaussian). | |
void | decimate (SurfaceMesh &mesh, unsigned int n_vertices, Scalar aspect_ratio=0.0, Scalar edge_length=0.0, unsigned int max_valence=0, Scalar normal_deviation=0.0, Scalar hausdorff_error=0.0, Scalar seam_threshold=1e-2, Scalar seam_angle_deviation=1) |
Mesh decimation based on approximation error and fairness criteria. | |
Scalar | clamp_cot (const Scalar v) |
Clamp cotangent values as if angles are in [3, 177]. | |
Scalar | clamp_cos (const Scalar v) |
Clamp cosine values as if angles are in [3, 177]. | |
Scalar | angle (const Point &v0, const Point &v1) |
Compute the angle between two (un-normalized) vectors. | |
Scalar | sin (const Point &v0, const Point &v1) |
Compute the sine of angle between two (un-normalized) vectors. | |
Scalar | cos (const Point &v0, const Point &v1) |
Compute the cosine of angle between two (un-normalized) vectors. | |
Scalar | cotan (const Point &v0, const Point &v1) |
Compute the cotangent of angle between two (un-normalized) vectors. | |
Scalar | triangle_area (const Point &p0, const Point &p1, const Point &p2) |
Compute the area of a triangle given by three points. | |
Scalar | face_area (const SurfaceMesh &mesh, Face f) |
Compute area of face f . | |
Scalar | surface_area (const SurfaceMesh &mesh) |
Compute the surface area of mesh as the sum of face areas. | |
Scalar | voronoi_area (const SurfaceMesh &mesh, Vertex v) |
Compute the (barycentric) Voronoi area of vertex v . | |
Scalar | voronoi_area_mixed (const SurfaceMesh &mesh, Vertex v) |
Compute mixed Voronoi area of a vertex. | |
Scalar | edge_area (const SurfaceMesh &mesh, Edge e) |
Compute the area assigned to edge e . | |
Scalar | volume (const SurfaceMesh &mesh) |
Compute the volume of a mesh. | |
Point | centroid (const SurfaceMesh &mesh, Face f) |
Compute the barycenter/centroid of face f . | |
Point | centroid (const SurfaceMesh &mesh) |
Compute the barycenter (centroid) of the mesh . | |
void | dual (SurfaceMesh &mesh) |
Compute dual of a mesh. | |
double | cotan_weight (const SurfaceMesh &mesh, Edge e) |
Compute the cotangent weight for edge e . | |
Point | laplace (const SurfaceMesh &mesh, Vertex v) |
Compute the Laplace vector for vertex v , normalized by Voronoi area. | |
Scalar | dist_point_line_segment (const Point &p, const Point &v0, const Point &v1, Point &nearest_point) |
Compute the distance of a point p to a line segment given by points v0 and v1 . | |
Scalar | dist_point_triangle (const Point &p, const Point &v0, const Point &v1, const Point &v2, Point &nearest_point) |
Compute the distance of a point p to the triangle given by points v0 , v1 , v2 . | |
void | minimize_area (SurfaceMesh &mesh) |
Minimize surface area. | |
void | minimize_curvature (SurfaceMesh &mesh) |
Minimize surface curvature. | |
void | fair (SurfaceMesh &mesh, unsigned int k=2) |
Implicit surface fairing. | |
size_t | detect_features (SurfaceMesh &mesh, Scalar angle) |
Mark edges with dihedral angle larger than angle as feature. | |
size_t | detect_boundary (SurfaceMesh &mesh) |
Mark all boundary edges as features. | |
void | clear_features (SurfaceMesh &mesh) |
Clear feature and boundary edges. | |
unsigned int | geodesics (SurfaceMesh &mesh, const std::vector< Vertex > &seeds, Scalar maxdist=std::numeric_limits< Scalar >::max(), unsigned int maxnum=std::numeric_limits< unsigned int >::max(), std::vector< Vertex > *neighbors=nullptr) |
Compute geodesic distance from a set of seed vertices. | |
void | geodesics_heat (SurfaceMesh &mesh, const std::vector< Vertex > &seeds) |
Compute geodesic distance from a set of seed vertices. | |
void | fill_hole (SurfaceMesh &mesh, Halfedge h) |
Fill the hole specified by halfedge h . | |
void | uniform_mass_matrix (const SurfaceMesh &mesh, DiagonalMatrix &M) |
Construct the mass matrix for the uniform Laplacian. | |
void | uniform_laplace_matrix (const SurfaceMesh &mesh, SparseMatrix &L) |
Construct the uniform Laplace matrix. | |
void | mass_matrix (const SurfaceMesh &mesh, DiagonalMatrix &M) |
Construct the (lumped) mass matrix for the cotangent Laplacian. | |
void | laplace_matrix (const SurfaceMesh &mesh, SparseMatrix &L, bool clamp=false) |
Construct the cotan Laplace matrix. | |
void | gradient_matrix (const SurfaceMesh &mesh, SparseMatrix &G) |
Construct the cotan gradient matrix. | |
void | divergence_matrix (const SurfaceMesh &mesh, SparseMatrix &D) |
Construct the cotan divergence matrix. | |
void | coordinates_to_matrix (const SurfaceMesh &mesh, DenseMatrix &X) |
For a mesh with N vertices, construct an Nx3 matrix containing the vertex coordinates in its rows. | |
void | matrix_to_coordinates (const DenseMatrix &X, SurfaceMesh &mesh) |
For a mesh with N vertices, set the vertex coordinates from the rows of an Nx3 matrix. | |
void | vertex_normals (SurfaceMesh &mesh) |
Compute vertex normals for the whole mesh . | |
void | face_normals (SurfaceMesh &mesh) |
Compute face normals for the whole mesh . | |
Normal | vertex_normal (const SurfaceMesh &mesh, Vertex v) |
Compute the normal vector of vertex v . | |
Normal | face_normal (const SurfaceMesh &mesh, Face f) |
Compute the normal vector of face f . | |
Normal | corner_normal (const SurfaceMesh &mesh, Halfedge h, Scalar crease_angle) |
Compute the normal vector of the polygon corner specified by the target vertex of halfedge h . | |
void | harmonic_parameterization (SurfaceMesh &mesh, bool use_uniform_weights=false) |
Compute discrete harmonic parameterization. | |
void | lscm_parameterization (SurfaceMesh &mesh) |
Compute parameterization based on least squares conformal mapping. | |
void | uniform_remeshing (SurfaceMesh &mesh, Scalar edge_length, unsigned int iterations=10, bool use_projection=true) |
Perform uniform remeshing. | |
void | adaptive_remeshing (SurfaceMesh &mesh, Scalar min_edge_length, Scalar max_edge_length, Scalar approx_error, unsigned int iterations=10, bool use_projection=true) |
Perform adaptive remeshing. | |
SurfaceMesh | tetrahedron () |
Generate tetrahedron. | |
SurfaceMesh | hexahedron () |
Generate hexahedron. | |
SurfaceMesh | octahedron () |
Generate octahedron. | |
SurfaceMesh | dodecahedron () |
Generate dodecahedron. | |
SurfaceMesh | icosahedron () |
Generate icosahedron. | |
SurfaceMesh | icosphere (size_t n_subdivisions=3) |
Generate icosphere refined by n_subdivisions . | |
SurfaceMesh | quad_sphere (size_t n_subdivisions=3) |
Generate quad sphere refined by n_subdivisions . | |
SurfaceMesh | uv_sphere (const Point ¢er=Point(0, 0, 0), Scalar radius=1.0, size_t n_slices=15, size_t n_stacks=15) |
Generate UV sphere with given center , radius , n_slices , and n_stacks . | |
SurfaceMesh | plane (size_t resolution=4) |
Generate a plane mesh. | |
SurfaceMesh | cone (size_t n_subdivisions=30, Scalar radius=1.0, Scalar height=2.5) |
Generate a cone mesh. | |
SurfaceMesh | cylinder (size_t n_subdivisions=30, Scalar radius=1.0, Scalar height=2.5) |
Generate a cylinder mesh. | |
SurfaceMesh | torus (size_t radial_resolution=20, size_t tubular_resolution=40, Scalar radius=1.0, Scalar thickness=0.4) |
Generate a torus mesh. | |
void | explicit_smoothing (SurfaceMesh &mesh, unsigned int iterations=10, bool use_uniform_laplace=false) |
Perform explicit Laplacian smoothing. | |
void | implicit_smoothing (SurfaceMesh &mesh, Scalar timestep=0.001, unsigned int iterations=1, bool use_uniform_laplace=false, bool rescale=true) |
Perform implicit Laplacian smoothing. | |
void | catmull_clark_subdivision (SurfaceMesh &mesh, BoundaryHandling boundary_handling=BoundaryHandling::Interpolate) |
Perform one step of Catmull-Clark subdivision. | |
void | loop_subdivision (SurfaceMesh &mesh, BoundaryHandling boundary_handling=BoundaryHandling::Interpolate) |
Perform one step of Loop subdivision. | |
void | quad_tri_subdivision (SurfaceMesh &mesh, BoundaryHandling boundary_handling=BoundaryHandling::Interpolate) |
Perform one step of quad-tri subdivision. | |
void | linear_subdivision (SurfaceMesh &mesh) |
Perform one step of linear quad-tri subdivision. | |
void | triangulate (SurfaceMesh &mesh) |
Triangulate all faces in mesh by applying triangulate(). | |
void | triangulate (SurfaceMesh &mesh, Face f) |
Triangulate the Face f . | |
BoundingBox | bounds (const SurfaceMesh &mesh) |
Compute the bounding box of mesh . | |
void | flip_faces (SurfaceMesh &mesh) |
Flip the orientation of all faces in mesh . | |
Scalar | min_face_area (const SurfaceMesh &mesh) |
Compute the minimum area of all faces in mesh . | |
Scalar | edge_length (const SurfaceMesh &mesh, Edge e) |
Compute length of an edge e in mesh . | |
Scalar | mean_edge_length (const SurfaceMesh &mesh) |
Compute mean edge length of mesh . | |
Mesh processing algorithms.
|
strong |
void adaptive_remeshing | ( | SurfaceMesh & | mesh, |
Scalar | min_edge_length, | ||
Scalar | max_edge_length, | ||
Scalar | approx_error, | ||
unsigned int | iterations = 10 , |
||
bool | use_projection = true |
||
) |
Perform adaptive remeshing.
Performs incremental remeshing based on edge collapse, split, flip, and tangential relaxation. See [2] and [10] for details.
mesh | The input mesh, modified in place. |
min_edge_length | The minimum edge length. |
max_edge_length | The maximum edge length. |
approx_error | The maximum approximation error. |
iterations | The number of iterations. |
use_projection | Use back-projection to the input surface. |
InvalidInputException | if the input precondition is violated. |
void catmull_clark_subdivision | ( | SurfaceMesh & | mesh, |
BoundaryHandling | boundary_handling = BoundaryHandling::Interpolate |
||
) |
Perform one step of Catmull-Clark subdivision.
See [5] for details.
mesh | The input mesh, modified in place. |
boundary_handling | Specify to interpolate or preserve boundary edges. |
Point centroid | ( | const SurfaceMesh & | mesh | ) |
Compute the barycenter (centroid) of the mesh
.
Computed as area-weighted mean of vertices.
void clear_features | ( | SurfaceMesh & | mesh | ) |
Clear feature and boundary edges.
Sets all "e:feature"
and "v:feature"
properties to false
.
SurfaceMesh cone | ( | size_t | n_subdivisions = 30 , |
Scalar | radius = 1.0 , |
||
Scalar | height = 2.5 |
||
) |
Generate a cone mesh.
Generates a polygonal mesh of a cone. The circular base lies in the x-y-plane and the tip points in positive z-direction.
n_subdivisions | Number of subdivisions of the base circle. Needs to be >= 3. Default: 30. |
radius | Radius of the base circle. Default: 1. |
height | Height of the the cone. Default: 2.5. |
void coordinates_to_matrix | ( | const SurfaceMesh & | mesh, |
DenseMatrix & | X | ||
) |
For a mesh with N vertices, construct an Nx3 matrix containing the vertex coordinates in its rows.
mesh | The input mesh. |
X | The output matrix. |
Normal corner_normal | ( | const SurfaceMesh & | mesh, |
Halfedge | h, | ||
Scalar | crease_angle | ||
) |
Compute the normal vector of the polygon corner specified by the target vertex of halfedge h
.
Averages incident corner normals if they are within crease_angle of the face normal. crease_angle
is in radians, not degrees.
double cotan_weight | ( | const SurfaceMesh & | mesh, |
Edge | e | ||
) |
Compute the cotangent weight for edge e
.
void curvature | ( | SurfaceMesh & | mesh, |
Curvature | c = Curvature::mean , |
||
int | smoothing_steps = 0 , |
||
bool | use_tensor = false , |
||
bool | use_two_ring = false |
||
) |
SurfaceMesh cylinder | ( | size_t | n_subdivisions = 30 , |
Scalar | radius = 1.0 , |
||
Scalar | height = 2.5 |
||
) |
Generate a cylinder mesh.
Generates a polygonal mesh of a cylinder. The cylinder is oriented in z-direction.
n_subdivisions | Number of subdivisions of the cylinder. Needs to be >= 3. Default: 30. |
radius | Radius of the cylinder. Default: 1. |
height | Height of the cylinder. Default: 2.5. |
void decimate | ( | SurfaceMesh & | mesh, |
unsigned int | n_vertices, | ||
Scalar | aspect_ratio = 0.0 , |
||
Scalar | edge_length = 0.0 , |
||
unsigned int | max_valence = 0 , |
||
Scalar | normal_deviation = 0.0 , |
||
Scalar | hausdorff_error = 0.0 , |
||
Scalar | seam_threshold = 1e-2 , |
||
Scalar | seam_angle_deviation = 1 |
||
) |
Mesh decimation based on approximation error and fairness criteria.
Performs incremental greedy mesh decimation based on halfedge collapses. See [14] and [11] for details.
mesh | Target mesh. Modified in place. |
n_vertices | Target number of vertices. |
aspect_ratio | Minimum aspect ratio of the triangles. |
edge_length | Minimum target edge length. |
max_valence | Maximum number of incident edges per vertex. |
normal_deviation | Maximum deviation of face normals. |
hausdorff_error | Maximum deviation from the original surface. |
seam_threshold | Threshold for texture seams. |
seam_angle_deviation | Maximum texture seam deviation. |
InvalidInputException | if the input precondition is violated. |
size_t detect_boundary | ( | SurfaceMesh & | mesh | ) |
Mark all boundary edges as features.
size_t detect_features | ( | SurfaceMesh & | mesh, |
Scalar | angle | ||
) |
Mark edges with dihedral angle larger than angle
as feature.
void divergence_matrix | ( | const SurfaceMesh & | mesh, |
SparseMatrix & | D | ||
) |
Construct the cotan divergence matrix.
Matrix is sparse and maps constant gradient vectors at non-boundary halfedges to values at vertices. The discrete operators are consistent, such that Laplacian is divergence of gradient. See [18] for details on triangle meshes and [4] for details on polygon meshes.
mesh | The input mesh. |
D | The output matrix. |
void dual | ( | SurfaceMesh & | mesh | ) |
Compute dual of a mesh.
Scalar edge_area | ( | const SurfaceMesh & | mesh, |
Edge | e | ||
) |
Compute the area assigned to edge e
.
A face with n edges assigns 1/n of its area to each edge.
void explicit_smoothing | ( | SurfaceMesh & | mesh, |
unsigned int | iterations = 10 , |
||
bool | use_uniform_laplace = false |
||
) |
Perform explicit Laplacian smoothing.
See [8] for details
mesh | The input mesh, modified in place. |
iterations | The number of iterations performed. |
use_uniform_laplace | Use uniform or cotan Laplacian. Default: cotan. |
Scalar face_area | ( | const SurfaceMesh & | mesh, |
Face | f | ||
) |
Compute area of face f
.
Computes standard area for triangles and norm of vector area for other polygons.
Normal face_normal | ( | const SurfaceMesh & | mesh, |
Face | f | ||
) |
Compute the normal vector of face f
.
Normal is computed as (normalized) sum of per-corner cross products of the two incident edges. This corresponds to the normalized vector area in [1]
void face_normals | ( | SurfaceMesh & | mesh | ) |
Compute face normals for the whole mesh
.
Calls face_normal() for each face and adds a new face property of type Normal named "f:normal".
void fair | ( | SurfaceMesh & | mesh, |
unsigned int | k = 2 |
||
) |
Implicit surface fairing.
Computes a surface by solving k-harmonic equation. See also [8] .
SolverException | in case of failure to solve the linear system |
InvalidInputException | in case of missing boundary constraints |
void fill_hole | ( | SurfaceMesh & | mesh, |
Halfedge | h | ||
) |
Fill the hole specified by halfedge h
.
Close simple holes (boundary loops of manifold vertices) by first filling the hole with an angle/area-minimizing triangulation, followed by isometric remeshing, and finished by curvature-minimizing fairing of the filled-in patch. See [15] for details.
InvalidInputException | in case on of the input preconditions is violated |
unsigned int geodesics | ( | SurfaceMesh & | mesh, |
const std::vector< Vertex > & | seeds, | ||
Scalar | maxdist = std::numeric_limits< Scalar >::max() , |
||
unsigned int | maxnum = std::numeric_limits< unsigned int >::max() , |
||
std::vector< Vertex > * | neighbors = nullptr |
||
) |
Compute geodesic distance from a set of seed vertices.
The method works by a Dijkstra-like breadth first traversal from the seed vertices, implemented by a heap structure. See [13] for details.
mesh | The input mesh, modified in place. | |
[in] | seeds | The vector of seed vertices. |
[in] | maxdist | The maximum distance up to which to compute the geodesic distances. |
[in] | maxnum | The maximum number of neighbors up to which to compute the geodesic distances. |
[out] | neighbors | The vector of neighbor vertices. |
void geodesics_heat | ( | SurfaceMesh & | mesh, |
const std::vector< Vertex > & | seeds | ||
) |
Compute geodesic distance from a set of seed vertices.
Compute geodesic distances based on the heat method, by solving two Poisson systems. Works on general polygon meshes. See [7] for details.
mesh | The input mesh, modified in place. |
seeds | The vector of seed vertices. |
void gradient_matrix | ( | const SurfaceMesh & | mesh, |
SparseMatrix & | G | ||
) |
Construct the cotan gradient matrix.
Matrix is sparse and maps values at vertices to constant gradient 3D-vectors at non-boundary halfedges. The discrete operators are consistent, such that Laplacian is divergence of gradient. See [18] for details on triangle meshes and [4] for details on polygon meshes.
mesh | The input mesh. |
G | The output matrix. |
void harmonic_parameterization | ( | SurfaceMesh & | mesh, |
bool | use_uniform_weights = false |
||
) |
Compute discrete harmonic parameterization.
See [9] for details.
InvalidInputException | if the input precondition is violated. |
SolverException | in case of failure to solve the linear system. |
SurfaceMesh icosphere | ( | size_t | n_subdivisions = 3 | ) |
Generate icosphere refined by n_subdivisions
.
Uses Loop subdivision to refine the initial icosahedron.
void implicit_smoothing | ( | SurfaceMesh & | mesh, |
Scalar | timestep = 0.001 , |
||
unsigned int | iterations = 1 , |
||
bool | use_uniform_laplace = false , |
||
bool | rescale = true |
||
) |
Perform implicit Laplacian smoothing.
mesh | The input mesh, modified in place. |
timestep | The time step taken. |
iterations | The number of iterations performed. |
use_uniform_laplace | Use uniform or cotan Laplacian. Default: cotan. |
rescale | Re-center and re-scale model after smoothing. Default: true. |
SolverException | in case of a failure to solve the linear system. |
Point laplace | ( | const SurfaceMesh & | mesh, |
Vertex | v | ||
) |
Compute the Laplace vector for vertex v
, normalized by Voronoi area.
void laplace_matrix | ( | const SurfaceMesh & | mesh, |
SparseMatrix & | L, | ||
bool | clamp = false |
||
) |
Construct the cotan Laplace matrix.
Matrix is sparse, symmetric and negative semi-definite. M(i,i) is the negative valence of vertex i. M(i,j) is cotangent weight of edge (i,j). M(i,i) is negative sum of off-diagonals. The discrete operators are consistent, such that Laplacian is divergence of gradient. See [18] for details on triangle meshes and [4] for details on polygon meshes.
mesh | The input mesh. |
L | The output matrix. |
clamp | Whether or not negative off-diagonal entries should be clamped to zero. |
void linear_subdivision | ( | SurfaceMesh & | mesh | ) |
Perform one step of linear quad-tri subdivision.
Suitable for mixed quad/triangle meshes.
mesh | The input mesh, modified in place. |
void loop_subdivision | ( | SurfaceMesh & | mesh, |
BoundaryHandling | boundary_handling = BoundaryHandling::Interpolate |
||
) |
Perform one step of Loop subdivision.
See [16] for details.
mesh | The input mesh, modified in place. |
boundary_handling | Specify to interpolate or preserve boundary edges. |
InvalidInputException | in case the input violates the precondition. |
void lscm_parameterization | ( | SurfaceMesh & | mesh | ) |
Compute parameterization based on least squares conformal mapping.
See [17] for details.
InvalidInputException | if the input precondition is violated. |
SolverException | in case of failure to solve the linear system. |
void mass_matrix | ( | const SurfaceMesh & | mesh, |
DiagonalMatrix & | M | ||
) |
void matrix_to_coordinates | ( | const DenseMatrix & | X, |
SurfaceMesh & | mesh | ||
) |
For a mesh with N vertices, set the vertex coordinates from the rows of an Nx3 matrix.
X | The input matrix. |
mesh | The mesh. |
void minimize_area | ( | SurfaceMesh & | mesh | ) |
void minimize_curvature | ( | SurfaceMesh & | mesh | ) |
SurfaceMesh plane | ( | size_t | resolution = 4 | ) |
Generate a plane mesh.
Generates a pure quad mesh in the x-y plane with origin (0,0,0) and side length 1.
resolution | Number of faces in each direction. Needs to be >= 1. Default: 4. |
SurfaceMesh quad_sphere | ( | size_t | n_subdivisions = 3 | ) |
Generate quad sphere refined by n_subdivisions
.
Uses Catmull-Clark subdivision to refine the initial hexahedron.
void quad_tri_subdivision | ( | SurfaceMesh & | mesh, |
BoundaryHandling | boundary_handling = BoundaryHandling::Interpolate |
||
) |
Perform one step of quad-tri subdivision.
Suitable for mixed quad/triangle meshes. See [22] for details.
mesh | The input mesh, modified in place. |
boundary_handling | Specify to interpolate or preserve boundary edges. |
SurfaceMesh torus | ( | size_t | radial_resolution = 20 , |
size_t | tubular_resolution = 40 , |
||
Scalar | radius = 1.0 , |
||
Scalar | thickness = 0.4 |
||
) |
Generate a torus mesh.
Generates a quad mesh of a torus with its major circle in the x-y plane.
radial_resolution | Number of subdivisions of the major circle. Needs to be >= 3. Default: 20. |
tubular_resolution | Number of subdivisions of along the tube. Needs to be >= 3. Default: 40. |
radius | Radius of the major circle. Default: 1. |
thickness | Thickness of the tube. Default: 0.4. |
void triangulate | ( | SurfaceMesh & | mesh, |
Face | f | ||
) |
Triangulate the Face f
.
Triangulate n-gons into n-2 triangles. Finds the triangulation that minimizes the sum of squared triangle areas. See [15] for details.
InvalidInputException | in case the input precondition is violated |
void uniform_laplace_matrix | ( | const SurfaceMesh & | mesh, |
SparseMatrix & | L | ||
) |
Construct the uniform Laplace matrix.
Matrix is sparse, symmetric and negative semi-definite. M(i,i) is the negative valence of vertex i. M(i,j) is +1 if vertex i and vertex j are neighbors.
mesh | The input mesh. |
L | The output matrix. |
void uniform_mass_matrix | ( | const SurfaceMesh & | mesh, |
DiagonalMatrix & | M | ||
) |
Construct the mass matrix for the uniform Laplacian.
Matrix is diagonal and positive definite. M(i,i) is the valence of vertex i.
mesh | The input mesh. |
M | The output matrix. |
void uniform_remeshing | ( | SurfaceMesh & | mesh, |
Scalar | edge_length, | ||
unsigned int | iterations = 10 , |
||
bool | use_projection = true |
||
) |
Perform uniform remeshing.
Performs incremental remeshing based on edge collapse, split, flip, and tangential relaxation. See [2] and [10] for details.
mesh | The input mesh, modified in place. |
edge_length | The target edge length. |
iterations | The number of iterations |
use_projection | Use back-projection to the input surface. |
InvalidInputException | if the input precondition is violated. |
Normal vertex_normal | ( | const SurfaceMesh & | mesh, |
Vertex | v | ||
) |
Compute the normal vector of vertex v
.
void vertex_normals | ( | SurfaceMesh & | mesh | ) |
Compute vertex normals for the whole mesh
.
Calls vertex_normal() for each vertex and adds a new vertex property of type Normal named "v:normal".
Scalar volume | ( | const SurfaceMesh & | mesh | ) |
Compute the volume of a mesh.
See [24] for details.
InvalidInputException | if the input precondition is violated. |
Scalar voronoi_area_mixed | ( | const SurfaceMesh & | mesh, |
Vertex | v | ||
) |
Compute mixed Voronoi area of a vertex.
This version is preferred for irregular triangles with obtuse angles. See [18] for details.