Polygon Mesh Processing Library
Loading...
Searching...
No Matches
algorithms

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 angle between two (un-normalized) vectors
 
Scalar sin (const Point &v0, const Point &v1)
 compute sine of angle between two (un-normalized) vectors
 
Scalar cos (const Point &v0, const Point &v1)
 compute cosine of angle between two (un-normalized) vectors
 
Scalar cotan (const Point &v0, const Point &v1)
 compute cotangent of angle between two (un-normalized) vectors
 
Scalar triangle_area (const Point &p0, const Point &p1, const Point &p2)
 compute 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 sum of face areas).
 
Scalar voronoi_area (const SurfaceMesh &mesh, Vertex v)
 compute (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 area assigned to edge e (a face with n edges assigns 1/n of its area to each edge)
 
Scalar volume (const SurfaceMesh &mesh)
 Compute the volume of a mesh.
 
Point centroid (const SurfaceMesh &mesh, Face f)
 compute barycenter/centroid of a face
 
Point centroid (const SurfaceMesh &mesh)
 barycenter/centroid of mesh, computed as area-weighted mean of vertices.
 
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 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,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 &center=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 .
 

Detailed Description

Mesh processing algorithms.

Enumeration Type Documentation

◆ Curvature

enum class Curvature
strong

Type of curvature to be computed.

Enumerator
min 

minimum curvature

max 

maximum curvature

mean 

mean curvature

gauss 

Gauss curvature.

max_abs 

maximum absolute curvature

Function Documentation

◆ adaptive_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.

Performs incremental remeshing based on edge collapse, split, flip, and tangential relaxation. See [2] and [10] for details.

Parameters
meshThe input mesh, modified in place.
min_edge_lengthThe minimum edge length.
max_edge_lengthThe maximum edge length.
approx_errorThe maximum approximation error.
iterationsThe number of iterations.
use_projectionUse back-projection to the input surface.
Precondition
Input mesh needs to be a triangle mesh.
Exceptions
InvalidInputExceptionif the input precondition is violated.

◆ catmull_clark_subdivision()

void catmull_clark_subdivision ( SurfaceMesh mesh,
BoundaryHandling  boundary_handling = BoundaryHandling::Interpolate 
)

Perform one step of Catmull-Clark subdivision.

See [5] for details.

Parameters
meshThe input mesh, modified in place.
boundary_handlingSpecify to interpolate or preserve boundary edges.

◆ clear_features()

void clear_features ( SurfaceMesh mesh)

Clear feature and boundary edges.

Sets all "e:feature" and "v:feature" properties to false.

Note
This does not remove the corresponding property arrays.

◆ cone()

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.

Parameters
n_subdivisionsNumber of subdivisions of the base circle. Needs to be >= 3. Default: 30.
radiusRadius of the base circle. Default: 1.
heightHeight of the the cone. Default: 2.5.

◆ coordinates_to_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.

Parameters
meshThe input mesh.
XThe output matrix.

◆ corner_normal()

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.

Note
This algorithm works on general polygon meshes.

◆ cotan_weight()

double cotan_weight ( const SurfaceMesh mesh,
Edge  e 
)

compute the cotangent weight for edge e

Precondition
Input mesh needs to be a triangle mesh.

◆ curvature()

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).

Curvature values for boundary vertices are interpolated from their interior neighbors. Curvature values can be smoothed. See [18] and [6] for details.

Note
This algorithm works on general polygon meshes.

◆ cylinder()

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.

Parameters
n_subdivisionsNumber of subdivisions of the cylinder. Needs to be >= 3. Default: 30.
radiusRadius of the cylinder. Default: 1.
heightHeight of the cylinder. Default: 2.5.

◆ decimate()

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.

Parameters
meshTarget mesh. Modified in place.
n_verticesTarget number of vertices.
aspect_ratioMinimum aspect ratio of the triangles.
edge_lengthMinimum target edge length.
max_valenceMaximum number of incident edges per vertex.
normal_deviationMaximum deviation of face normals.
hausdorff_errorMaximum deviation from the original surface.
seam_thresholdThreshold for texture seams.
seam_angle_deviationMaximum texture seam deviation.
Precondition
Input mesh needs to be a triangle mesh.
Exceptions
InvalidInputExceptionif the input precondition is violated.

◆ detect_boundary()

size_t detect_boundary ( SurfaceMesh mesh)

Mark all boundary edges as features.

Returns
The number of boundary edges detected.

◆ detect_features()

size_t detect_features ( SurfaceMesh mesh,
Scalar  angle 
)

Mark edges with dihedral angle larger than angle as feature.

Returns
The number of feature edges detected.

◆ divergence_matrix()

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.

Parameters
meshThe input mesh.
DThe output matrix.
See also
laplace_matrix
gradient_matrix

◆ dual()

void dual ( SurfaceMesh mesh)

Compute dual of a mesh.

Warning
Changes the mesh in place. All properties are cleared.

◆ explicit_smoothing()

void explicit_smoothing ( SurfaceMesh mesh,
unsigned int  iterations = 10,
bool  use_uniform_laplace = false 
)

Perform explicit Laplacian smoothing.

See [8] for details

Note
This algorithm works on general polygon meshes.
Parameters
meshThe input mesh, modified in place.
iterationsThe number of iterations performed.
use_uniform_laplaceUse uniform or cotan Laplacian. Default: cotan.

◆ face_area()

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.

◆ face_normal()

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]

Note
This algorithm works on general polygon meshes.

◆ face_normals()

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".

Note
This algorithm works on general polygon meshes.

◆ fair()

void fair ( SurfaceMesh mesh,
unsigned int  k = 2 
)

Implicit surface fairing.

Computes a surface by solving k-harmonic equation. See also [8] .

Note
This algorithm works on general polygon meshes.
Exceptions
SolverExceptionin case of failure to solve the linear system
InvalidInputExceptionin case of missing boundary constraints

◆ fill_hole()

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.

Precondition
The specified halfedge is valid.
The specified halfedge is a boundary halfedge.
The specified halfedge is not adjacent to a non-manifold hole.
Exceptions
InvalidInputExceptionin case on of the input preconditions is violated
Note
This algorithm works on general polygon meshes.

◆ geodesics()

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.

Parameters
meshThe input mesh, modified in place.
[in]seedsThe vector of seed vertices.
[in]maxdistThe maximum distance up to which to compute the geodesic distances.
[in]maxnumThe maximum number of neighbors up to which to compute the geodesic distances.
[out]neighborsThe vector of neighbor vertices.
Returns
The number of neighbors that have been found.
Precondition
Input mesh needs to be a triangle mesh.

◆ geodesics_heat()

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.

Parameters
meshThe input mesh, modified in place.
seedsThe vector of seed vertices.
Note
This algorithm works on general polygon meshes.

◆ gradient_matrix()

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.

Parameters
meshThe input mesh.
GThe output matrix.
See also
laplace_matrix
divergence_matrix

◆ harmonic_parameterization()

void harmonic_parameterization ( SurfaceMesh mesh,
bool  use_uniform_weights = false 
)

Compute discrete harmonic parameterization.

See [9] for details.

Note
This algorithm works on general polygon meshes.
Precondition
The mesh has a boundary.
Exceptions
InvalidInputExceptionif the input precondition is violated.
SolverExceptionin case of failure to solve the linear system.
Note
This algorithm works on general polygon meshes.

◆ icosphere()

SurfaceMesh icosphere ( size_t  n_subdivisions = 3)

Generate icosphere refined by n_subdivisions .

Uses Loop subdivision to refine the initial icosahedron.

◆ implicit_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.

See [8] and [12] .

Note
This algorithm works on general polygon meshes.
Parameters
meshThe input mesh, modified in place.
timestepThe time step taken.
iterationsThe number of iterations performed.
use_uniform_laplaceUse uniform or cotan Laplacian. Default: cotan.
rescaleRe-center and re-scale model after smoothing. Default: true.
Exceptions
SolverExceptionin case of a failure to solve the linear system.

◆ laplace()

Point laplace ( const SurfaceMesh mesh,
Vertex  v 
)

compute Laplace vector for vertex v (normalized by Voronoi area)

Precondition
Input mesh needs to be a triangle mesh.

◆ laplace_matrix()

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.

Parameters
meshThe input mesh.
LThe output matrix.
clampWhether or not negative off-diagonal entries should be clamped to zero.
See also
gradient_matrix
divergence_matrix

◆ linear_subdivision()

void linear_subdivision ( SurfaceMesh mesh)

Perform one step of linear quad-tri subdivision.

Suitable for mixed quad/triangle meshes.

Parameters
meshThe input mesh, modified in place.

◆ loop_subdivision()

void loop_subdivision ( SurfaceMesh mesh,
BoundaryHandling  boundary_handling = BoundaryHandling::Interpolate 
)

Perform one step of Loop subdivision.

See [16] for details.

Parameters
meshThe input mesh, modified in place.
boundary_handlingSpecify to interpolate or preserve boundary edges.
Precondition
Requires a triangle mesh as input.
Exceptions
InvalidInputExceptionin case the input violates the precondition.

◆ lscm_parameterization()

void lscm_parameterization ( SurfaceMesh mesh)

Compute parameterization based on least squares conformal mapping.

See [17] for details.

Precondition
The mesh has a boundary.
Input mesh needs to be a triangle mesh.
Exceptions
InvalidInputExceptionif the input precondition is violated.
SolverExceptionin case of failure to solve the linear system.

◆ mass_matrix()

void mass_matrix ( const SurfaceMesh mesh,
DiagonalMatrix M 
)

Construct the (lumped) mass matrix for the cotangent Laplacian.

Matrix is diagonal and positive definite. M(i,i) is the (mixed) Voronoi area of vertex i. See [18] for details on triangle meshes and [4] for details on polygon meshes.

Parameters
meshThe input mesh.
MThe output matrix.

◆ matrix_to_coordinates()

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.

Parameters
XThe input matrix.
meshThe mesh.

◆ minimize_area()

void minimize_area ( SurfaceMesh mesh)

Minimize surface area.

Note
This algorithm works on general polygon meshes.
See also
fair()

◆ minimize_curvature()

void minimize_curvature ( SurfaceMesh mesh)

Minimize surface curvature.

Note
This algorithm works on general polygon meshes.
See also
fair()

◆ plane()

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.

Parameters
resolutionNumber of faces in each direction. Needs to be >= 1. Default: 4.

◆ quad_sphere()

SurfaceMesh quad_sphere ( size_t  n_subdivisions = 3)

Generate quad sphere refined by n_subdivisions .

Uses Catmull-Clark subdivision to refine the initial hexahedron.

◆ quad_tri_subdivision()

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.

Parameters
meshThe input mesh, modified in place.
boundary_handlingSpecify to interpolate or preserve boundary edges.

◆ torus()

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.

Parameters
radial_resolutionNumber of subdivisions of the major circle. Needs to be >= 3. Default: 20.
tubular_resolutionNumber of subdivisions of along the tube. Needs to be >= 3. Default: 40.
radiusRadius of the major circle. Default: 1.
thicknessThickness of the tube. Default: 0.4.

◆ triangulate()

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.

Precondition
The input face is manifold
Note
This algorithm works on general polygon meshes.
Exceptions
InvalidInputExceptionin case the input precondition is violated

◆ uniform_laplace_matrix()

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.

Parameters
meshThe input mesh.
LThe output matrix.

◆ uniform_mass_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.

Parameters
meshThe input mesh.
MThe output matrix.

◆ uniform_remeshing()

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.

Parameters
meshThe input mesh, modified in place.
edge_lengthThe target edge length.
iterationsThe number of iterations
use_projectionUse back-projection to the input surface.
Precondition
Input mesh needs to be a triangle mesh.
Exceptions
InvalidInputExceptionif the input precondition is violated.

◆ vertex_normal()

Normal vertex_normal ( const SurfaceMesh mesh,
Vertex  v 
)

Compute the normal vector of vertex v.

Note
This algorithm works on general polygon meshes.

◆ vertex_normals()

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".

Note
This algorithm works on general polygon meshes.

◆ volume()

Scalar volume ( const SurfaceMesh mesh)

Compute the volume of a mesh.

See [24] for details.

Precondition
Input mesh needs to be a triangle mesh.
Exceptions
InvalidInputExceptionif the input precondition is violated.

◆ voronoi_area_mixed()

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.

Precondition
Input mesh needs to be a triangle mesh.