37#include <unordered_map>
38#include "DGtal/base/ConstAlias.h"
39#include "DGtal/base/Common.h"
40#include "DGtal/shapes/SurfaceMesh.h"
41#include "DGtal/math/linalg/EigenSupport.h"
60 template <
typename TRealPo
int,
typename TRealVector>
91 const auto nn = (*myNormals)[f];
93 return p - nn.dot(p)*nn;
127template <
typename TRealPo
int,
typename TRealVector>
175 bool globalInternalCacheEnabled =
false):
191 bool globalInternalCacheEnabled =
false):
206 bool globalInternalCacheEnabled =
false):
227 bool globalInternalCacheEnabled =
false) :
303 for(
auto v: vertices)
326 for(
auto i=0u; i < nf; ++i)
360 for(
auto i=0u; i < nf; ++i)
362 a(i, (i+1)%nf) = 0.5;
378 auto it = vertices.cbegin();
379 auto itnext = vertices.cbegin();
381 while (it != vertices.cend())
385 af += xi.crossProduct(xip);
388 if (itnext == vertices.cend())
389 itnext = vertices.cbegin();
391 Eigen::Vector3d output = {af[0],af[1],af[2]};
419 return {v(0),v(1),v(2)};
439 brack << 0.0 , -n(2), n(1),
467 DenseMatrix op =
E(f)*( DenseMatrix::Identity(3,3) - n*n.transpose());
489 return 1.0/(double)nf *
X(f).transpose() * Vector::Ones(nf);
497 return {c(0),c(1),c(2)};
510 (
B(f).transpose() -
centroid(f)* Vector::Ones(nf).transpose() );
609 DenseMatrix op = -1.0 * Df.transpose() *
M(f,lambda) * Df;
626 Eigen::Vector3d nv =
n_v(v);
627 ASSERT(std::abs(nv.norm() - 1.0) < 0.001);
629 auto neighbor = *N.begin();
632 Eigen::Vector3d w =
toVec3(tangentVector);
633 Eigen::Vector3d uu =
project(w,nv).normalized();
634 Eigen::Vector3d vv = nv.cross(uu);
647 ASSERT(std::abs(nf.norm() - 1.0) < 0.001);
649 auto v1 = *(N.begin());
650 auto v2 = *(N.begin() + 1);
653 Eigen::Vector3d w =
toVec3(tangentVector);
654 Eigen::Vector3d uu =
project(w,nf).normalized();
655 Eigen::Vector3d vv = nf.cross(uu);
671 return T.col(0) * I(0) + T.col(1) * I(1);
691 Eigen::Vector3d nv =
n_v(v);
692 double c = nv.dot(nf);
693 ASSERT(std::abs( c + 1.0) > 0.0001);
695 if (std::abs( c + 1.0) < 0.00001)
696 return -Eigen::Matrix3d::Identity();
698 auto vv = nv.cross(nf);
700 return Eigen::Matrix3d::Identity() + skew +
701 1.0 / (1.0 + c) * skew * skew;
708 return Tf(f).transpose() *
Qvf(v,f) *
Tv(v);
719 N.block(cpt,0,3,1) =
n_v(v).transpose();
724 return 0.5 *
Tf.transpose() * (GN + GN.transpose()) *
Tf;
740 uf_nabla.block(cpt,0,1,2) =
741 (
Rvf(v,f) * uf.block(2 * cpt,0,2,1)).transpose();
858 std::vector<Triplet> triplets;
864 for(
auto i=0u; i < nf; ++i)
865 for(
auto j=0u; j < nf; ++j)
869 triplets.emplace_back(
Triplet( (SparseMatrix::StorageIndex)vertices[ i ], (SparseMatrix::StorageIndex)vertices[ j ],
873 lapGlobal.setFromTriplets(triplets.begin(), triplets.end());
885 std::vector<Triplet> triplets;
892 triplets.emplace_back(
Triplet(v,v,varea));
894 M.setFromTriplets(triplets.begin(),triplets.end());
905 for (
int k = 0; k < iM0.outerSize(); ++k )
906 for (
typename SparseMatrix::InnerIterator it( iM0, k ); it; ++it )
907 it.valueRef() = 1.0 / it.value();
933 std::vector<Triplet> triplets;
939 for (
auto i = 0u; i < nf; ++i)
940 for (
auto j = 0u; j < nf; ++j)
941 for (
short k1 = 0; k1 < 2; k1++)
942 for (
short k2 = 0; k2 < 2; k2++)
944 auto v = Lap(2 * i + k1, 2 * j + k2);
946 triplets.emplace_back(
Triplet(2 * vertices[i] + k1,
947 2 * vertices[j] + k2, v));
950 lapGlobal.setFromTriplets(triplets.begin(), triplets.end());
964 std::vector<Triplet> triplets;
971 triplets.emplace_back(
Triplet(2 * v, 2 * v, varea));
972 triplets.emplace_back(
Triplet(2 * v + 1, 2 * v + 1, varea));
974 M.setFromTriplets(triplets.begin(), triplets.end());
1001 std::vector<DenseMatrix> cache;
1003 cache.push_back(perFaceOperator(f));
1025 std::vector<Vector> cache;
1027 cache.push_back(perFaceVectorOperator(f));
1099 out <<
"[PolygonalCalculus]: ";
1101 out<<
"internal cache enabled, ";
1103 out<<
"internal cache disabled, ";
1126 enum OPERATOR {
X_,
D_,
E_,
A_,
COGRAD_,
GRAD_,
FLAT_,
B_,
SHARP_,
P_,
M_,
DIVERGENCE_,
CURL_,
L_,
CON_L_};
1135 auto nf = vertices.size();
1169 return u - (u.dot(n) / n.squaredNorm()) * n;
1179 for (
int i = 0; i < 3; i++)
1189 return Eigen::Vector3d(x(0),x(1),x(2));
1198 return { x(0), x(1), x(2)};
1218 for (
auto f : faces)
1221 if (fabs(n.norm() - 0.0) < 0.00001)
1224 trace.
warning()<<
"[PolygonalCalculus] Trying to compute the normal vector at a boundary vertex incident to pnon-manifold edge, we return a random vector."<<std::endl;
1225 n << Vector::Random(3);
1243 DenseMatrix RU_fO = DenseMatrix::Zero(nf * 2,nf * 2);
1248 RU_fO.block<2,2>(2 * cpt,2 * cpt) = Rv;
1257 size_t h =
M.rows();
1258 size_t w =
M.cols();
1260 for (
size_t j = 0; j < h; j++)
1261 for (
size_t i = 0; i < w; i++)
1263 MK(2 * j, 2 * i) =
M(j, i);
1264 MK(2 * j + 1, 2 * i + 1) =
M(j, i);
1301template <
typename TP,
typename TV>
1305 object.selfDisplay( out );
Aim: This class encapsulates its parameter class so that to indicate to the user that the object/poin...
Implements differential operators on polygonal surfaces from degoes2020discrete.
Eigen::Vector3d n_v(const Vertex &v) const
DenseMatrix Rvf(const Vertex &v, const Face &f) const
PolygonalCalculus(const ConstAlias< MySurfaceMesh > surf, const std::function< Real3dPoint(Face, Vertex)> &pos_embedder, const std::function< Vector(Vertex)> &normal_embedder, bool globalInternalCacheEnabled=false)
EigenLinearAlgebraBackend LinAlg
Linear Algebra Backend from Eigen.
DenseMatrix transportAndFormatVectorField(const Face f, const Vector &uf)
SurfaceMesh< TRealPoint, TRealVector > MySurfaceMesh
Type of SurfaceMesh.
std::function< Real3dPoint(Face, Vertex)> myEmbedder
Embedding function (face,vertex)->R^3 for the vertex position wrt. the face.
DenseMatrix D(const Face f) const
DenseMatrix laplaceBeltrami(const Face f, const double lambda=1.0) const
static Eigen::Vector3d toVec3(const Real3dPoint &x)
toVec3 convert Real3dPoint to Eigen::Vector3d
DenseMatrix Tv(const Vertex &v) const
LinAlg::SparseMatrix SparseMatrix
Type of sparse matrix.
OPERATOR
Enum for operators in the internal cache strategy.
std::vector< size_t > myFaceDegree
Cache containing the face degree.
DenseMatrix X(const Face f) const
void selfDisplay(std::ostream &out) const
DenseMatrix Tf(const Face &f) const
PolygonalCalculus(PolygonalCalculus &&other)=delete
size_t nbVertices() const
std::vector< DenseMatrix > getOperatorCacheMatrix(const std::function< DenseMatrix(Face)> &perFaceOperator) const
Vector faceNormal(const Face f) const
Vector centroid(const Face f) const
const MySurfaceMesh * mySurfaceMesh
Underlying SurfaceMesh.
MySurfaceMesh::RealVector Real3dVector
Real vector type.
DenseMatrix curl(const Face f) const
static const Dimension dimension
Concept checking.
Real3dVector faceNormalAsDGtalVector(const Face f) const
PolygonalCalculus(const ConstAlias< MySurfaceMesh > surf, const std::function< Real3dPoint(Face, Vertex)> &embedder, bool globalInternalCacheEnabled=false)
Real3dPoint centroidAsDGtalPoint(const Face f) const
DenseMatrix covariantProjection(const Face f, const Vector &uf)
PolygonalCalculus & operator=(const PolygonalCalculus &other)=delete
Vector vectorArea(const Face f) const
DenseMatrix E(const Face f) const
DenseMatrix blockConnection(const Face &f) const
MySurfaceMesh::Vertex Vertex
Vertex type.
SparseMatrix globalInverseLumpedMassMatrix() const
Vector computeVertexNormal(const Vertex &v) const
std::array< std::unordered_map< Face, DenseMatrix >, 15 > myGlobalCache
size_t faceDegree(Face f) const
Vector toExtrinsicVector(const Vertex v, const Vector &I) const
toExtrinsicVector
static Real3dVector toReal3dVector(const Eigen::Vector3d &x)
toReal3dVector converts Eigen::Vector3d to Real3dVector.
SparseMatrix doubledGlobalLumpedMassMatrix() const
bool myGlobalCacheEnabled
Global cache.
DenseMatrix bracket(const Vector &n) const
PolygonalCalculus(const PolygonalCalculus &other)=delete
DenseMatrix kroneckerWithI2(const DenseMatrix &M) const
DenseMatrix P(const Face f) const
SparseMatrix globalLaplaceBeltrami(const double lambda=1.0) const
PolygonalCalculus(const ConstAlias< MySurfaceMesh > surf, const std::function< Real3dVector(Vertex)> &embedder, bool globalInternalCacheEnabled=false)
DenseMatrix Qvf(const Vertex &v, const Face &f) const
void setEmbedder(const std::function< Real3dPoint(Face, Vertex)> &externalFunctor)
LinAlg::SolverSimplicialLDLT Solver
Type of a sparse matrix solver.
DenseMatrix sharp(const Face f) const
bool checkCache(OPERATOR key, const Face f) const
double faceArea(const Face f) const
DenseMatrix covariantGradient_f(const Face &f) const
PolygonalCalculus()=delete
std::vector< Vector > toExtrinsicVectors(const std::vector< Vector > &I) const
std::vector< Vector > getOperatorCacheVector(const std::function< Vector(Face)> &perFaceVectorOperator) const
DenseMatrix coGradient(const Face f) const
DenseMatrix divergence(const Face f) const
DenseMatrix M(const Face f, const double lambda=1.0) const
DenseMatrix shapeOperator(const Face f) const
~PolygonalCalculus()=default
DenseMatrix connectionLaplacian(const Face &f, double lambda=1.0) const
size_t degree(const Face f) const
void updateFaceDegree()
Update the face degree cache.
SparseMatrix globalLumpedMassMatrix() const
DenseMatrix gradient(const Face f) const
MySurfaceMesh::Face Face
Face type.
PolygonalCalculus< TRealPoint, TRealVector > Self
Self type.
LinAlg::DenseMatrix DenseMatrix
Type of dense matrix.
MySurfaceMesh::RealPoint Real3dPoint
Position type.
static Vector toVector(const Eigen::Vector3d &x)
toVector convert Real3dPoint to Eigen::VectorXd
void setInCache(OPERATOR key, const Face f, const DenseMatrix &ope) const
LinAlg::Triplet Triplet
Type of sparse matrix triplet.
std::function< Real3dVector(Vertex)> myVertexNormalEmbedder
Embedding function (vertex)->R^3 for the vertex normal.
PolygonalCalculus(const ConstAlias< MySurfaceMesh > surf, bool globalInternalCacheEnabled=false)
void disableInternalGlobalCache()
DenseMatrix B(const Face f) const
const MySurfaceMesh * getSurfaceMeshPtr() const
Vector Form
Global 0-form, 1-form, 2-form are Vector.
DenseMatrix covariantProjection_f(const Face &f) const
void enableInternalGlobalCache()
LinAlg::DenseVector Vector
Type of Vector.
SparseMatrix globalConnectionLaplace(const double lambda=1.0) const
static Vector project(const Vector &u, const Vector &n)
DenseMatrix A(const Face f) const
DenseMatrix covariantGradient(const Face f, const Vector &uf)
BOOST_STATIC_ASSERT((dimension==3))
DenseMatrix flat(const Face f) const
SparseMatrix identity0() const
DGtal is the top-level namespace which contains all DGtal functions and types.
std::ostream & operator<<(std::ostream &out, const ClosedIntegerHalfPlane< TSpace > &object)
DGtal::uint32_t Dimension
Aim: Provide linear algebra backend using Eigen dense and sparse matrix as well as dense vector....
Eigen::VectorXd DenseVector
Eigen::SimplicialLDLT< SparseMatrix > SolverSimplicialLDLT
Eigen::Triplet< double, SparseMatrix::StorageIndex > Triplet
Eigen::MatrixXd DenseMatrix
Eigen::SparseMatrix< DenseVector::Scalar, Eigen::ColMajor, DenseVector::Index > SparseMatrix
Aim: Represents an embedded mesh as faces and a list of vertices. Vertices may be shared among faces ...
const Faces & incidentFaces(Vertex v) const
std::size_t Index
The type used for numbering vertices and faces.
const Vertices & neighborVertices(Vertex v) const
const Vertices & incidentVertices(Face f) const
RealPoint & position(Vertex v)
Functor that projects a face vertex of a surface mesh onto the tangent plane given by a per-face norm...
const MySurfaceMesh * mySurfaceMesh
Alias to the surface mesh.
EmbedderFromNormalVectors(ConstAlias< std::vector< Real3dPoint > > normals, ConstAlias< MySurfaceMesh > surfmesh)
const std::vector< Real3dPoint > * myNormals
Alias to the normal vectors.
SurfaceMesh< TRealPoint, TRealVector > MySurfaceMesh
Type of SurfaceMesh.
MySurfaceMesh::RealPoint Real3dPoint
Position type.
MySurfaceMesh::Face Face
Face type.
MySurfaceMesh::Vertex Vertex
Vertex type.
Real3dPoint operator()(const Face &f, const Vertex &v)
EmbedderFromNormalVectors()=delete