The pmp-library namespace. More...
Classes | |
class | AllocationException |
Exception indicating failure to allocate a new resource. More... | |
class | BoundingBox |
Simple class for representing a bounding box. More... | |
class | Edge |
this type represents an edge (internally it is basically an index) More... | |
class | EdgeProperty |
Edge property of type T. More... | |
class | Face |
this type represents a face (internally it is basically an index) More... | |
class | FaceProperty |
Face property of type T. More... | |
class | GLException |
Exception indicating an OpenGL error. More... | |
class | Halfedge |
this type represents a halfedge (internally it is basically an index) More... | |
class | HalfedgeProperty |
Halfedge property of type T. More... | |
class | Handle |
Base class for all entity handles types. More... | |
class | InvalidInputException |
Exception indicating invalid input passed to a function. More... | |
class | IOException |
Exception indicating an error occurred while performing IO. More... | |
struct | IOFlags |
Flags to control reading and writing. More... | |
class | Matrix |
Generic class for \( M \times N \) matrices. More... | |
class | MemoryUsage |
A simple class to retrieve memory usage information. More... | |
class | MeshViewer |
Simple viewer for a SurfaceMesh. More... | |
class | Renderer |
Class for rendering surface meshes using OpenGL. More... | |
class | Shader |
Class for handling shaders. More... | |
class | SolverException |
Exception indicating failure so solve an equation system. More... | |
class | StopWatch |
A simple stop watch class. More... | |
class | SurfaceMesh |
A class for representing polygon surface meshes. More... | |
class | TopologyException |
Exception indicating a topological error has occurred. More... | |
class | TrackballViewer |
A simple GLFW viewer with trackball user interface. More... | |
class | Vertex |
this type represents a vertex (internally it is basically an index) More... | |
class | VertexProperty |
Vertex property of type T. More... | |
class | Window |
A window provided by GLFW. More... | |
Typedefs | |
using | SparseMatrix = Eigen::SparseMatrix< double > |
PMP uses Eigen's double-precision sparse matrices. | |
using | DiagonalMatrix = Eigen::DiagonalMatrix< double, Eigen::Dynamic > |
PMP uses Eigen's double-precision diagonal matrices. | |
using | DenseMatrix = Eigen::MatrixXd |
PMP uses Eigen's double-precision dense matrices. | |
using | Triplet = Eigen::Triplet< double > |
PMP uses Eigen's double-precision triplets. | |
template<typename Scalar , int M> | |
using | Vector = Matrix< Scalar, M, 1 > |
template specialization for Vector as Nx1 matrix | |
template<typename Scalar > | |
using | Mat4 = Matrix< Scalar, 4, 4 > |
template specialization for 4x4 matrices | |
template<typename Scalar > | |
using | Mat3 = Matrix< Scalar, 3, 3 > |
template specialization for 3x3 matrices | |
template<typename Scalar > | |
using | Mat2 = Matrix< Scalar, 2, 2 > |
template specialization for 2x2 matrices | |
using | vec2 = Vector< float, 2 > |
template specialization for a vector of two float values | |
using | dvec2 = Vector< double, 2 > |
template specialization for a vector of two double values | |
using | bvec2 = Vector< bool, 2 > |
template specialization for a vector of two bool values | |
using | ivec2 = Vector< int, 2 > |
template specialization for a vector of two int values | |
using | uvec2 = Vector< unsigned int, 2 > |
template specialization for a vector of two unsigned int values | |
using | vec3 = Vector< float, 3 > |
template specialization for a vector of three float values | |
using | dvec3 = Vector< double, 3 > |
template specialization for a vector of three double values | |
using | bvec3 = Vector< bool, 3 > |
template specialization for a vector of three bool values | |
using | ivec3 = Vector< int, 3 > |
template specialization for a vector of three int values | |
using | uvec3 = Vector< unsigned int, 3 > |
template specialization for a vector of three unsigned int values | |
using | vec4 = Vector< float, 4 > |
template specialization for a vector of four float values | |
using | dvec4 = Vector< double, 4 > |
template specialization for a vector of four double values | |
using | bvec4 = Vector< bool, 4 > |
template specialization for a vector of four bool values | |
using | ivec4 = Vector< int, 4 > |
template specialization for a vector of four int values | |
using | uvec4 = Vector< unsigned int, 4 > |
template specialization for a vector of four unsigned int values | |
using | mat2 = Mat2< float > |
template specialization for a 2x2 matrix of float values | |
using | dmat2 = Mat2< double > |
template specialization for a 2x2 matrix of double values | |
using | mat3 = Mat3< float > |
template specialization for a 3x3 matrix of float values | |
using | dmat3 = Mat3< double > |
template specialization for a 3x3 matrix of double values | |
using | mat4 = Mat4< float > |
template specialization for a 4x4 matrix of float values | |
using | dmat4 = Mat4< double > |
template specialization for a 4x4 matrix of double values | |
using | Scalar = float |
Scalar type. | |
using | Point = Vector< Scalar, 3 > |
Point type. | |
using | Normal = Vector< Scalar, 3 > |
Normal type. | |
using | Color = Vector< Scalar, 3 > |
Color type. | |
using | TexCoord = Vector< Scalar, 2 > |
Texture coordinate type. | |
Enumerations | |
enum class | Curvature { min , max , mean , gauss , max_abs } |
Type of curvature to be computed. More... | |
Functions | |
void | curvature_to_texture_coordinates (SurfaceMesh &mesh) |
convert curvature values "v:curv" to 1D texture coordinates stored in vertex property "v:tex" | |
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 | 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 | edge_area (const SurfaceMesh &mesh, Edge e) |
Compute the area assigned to edge e . | |
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 | 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 | 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 | 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. | |
void | distance_to_texture_coordinates (SurfaceMesh &mesh) |
Use the normalized distances as texture coordinates. | |
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 | mass_matrix (const SurfaceMesh &mesh, DiagonalMatrix &M) |
Construct the (lumped) mass matrix for the cotangent Laplacian. | |
void | uniform_laplace_matrix (const SurfaceMesh &mesh, SparseMatrix &L) |
Construct the uniform Laplace matrix. | |
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. | |
Normal | face_normal (const SurfaceMesh &mesh, Face f) |
Compute the normal vector of face f . | |
Normal | vertex_normal (const SurfaceMesh &mesh, Vertex v) |
Compute the normal vector of vertex v . | |
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 | vertex_normals (SurfaceMesh &mesh) |
Compute vertex normals for the whole mesh . | |
void | face_normals (SurfaceMesh &mesh) |
Compute face normals for the whole mesh . | |
DenseMatrix | cholesky_solve (const SparseMatrix &A, const DenseMatrix &B) |
Solve the linear system A*X=B using sparse Cholesky decomposition. | |
DenseMatrix | cholesky_solve (const SparseMatrix &A, const DenseMatrix &B, const std::function< bool(unsigned int)> &is_constrained, const DenseMatrix &C) |
Solve the linear system A*X=B with given hard constraints using sparse Cholesky decomposition. | |
void | selector_matrix (const SurfaceMesh &mesh, const std::function< bool(Vertex)> &is_selected, SparseMatrix &S) |
Constructs a selector matrix for a mesh with N vertices. | |
void | matrices_to_mesh (const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, SurfaceMesh &mesh) |
Build SurfaceMesh from Eigen matrices containing vertex coordinates and triangle indices. | |
void | mesh_to_matrices (const SurfaceMesh &mesh, Eigen::MatrixXd &V, Eigen::MatrixXi &F) |
Convert SurfaceMesh to Eigen matrices of vertex coordinates and triangle indices. | |
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 | mean_edge_length (const SurfaceMesh &mesh) |
Compute mean edge length of mesh . | |
Scalar | edge_length (const SurfaceMesh &mesh, Edge e) |
Compute length of an edge e in mesh . | |
void | read (SurfaceMesh &mesh, const std::filesystem::path &file) |
Read into mesh from file . | |
void | write (const SurfaceMesh &mesh, const std::filesystem::path &file, const IOFlags &flags=IOFlags()) |
Write mesh to file controlled by flags . | |
template<typename Scalar , int M, int N> | |
std::ostream & | operator<< (std::ostream &os, const Matrix< Scalar, M, N > &m) |
output a matrix by printing its space-separated components | |
template<typename Scalar , int M, int N, int K> | |
Matrix< Scalar, M, N > | operator* (const Matrix< Scalar, M, K > &m1, const Matrix< Scalar, K, N > &m2) |
matrix-matrix multiplication | |
template<typename Scalar , int M, int N> | |
Matrix< Scalar, M, N > | cmult (const Matrix< Scalar, M, N > &m1, const Matrix< Scalar, M, N > &m2) |
component-wise multiplication | |
template<typename Scalar , int M, int N> | |
Matrix< Scalar, N, M > | transpose (const Matrix< Scalar, M, N > &m) |
transpose MxN matrix to NxM matrix | |
template<typename Scalar , int M, int N> | |
Matrix< Scalar, M, N > | operator+ (const Matrix< Scalar, M, N > &m1, const Matrix< Scalar, M, N > &m2) |
matrix addition: m1 + m2 | |
template<typename Scalar , int M, int N> | |
Matrix< Scalar, M, N > | operator- (const Matrix< Scalar, M, N > &m1, const Matrix< Scalar, M, N > &m2) |
matrix subtraction: m1 - m2 | |
template<typename Scalar , int M, int N> | |
Matrix< Scalar, M, N > | operator- (const Matrix< Scalar, M, N > &m) |
matrix negation: -m | |
template<typename Scalar , typename Scalar2 , int M, int N> | |
Matrix< Scalar, M, N > | operator* (const Scalar2 s, const Matrix< Scalar, M, N > &m) |
scalar multiplication of matrix: s*m | |
template<typename Scalar , typename Scalar2 , int M, int N> | |
Matrix< Scalar, M, N > | operator* (const Matrix< Scalar, M, N > &m, const Scalar2 s) |
scalar multiplication of matrix: m*s | |
template<typename Scalar , typename Scalar2 , int M, int N> | |
Matrix< Scalar, M, N > | operator/ (const Matrix< Scalar, M, N > &m, const Scalar2 s) |
divide matrix by scalar: m/s | |
template<typename Scalar , int M, int N> | |
Scalar | norm (const Matrix< Scalar, M, N > &m) |
compute the Frobenius norm of a matrix (or Euclidean norm of a vector) | |
template<typename Scalar , int M, int N> | |
Scalar | sqrnorm (const Matrix< Scalar, M, N > &m) |
compute the squared Frobenius norm of a matrix (or squared Euclidean norm of a vector) | |
template<typename Scalar , int M, int N> | |
Matrix< Scalar, M, N > | normalize (const Matrix< Scalar, M, N > &m) |
return a normalized copy of a matrix or a vector | |
template<typename Scalar , int M, int N> | |
Matrix< Scalar, M, N > | min (const Matrix< Scalar, M, N > &m1, const Matrix< Scalar, M, N > &m2) |
return component-wise minimum | |
template<typename Scalar , int M, int N> | |
Matrix< Scalar, M, N > | max (const Matrix< Scalar, M, N > &m1, const Matrix< Scalar, M, N > &m2) |
return component-wise maximum | |
template<typename Scalar > | |
Mat4< Scalar > | viewport_matrix (Scalar l, Scalar b, Scalar w, Scalar h) |
OpenGL viewport matrix with parameters left, bottom, width, height. | |
template<typename Scalar > | |
Mat4< Scalar > | inverse_viewport_matrix (Scalar l, Scalar b, Scalar w, Scalar h) |
inverse of OpenGL viewport matrix with parameters left, bottom, width, height | |
template<typename Scalar > | |
Mat4< Scalar > | frustum_matrix (Scalar l, Scalar r, Scalar b, Scalar t, Scalar n, Scalar f) |
OpenGL frustum matrix with parameters left, right, bottom, top, near, far. | |
template<typename Scalar > | |
Mat4< Scalar > | inverse_frustum_matrix (Scalar l, Scalar r, Scalar b, Scalar t, Scalar n, Scalar f) |
inverse of OpenGL frustum matrix with parameters left, right, bottom, top, near, far | |
template<typename Scalar > | |
Mat4< Scalar > | perspective_matrix (Scalar fovy, Scalar aspect, Scalar zNear, Scalar zFar) |
OpenGL perspective matrix with parameters field of view in y-direction, aspect ratio, and distance of near and far planes. | |
template<typename Scalar > | |
Mat4< Scalar > | inverse_perspective_matrix (Scalar fovy, Scalar aspect, Scalar zNear, Scalar zFar) |
inverse of perspective matrix | |
template<typename Scalar > | |
Mat4< Scalar > | ortho_matrix (Scalar left, Scalar right, Scalar bottom, Scalar top, Scalar zNear=-1, Scalar zFar=1) |
OpenGL orthogonal projection matrix with parameters left, right, bottom, top, near, far. | |
template<typename Scalar > | |
Mat4< Scalar > | look_at_matrix (const Vector< Scalar, 3 > &eye, const Vector< Scalar, 3 > ¢er, const Vector< Scalar, 3 > &up) |
OpenGL look-at camera matrix with parameters eye position, scene center, up-direction. | |
template<typename Scalar > | |
Mat4< Scalar > | translation_matrix (const Vector< Scalar, 3 > &t) |
OpenGL matrix for translation by vector t. | |
template<typename Scalar > | |
Mat4< Scalar > | scaling_matrix (const Scalar s) |
OpenGL matrix for scaling x/y/z by s. | |
template<typename Scalar > | |
Mat4< Scalar > | scaling_matrix (const Vector< Scalar, 3 > &s) |
OpenGL matrix for scaling x/y/z by the components of s. | |
template<typename Scalar > | |
Mat4< Scalar > | rotation_matrix_x (Scalar angle) |
OpenGL matrix for rotation around x-axis by given angle (in degrees) | |
template<typename Scalar > | |
Mat4< Scalar > | rotation_matrix_y (Scalar angle) |
OpenGL matrix for rotation around y-axis by given angle (in degrees) | |
template<typename Scalar > | |
Mat4< Scalar > | rotation_matrix_z (Scalar angle) |
OpenGL matrix for rotation around z-axis by given angle (in degrees) | |
template<typename Scalar > | |
Mat4< Scalar > | rotation_matrix (const Vector< Scalar, 3 > &axis, Scalar angle) |
OpenGL matrix for rotation around given axis by given angle (in degrees) | |
template<typename Scalar > | |
Mat4< Scalar > | rotation_matrix (const Vector< Scalar, 4 > &quat) |
OpenGL matrix for rotation specified by unit quaternion. | |
template<typename Scalar > | |
Mat3< Scalar > | linear_part (const Mat4< Scalar > &m) |
return upper 3x3 matrix from given 4x4 matrix, corresponding to the linear part of an affine transformation | |
template<typename Scalar > | |
Vector< Scalar, 3 > | projective_transform (const Mat4< Scalar > &m, const Vector< Scalar, 3 > &v) |
projective transformation of 3D vector v by a 4x4 matrix m: add 1 as 4th component of v, multiply m*v, divide by 4th component | |
template<typename Scalar > | |
Vector< Scalar, 3 > | affine_transform (const Mat4< Scalar > &m, const Vector< Scalar, 3 > &v) |
affine transformation of 3D vector v by a 4x4 matrix m: add 1 as 4th component of v, multiply m*v, do NOT divide by 4th component | |
template<typename Scalar > | |
Vector< Scalar, 3 > | linear_transform (const Mat4< Scalar > &m, const Vector< Scalar, 3 > &v) |
linear transformation of 3D vector v by a 4x4 matrix m: transform vector by upper-left 3x3 submatrix of m | |
template<typename Scalar > | |
Mat4< Scalar > | inverse (const Mat4< Scalar > &m) |
return the inverse of a 4x4 matrix | |
template<typename Scalar > | |
Scalar | determinant (const Mat3< Scalar > &m) |
return determinant of 3x3 matrix | |
template<typename Scalar > | |
Mat3< Scalar > | inverse (const Mat3< Scalar > &m) |
return the inverse of a 3x3 matrix | |
template<typename Scalar > | |
bool | symmetric_eigendecomposition (const Mat3< Scalar > &m, Scalar &eval1, Scalar &eval2, Scalar &eval3, Vector< Scalar, 3 > &evec1, Vector< Scalar, 3 > &evec2, Vector< Scalar, 3 > &evec3) |
compute eigenvector/eigenvalue decomposition of a 3x3 matrix | |
template<typename Scalar , int N> | |
std::istream & | operator>> (std::istream &is, Vector< Scalar, N > &vec) |
read the space-separated components of a vector from a stream | |
template<typename Scalar , int N> | |
std::ostream & | operator<< (std::ostream &os, const Vector< Scalar, N > &vec) |
output a vector by printing its space-separated components | |
template<typename Scalar , int N> | |
Scalar | dot (const Vector< Scalar, N > &v0, const Vector< Scalar, N > &v1) |
compute the dot product of two vectors | |
template<typename Scalar , int N> | |
Scalar | distance (const Vector< Scalar, N > &v0, const Vector< Scalar, N > &v1) |
compute the Euclidean distance between two points | |
template<typename Scalar > | |
Vector< Scalar, 2 > | perp (const Vector< Scalar, 2 > &v) |
compute perpendicular vector (rotate vector counter-clockwise by 90 degrees) | |
template<typename Scalar > | |
Vector< Scalar, 3 > | cross (const Vector< Scalar, 3 > &v0, const Vector< Scalar, 3 > &v1) |
compute the cross product of two vectors (only valid for 3D vectors) | |
std::ostream & | operator<< (std::ostream &os, const StopWatch &watch) |
output a elapsed time to a stream | |
std::ostream & | operator<< (std::ostream &os, Vertex v) |
output a Vertex to a stream | |
std::ostream & | operator<< (std::ostream &os, Halfedge h) |
output a Halfedge to a stream | |
std::ostream & | operator<< (std::ostream &os, Edge e) |
output an Edge to a stream | |
std::ostream & | operator<< (std::ostream &os, Face f) |
output a Face to a stream | |
The pmp-library namespace.
DenseMatrix cholesky_solve | ( | const SparseMatrix & | A, |
const DenseMatrix & | B | ||
) |
Solve the linear system A*X=B using sparse Cholesky decomposition.
Returns the solution vector/matrix X.
A | The system matrix. |
B | The right hand side. |
DenseMatrix cholesky_solve | ( | const SparseMatrix & | A, |
const DenseMatrix & | B, | ||
const std::function< bool(unsigned int)> & | is_constrained, | ||
const DenseMatrix & | C | ||
) |
Solve the linear system A*X=B with given hard constraints using sparse Cholesky decomposition.
Returns the solution vector or matrix X.
A | The system matrix. |
B | The right hand side. |
is_constrained | A function returning whether or not X(i) is constrained or not. |
C | A matrix storing the Dirichlet constraints: X(i) should be C(i) is entry i is constrained. |
void distance_to_texture_coordinates | ( | SurfaceMesh & | mesh | ) |
Use the normalized distances as texture coordinates.
Stores the normalized distances in a vertex property of type TexCoord named "v:tex". Re-uses any existing vertex property of the same type and name.
void matrices_to_mesh | ( | const Eigen::MatrixXd & | V, |
const Eigen::MatrixXi & | F, | ||
SurfaceMesh & | mesh | ||
) |
Build SurfaceMesh from Eigen matrices containing vertex coordinates and triangle indices.
V | \(n\times 3\) matrix of double precision vertex coordinates. |
F | \(m\times 3\) matrix of integer triangle indices. |
mesh | The mesh to be build from V and F . The mesh will be cleared. |
void mesh_to_matrices | ( | const SurfaceMesh & | mesh, |
Eigen::MatrixXd & | V, | ||
Eigen::MatrixXi & | F | ||
) |
Convert SurfaceMesh to Eigen matrices of vertex coordinates and triangle indices.
mesh | The mesh used to fill V and F . |
V | The resulting \(n\times 3\) matrix of double precision vertex coordinates. |
F | The resulting \(m\times 3\) matrix of integer triangle indices. |
void selector_matrix | ( | const SurfaceMesh & | mesh, |
const std::function< bool(Vertex)> & | is_selected, | ||
SparseMatrix & | S | ||
) |
Constructs a selector matrix for a mesh with N vertices.
Returns a matrix built from the rows of the NxN identity matrix that belong to selected vertices.
mesh | The input mesh. |
is_selected | A function returning whether or not vertex i is selected or not. |
S | The output matrix. |