CortidQCT  1.2.2.52
MeshHelpers.h
Go to the documentation of this file.
1 
12 #pragma once
13 
14 #include "Mesh.h"
15 
16 #include <Eigen/Core>
17 #include <gsl/gsl>
18 #include <igl/cotmatrix.h>
19 #include <igl/per_vertex_normals.h>
20 
21 namespace CortidQCT {
22 
23 namespace Internal {
24 
25 template <class T = float>
26 using VertexMatrix = Eigen::Matrix<T, Eigen::Dynamic, 3>;
27 template <class T = float> using NormalMatrix = VertexMatrix<T>;
28 using FacetMatrix =
29  Eigen::Matrix<typename Mesh<float>::Index, Eigen::Dynamic, 3>;
30 using LabelVector =
31  Eigen::Matrix<typename Mesh<float>::Label, Eigen::Dynamic, 1>;
32 template <class T = float> using LaplacianMatrix = Eigen::SparseMatrix<T>;
33 
35 template <class T> inline VertexMatrix<T> vertexMatrix(Mesh<T> const &mesh) {
36  return mesh
37  .withUnsafeVertexPointer([&mesh](auto const *ptr) {
38  return Eigen::Map<Eigen::Matrix<T, 3, Eigen::Dynamic> const>{
39  ptr, 3, gsl::narrow<Eigen::Index>(mesh.vertexCount())};
40  })
41  .transpose();
42 }
43 
45 template <class T> inline FacetMatrix facetMatrix(Mesh<T> const &mesh) {
46  using Index = typename Mesh<T>::Index;
47  return mesh
48  .withUnsafeIndexPointer([&mesh](auto const *ptr) {
49  return Eigen::Map<Eigen::Matrix<Index, 3, Eigen::Dynamic> const>{
50  ptr, 3, gsl::narrow<Eigen::Index>(mesh.triangleCount())};
51  })
52  .transpose();
53 }
54 
56 template <class T> inline LabelVector labelVector(Mesh<T> const &mesh) {
57  using Label = typename Mesh<T>::Label;
58  return mesh.withUnsafeLabelPointer([&mesh](auto const *ptr) {
59  return Eigen::Map<Eigen::Matrix<Label, Eigen::Dynamic, 1> const>{
60  ptr, gsl::narrow<Eigen::Index>(mesh.vertexCount())};
61  });
62 }
63 
65 template <class DerivedV, class DerivedF>
66 inline NormalMatrix<typename DerivedV::Scalar>
67 perVertexNormalMatrix(Eigen::MatrixBase<DerivedV> const &V,
68  Eigen::MatrixBase<DerivedF> const &F) {
69  // output matrix
70  NormalMatrix<typename DerivedV::Scalar> normals;
71 
72  igl::per_vertex_normals(V, F, igl::PER_VERTEX_NORMALS_WEIGHTING_TYPE_ANGLE,
73  normals);
74 
75  return normals;
76 }
77 
79 template <class T>
80 inline NormalMatrix<T> perVertexNormalMatrix(Mesh<T> const &mesh) {
81  return perVertexNormalMatrix(vertexMatrix(mesh), facetMatrix(mesh));
82 }
83 
85 template <class DerivedV, class DerivedF>
86 inline LaplacianMatrix<typename DerivedV::Scalar>
87 laplacianMatrix(Eigen::MatrixBase<DerivedV> const &V,
88  Eigen::MatrixBase<DerivedF> const &F) {
89  // output matrix
90  Eigen::SparseMatrix<typename DerivedV::Scalar> laplacian;
91  Expects(V.rows() > 0 && V.cols() == 3);
92  Expects(F.rows() > 0 && F.cols() == 3);
93  igl::cotmatrix(V, F, laplacian);
94 
95  return laplacian;
96 }
97 
99 template <class T>
100 inline LaplacianMatrix<T> laplacianMatrix(Mesh<T> const &mesh) {
101  return laplacianMatrix(vertexMatrix(mesh), facetMatrix(mesh));
102 }
103 
104 } // namespace Internal
105 } // namespace CortidQCT
Name namespace for CortidQCT library.
Definition: CortidQCT.h:23
auto withUnsafeIndexPointer(F &&f) const noexcept(noexcept(f(std::declval< const IndexData >().data())))
Calls the given functional with an unsafe pointer to the raw index storage.
Definition: Mesh.h:372
VertexMatrix< T > vertexMatrix(Mesh< T > const &mesh)
Returns the Nx3 vertex matrix of the mesh.
Definition: MeshHelpers.h:35
std::ptrdiff_t Index
Index type.
Definition: Mesh.h:45
auto withUnsafeVertexPointer(F &&f) const noexcept(noexcept(f(std::declval< const VertexData >().data())))
Calls the given functional with an unsafe pointer to the raw vertex storage.
Definition: Mesh.h:349
NormalMatrix< typename DerivedV::Scalar > perVertexNormalMatrix(Eigen::MatrixBase< DerivedV > const &V, Eigen::MatrixBase< DerivedF > const &F)
Returns a Nx3 matrix with per-vertex normals.
Definition: MeshHelpers.h:67
A triangle mesh class.
Definition: BarycentricPoint.h:20
Include file defining the Mesh data type.
LabelVector labelVector(Mesh< T > const &mesh)
Returns a N-vector containing the per vertex labels.
Definition: MeshHelpers.h:56
unsigned int Label
Label type.
Definition: Mesh.h:51
Size triangleCount() const noexcept
Number of triangles.
Definition: Mesh.h:89
FacetMatrix facetMatrix(Mesh< T > const &mesh)
Returns the Mx3 index matrix of the mesh.
Definition: MeshHelpers.h:45
Size vertexCount() const noexcept
Number of vertices.
Definition: Mesh.h:86
auto withUnsafeLabelPointer(F &&f) const noexcept(noexcept(f(LabelData().data())))
Calls the given functional with an unsafe pointer to the raw label storage.
Definition: Mesh.h:395
LaplacianMatrix< typename DerivedV::Scalar > laplacianMatrix(Eigen::MatrixBase< DerivedV > const &V, Eigen::MatrixBase< DerivedF > const &F)
Returns the NxN sparse laplacian matrix (using cotangent weights)
Definition: MeshHelpers.h:87