DGtal  1.4.beta
DGtal::PolygonalCalculus< TRealPoint, TRealVector > Class Template Reference

Implements differential operators on polygonal surfaces from [44]. More...

#include <DGtal/dec/PolygonalCalculus.h>

Public Types

typedef PolygonalCalculus< TRealPoint, TRealVector > Self
 Self type. More...
 
typedef SurfaceMesh< TRealPoint, TRealVector > MySurfaceMesh
 Type of SurfaceMesh. More...
 
typedef MySurfaceMesh::Vertex Vertex
 Vertex type. More...
 
typedef MySurfaceMesh::Face Face
 Face type. More...
 
typedef MySurfaceMesh::RealPoint Real3dPoint
 Position type. More...
 
typedef MySurfaceMesh::RealVector Real3dVector
 Real vector type. More...
 
typedef EigenLinearAlgebraBackend LinAlg
 Linear Algebra Backend from Eigen. More...
 
typedef LinAlg::DenseVector Vector
 Type of Vector. More...
 
typedef Vector Form
 Global 0-form, 1-form, 2-form are Vector. More...
 
typedef LinAlg::DenseMatrix DenseMatrix
 Type of dense matrix. More...
 
typedef LinAlg::SparseMatrix SparseMatrix
 Type of sparse matrix. More...
 
typedef LinAlg::Triplet Triplet
 Type of sparse matrix triplet. More...
 
typedef LinAlg::SolverSimplicialLDLT Solver
 Type of a sparse matrix solver. More...
 

Public Member Functions

 BOOST_STATIC_ASSERT ((dimension==3))
 
Standard services
 PolygonalCalculus (const ConstAlias< MySurfaceMesh > surf, bool globalInternalCacheEnabled=false)
 
 PolygonalCalculus (const ConstAlias< MySurfaceMesh > surf, const std::function< Real3dPoint(Face, Vertex)> &embedder, bool globalInternalCacheEnabled=false)
 
 PolygonalCalculus (const ConstAlias< MySurfaceMesh > surf, const std::function< Real3dVector(Vertex)> &embedder, bool globalInternalCacheEnabled=false)
 
 PolygonalCalculus (const ConstAlias< MySurfaceMesh > surf, const std::function< Real3dPoint(Face, Vertex)> &pos_embedder, const std::function< Vector(Vertex)> &normal_embedder, bool globalInternalCacheEnabled=false)
 
 PolygonalCalculus ()=delete
 
 ~PolygonalCalculus ()=default
 
 PolygonalCalculus (const PolygonalCalculus &other)=delete
 
 PolygonalCalculus (PolygonalCalculus &&other)=delete
 
PolygonalCalculusoperator= (const PolygonalCalculus &other)=delete
 
PolygonalCalculusoperator= (PolygonalCalculus &&other)=delete
 
Embedding services
void setEmbedder (const std::function< Real3dPoint(Face, Vertex)> &externalFunctor)
 
Per face operators on scalars
DenseMatrix X (const Face f) const
 
DenseMatrix D (const Face f) const
 
DenseMatrix E (const Face f) const
 
DenseMatrix A (const Face f) const
 
Vector vectorArea (const Face f) const
 
double faceArea (const Face f) const
 
Vector faceNormal (const Face f) const
 
Real3dVector faceNormalAsDGtalVector (const Face f) const
 
DenseMatrix coGradient (const Face f) const
 
DenseMatrix bracket (const Vector &n) const
 
DenseMatrix gradient (const Face f) const
 
DenseMatrix flat (const Face f) const
 
DenseMatrix B (const Face f) const
 
Vector centroid (const Face f) const
 
Real3dPoint centroidAsDGtalPoint (const Face f) const
 
DenseMatrix sharp (const Face f) const
 
DenseMatrix P (const Face f) const
 
DenseMatrix M (const Face f, const double lambda=1.0) const
 
DenseMatrix divergence (const Face f) const
 
DenseMatrix curl (const Face f) const
 
DenseMatrix laplaceBeltrami (const Face f, const double lambda=1.0) const
 
Per face operators on vector fields
DenseMatrix Tv (const Vertex &v) const
 
DenseMatrix Tf (const Face &f) const
 
Vector toExtrinsicVector (const Vertex v, const Vector &I) const
 toExtrinsicVector More...
 
std::vector< VectortoExtrinsicVectors (const std::vector< Vector > &I) const
 
DenseMatrix Qvf (const Vertex &v, const Face &f) const
 
DenseMatrix Rvf (const Vertex &v, const Face &f) const
 
DenseMatrix shapeOperator (const Face f) const
 
DenseMatrix transportAndFormatVectorField (const Face f, const Vector &uf)
 
DenseMatrix covariantGradient (const Face f, const Vector &uf)
 
DenseMatrix covariantProjection (const Face f, const Vector &uf)
 
DenseMatrix covariantGradient_f (const Face &f) const
 
DenseMatrix covariantProjection_f (const Face &f) const
 
DenseMatrix connectionLaplacian (const Face &f, double lambda=1.0) const
 
Global operators
Form form0 () const
 
SparseMatrix identity0 () const
 
SparseMatrix globalLaplaceBeltrami (const double lambda=1.0) const
 
SparseMatrix globalLumpedMassMatrix () const
 
SparseMatrix globalInverseLumpedMassMatrix () const
 
SparseMatrix globalConnectionLaplace (const double lambda=1.0) const
 
SparseMatrix doubledGlobalLumpedMassMatrix () const
 
Cache mechanism
std::vector< DenseMatrixgetOperatorCacheMatrix (const std::function< DenseMatrix(Face)> &perFaceOperator) const
 
std::vector< VectorgetOperatorCacheVector (const std::function< Vector(Face)> &perFaceVectorOperator) const
 
void enableInternalGlobalCache ()
 
void disableInternalGlobalCache ()
 
Common services
void init ()
 
size_t faceDegree (Face f) const
 
size_t nbVertices () const
 
size_t nbFaces () const
 
size_t degree (const Face f) const
 
const MySurfaceMeshgetSurfaceMeshPtr () const
 
void selfDisplay (std::ostream &out) const
 
bool isValid () const
 

Static Public Attributes

static const Dimension dimension = TRealPoint::dimension
 Concept checking. More...
 

Private Attributes

const MySurfaceMeshmySurfaceMesh
 Underlying SurfaceMesh. More...
 
std::function< Real3dPoint(Face, Vertex)> myEmbedder
 Embedding function (face,vertex)->R^3 for the vertex position wrt. the face. More...
 
std::function< Real3dVector(Vertex)> myVertexNormalEmbedder
 Embedding function (vertex)->R^3 for the vertex normal. More...
 
std::vector< size_t > myFaceDegree
 Cache containing the face degree. More...
 
bool myGlobalCacheEnabled
 Global cache. More...
 
std::array< std::unordered_map< Face, DenseMatrix >, 15 > myGlobalCache
 

Protected services and types

enum  OPERATOR {
  X_ , D_ , E_ , A_ ,
  COGRAD_ , GRAD_ , FLAT_ , B_ ,
  SHARP_ , P_ , M_ , DIVERGENCE_ ,
  CURL_ , L_ , CON_L_
}
 Enum for operators in the internal cache strategy. More...
 
void updateFaceDegree ()
 Update the face degree cache. More...
 
bool checkCache (OPERATOR key, const Face f) const
 
void setInCache (OPERATOR key, const Face f, const DenseMatrix &ope) const
 
Vector computeVertexNormal (const Vertex &v) const
 
Eigen::Vector3d n_v (const Vertex &v) const
 
DenseMatrix blockConnection (const Face &f) const
 
DenseMatrix kroneckerWithI2 (const DenseMatrix &M) const
 
static Vector project (const Vector &u, const Vector &n)
 
static Vector toVector (const Eigen::Vector3d &x)
 toVector convert Real3dPoint to Eigen::VectorXd More...
 
static Eigen::Vector3d toVec3 (const Real3dPoint &x)
 toVec3 convert Real3dPoint to Eigen::Vector3d More...
 
static Real3dVector toReal3dVector (const Eigen::Vector3d &x)
 toReal3dVector converts Eigen::Vector3d to Real3dVector. More...
 

Detailed Description

template<typename TRealPoint, typename TRealVector>
class DGtal::PolygonalCalculus< TRealPoint, TRealVector >

Implements differential operators on polygonal surfaces from [44].

Description of template class 'PolygonalCalculus'

See Discrete differential calculus on polygonal surfaces for details.

Note
The sign convention for the divergence and the Laplacian operator is opposite to the one of [44]. This is to match the usual mathematical convention that the Laplacian (and the Laplacian-Beltrami) has negative eigenvalues (and is the sum of second derivatives in the cartesian grid). It also follows the formal adjointness of exterior derivative and opposite of divergence as relation \( \langle \mathrm{d} u, v \rangle = - \langle u, \mathrm{div} v \rangle \). See also https://en.wikipedia.org/wiki/Laplace–Beltrami_operator
Template Parameters
TRealPointa model of points \(\mathbb{R}^3\) (e.g. PointVector).
TRealVectora model of vectors in \(\mathbb{R}^3\) (e.g. PointVector).

Definition at line 70 of file PolygonalCalculus.h.

Member Typedef Documentation

◆ DenseMatrix

template<typename TRealPoint , typename TRealVector >
typedef LinAlg::DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::DenseMatrix

Type of dense matrix.

Definition at line 100 of file PolygonalCalculus.h.

◆ Face

template<typename TRealPoint , typename TRealVector >
typedef MySurfaceMesh::Face DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Face

Face type.

Definition at line 87 of file PolygonalCalculus.h.

◆ Form

template<typename TRealPoint , typename TRealVector >
typedef Vector DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Form

Global 0-form, 1-form, 2-form are Vector.

Definition at line 98 of file PolygonalCalculus.h.

◆ LinAlg

template<typename TRealPoint , typename TRealVector >
typedef EigenLinearAlgebraBackend DGtal::PolygonalCalculus< TRealPoint, TRealVector >::LinAlg

Linear Algebra Backend from Eigen.

Definition at line 94 of file PolygonalCalculus.h.

◆ MySurfaceMesh

template<typename TRealPoint , typename TRealVector >
typedef SurfaceMesh<TRealPoint, TRealVector> DGtal::PolygonalCalculus< TRealPoint, TRealVector >::MySurfaceMesh

Type of SurfaceMesh.

Definition at line 83 of file PolygonalCalculus.h.

◆ Real3dPoint

template<typename TRealPoint , typename TRealVector >
typedef MySurfaceMesh::RealPoint DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Real3dPoint

Position type.

Definition at line 89 of file PolygonalCalculus.h.

◆ Real3dVector

template<typename TRealPoint , typename TRealVector >
typedef MySurfaceMesh::RealVector DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Real3dVector

Real vector type.

Definition at line 91 of file PolygonalCalculus.h.

◆ Self

template<typename TRealPoint , typename TRealVector >
typedef PolygonalCalculus<TRealPoint, TRealVector> DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Self

Self type.

Definition at line 80 of file PolygonalCalculus.h.

◆ Solver

template<typename TRealPoint , typename TRealVector >
typedef LinAlg::SolverSimplicialLDLT DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Solver

Type of a sparse matrix solver.

Definition at line 107 of file PolygonalCalculus.h.

◆ SparseMatrix

template<typename TRealPoint , typename TRealVector >
typedef LinAlg::SparseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::SparseMatrix

Type of sparse matrix.

Definition at line 102 of file PolygonalCalculus.h.

◆ Triplet

template<typename TRealPoint , typename TRealVector >
typedef LinAlg::Triplet DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Triplet

Type of sparse matrix triplet.

Definition at line 104 of file PolygonalCalculus.h.

◆ Vector

template<typename TRealPoint , typename TRealVector >
typedef LinAlg::DenseVector DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Vector

Type of Vector.

Definition at line 96 of file PolygonalCalculus.h.

◆ Vertex

template<typename TRealPoint , typename TRealVector >
typedef MySurfaceMesh::Vertex DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Vertex

Vertex type.

Definition at line 85 of file PolygonalCalculus.h.

Member Enumeration Documentation

◆ OPERATOR

template<typename TRealPoint , typename TRealVector >
enum DGtal::PolygonalCalculus::OPERATOR
protected

Constructor & Destructor Documentation

◆ PolygonalCalculus() [1/7]

template<typename TRealPoint , typename TRealVector >
DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus ( const ConstAlias< MySurfaceMesh surf,
bool  globalInternalCacheEnabled = false 
)
inline

Create a Polygonal DEC structure from a surface mesh (surf) using an default identity embedder.

Parameters
surfan instance of SurfaceMesh
globalInternalCacheEnabledenable the internal cache for all operators (default: false)

Definition at line 116 of file PolygonalCalculus.h.

117  :
118  mySurfaceMesh(&surf), myGlobalCacheEnabled(globalInternalCacheEnabled)
119  {
120  myEmbedder =[&](Face f,Vertex v){ (void)f; return mySurfaceMesh->position(v);};
122  init();
123  };
std::function< Real3dPoint(Face, Vertex)> myEmbedder
Embedding function (face,vertex)->R^3 for the vertex position wrt. the face.
const MySurfaceMesh * mySurfaceMesh
Underlying SurfaceMesh.
Vector computeVertexNormal(const Vertex &v) const
static Real3dVector toReal3dVector(const Eigen::Vector3d &x)
toReal3dVector converts Eigen::Vector3d to Real3dVector.
bool myGlobalCacheEnabled
Global cache.
std::function< Real3dVector(Vertex)> myVertexNormalEmbedder
Embedding function (vertex)->R^3 for the vertex normal.
RealPoint & position(Vertex v)
Definition: SurfaceMesh.h:647
TriMesh::Face Face
TriMesh::Vertex Vertex

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::computeVertexNormal(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::init(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myEmbedder, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myVertexNormalEmbedder, DGtal::SurfaceMesh< TRealPoint, TRealVector >::position(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::toReal3dVector().

◆ PolygonalCalculus() [2/7]

template<typename TRealPoint , typename TRealVector >
DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus ( const ConstAlias< MySurfaceMesh surf,
const std::function< Real3dPoint(Face, Vertex)> &  embedder,
bool  globalInternalCacheEnabled = false 
)
inline

Create a Polygonal DEC structure from a surface mesh (surf) and an embedder for the vertex position: function with two parameters, a face and a vertex which outputs the embedding in R^3 of the vertex w.r.t. to the face.

Parameters
surfan instance of SurfaceMesh
embedderan embedder for the vertex position
globalInternalCacheEnabledif true, the global operator cache is enabled

Definition at line 131 of file PolygonalCalculus.h.

133  :
134  mySurfaceMesh(&surf), myEmbedder(embedder), myGlobalCacheEnabled(globalInternalCacheEnabled)
135  {
137  init();
138  };

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::computeVertexNormal(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::init(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myVertexNormalEmbedder, and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::toReal3dVector().

◆ PolygonalCalculus() [3/7]

template<typename TRealPoint , typename TRealVector >
DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus ( const ConstAlias< MySurfaceMesh surf,
const std::function< Real3dVector(Vertex)> &  embedder,
bool  globalInternalCacheEnabled = false 
)
inline

Create a Polygonal DEC structure from a surface mesh (surf) and an embedder for the vertex normal: function with a vertex as parameter which outputs the embedding in R^3 of the vertex normal.

Parameters
surfan instance of SurfaceMesh
embedderan embedder for the vertex position
globalInternalCacheEnabledif true, the global operator cache is enabled

Definition at line 146 of file PolygonalCalculus.h.

148  :
149  mySurfaceMesh(&surf), myVertexNormalEmbedder(embedder),
150  myGlobalCacheEnabled(globalInternalCacheEnabled)
151  {
152  myEmbedder = [&](Face f,Vertex v){(void)f; return mySurfaceMesh->position(v); };
153  init();
154  };

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::init(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myEmbedder, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh, and DGtal::SurfaceMesh< TRealPoint, TRealVector >::position().

◆ PolygonalCalculus() [4/7]

template<typename TRealPoint , typename TRealVector >
DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus ( const ConstAlias< MySurfaceMesh surf,
const std::function< Real3dPoint(Face, Vertex)> &  pos_embedder,
const std::function< Vector(Vertex)> &  normal_embedder,
bool  globalInternalCacheEnabled = false 
)
inline

Create a Polygonal DEC structure from a surface mesh (surf) and an embedder for the vertex position: function with two parameters, a face and a vertex which outputs the embedding in R^3 of the vertex w.r.t. to the face. and an embedder for the vertex normal: function with a vertex as parameter which outputs the embedding in R^3 of the vertex normal.

Parameters
surfan instance of SurfaceMesh
pos_embedderan embedder for the position
normal_embedderan embedder for the position
globalInternalCacheEnabled

Definition at line 166 of file PolygonalCalculus.h.

169  :
170  mySurfaceMesh(&surf), myEmbedder(pos_embedder),
171  myVertexNormalEmbedder(normal_embedder),
172  myGlobalCacheEnabled(globalInternalCacheEnabled)
173  {
174  init();
175  };

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::init().

◆ PolygonalCalculus() [5/7]

template<typename TRealPoint , typename TRealVector >
DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus ( )
delete

Deleted default constructor.

◆ ~PolygonalCalculus()

template<typename TRealPoint , typename TRealVector >
DGtal::PolygonalCalculus< TRealPoint, TRealVector >::~PolygonalCalculus ( )
default

Destructor (default).

◆ PolygonalCalculus() [6/7]

template<typename TRealPoint , typename TRealVector >
DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus ( const PolygonalCalculus< TRealPoint, TRealVector > &  other)
delete

Deleted copy constructor.

Parameters
otherthe object to clone.

◆ PolygonalCalculus() [7/7]

template<typename TRealPoint , typename TRealVector >
DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus ( PolygonalCalculus< TRealPoint, TRealVector > &&  other)
delete

Deleted move constructor.

Parameters
otherthe object to move.

Member Function Documentation

◆ A()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::A ( const Face  f) const
inline

Average operator to average, per edge, its vertex values.

Parameters
fthe face
Returns
a degree x degree matrix

Definition at line 295 of file PolygonalCalculus.h.

296  {
297  if (checkCache(A_,f))
298  return myGlobalCache[A_][f];
299 
300  const auto nf = myFaceDegree[f];
301  DenseMatrix a = DenseMatrix::Zero(nf ,nf);
302  for(auto i=0u; i < nf; ++i)
303  {
304  a(i, (i+1)%nf) = 0.5;
305  a(i,i) = 0.5;
306  }
307 
308  setInCache(A_,f,a);
309  return a;
310  }
std::vector< size_t > myFaceDegree
Cache containing the face degree.
std::array< std::unordered_map< Face, DenseMatrix >, 15 > myGlobalCache
bool checkCache(OPERATOR key, const Face f) const
void setInCache(OPERATOR key, const Face f, const DenseMatrix &ope) const
EigenLinearAlgebraBackend::DenseMatrix DenseMatrix

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::A_, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::checkCache(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myFaceDegree, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCache, and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::setInCache().

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::B(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::coGradient(), and TEST_CASE().

◆ B()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::B ( const Face  f) const
inline

◆ blockConnection()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::blockConnection ( const Face f) const
inlineprotected
Returns
Block Diagonal matrix with Rvf for each vertex v in face f

Definition at line 1182 of file PolygonalCalculus.h.

1183  {
1184  auto nf = degree(f);
1185  DenseMatrix RU_fO = DenseMatrix::Zero(nf * 2,nf * 2);
1186  size_t cpt = 0;
1187  for (auto v : getSurfaceMeshPtr()->incidentVertices(f))
1188  {
1189  auto Rv = Rvf(v,f);
1190  RU_fO.block<2,2>(2 * cpt,2 * cpt) = Rv;
1191  ++cpt;
1192  }
1193  return RU_fO;
1194  }
DenseMatrix Rvf(const Vertex &v, const Face &f) const
const MySurfaceMesh * getSurfaceMeshPtr() const
size_t degree(const Face f) const

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::degree(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::getSurfaceMeshPtr(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Rvf().

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantGradient_f(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantProjection_f().

◆ BOOST_STATIC_ASSERT()

template<typename TRealPoint , typename TRealVector >
DGtal::PolygonalCalculus< TRealPoint, TRealVector >::BOOST_STATIC_ASSERT ( (dimension==3)  )

◆ bracket()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::bracket ( const Vector n) const
inline

Return [n] as the 3x3 operator such that [n]q = n x q

Parameters
na vector

Definition at line 378 of file PolygonalCalculus.h.

379  {
380  DenseMatrix brack(3,3);
381  brack << 0.0 , -n(2), n(1),
382  n(2), 0.0 , -n(0),
383  -n(1) , n(0),0.0 ;
384  return brack;
385  }

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::gradient(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Qvf(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::sharp().

◆ centroid()

template<typename TRealPoint , typename TRealVector >
Vector DGtal::PolygonalCalculus< TRealPoint, TRealVector >::centroid ( const Face  f) const
inline
Returns
the centroid of the face
Parameters
fthe face

Definition at line 428 of file PolygonalCalculus.h.

429  {
430  const auto nf = myFaceDegree[f];
431  return 1.0/(double)nf * X(f).transpose() * Vector::Ones(nf);
432  }

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myFaceDegree, and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::X().

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::centroidAsDGtalPoint(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::sharp(), and TEST_CASE().

◆ centroidAsDGtalPoint()

template<typename TRealPoint , typename TRealVector >
Real3dPoint DGtal::PolygonalCalculus< TRealPoint, TRealVector >::centroidAsDGtalPoint ( const Face  f) const
inline
Returns
the centroid of the face as a DGtal RealPoint
Parameters
fthe face

Definition at line 436 of file PolygonalCalculus.h.

437  {
438  const Vector c = centroid(f);
439  return {c(0),c(1),c(2)};
440  }
Vector centroid(const Face f) const
DigitalPlane::Point Vector

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::centroid().

Referenced by TEST_CASE().

◆ checkCache()

template<typename TRealPoint , typename TRealVector >
bool DGtal::PolygonalCalculus< TRealPoint, TRealVector >::checkCache ( OPERATOR  key,
const Face  f 
) const
inlineprotected

◆ coGradient()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::coGradient ( const Face  f) const
inline

◆ computeVertexNormal()

template<typename TRealPoint , typename TRealVector >
Vector DGtal::PolygonalCalculus< TRealPoint, TRealVector >::computeVertexNormal ( const Vertex v) const
inlineprotected

Compute the (normalized) normal vector at a Vertex by averaging the adjacent face normal vectors.

Parameters
vthe vertex to compute the normal from
Returns
3D normal vector at vertex v

Definition at line 1149 of file PolygonalCalculus.h.

1150  {
1151  Vector n(3);
1152  n(0) = 0.;
1153  n(1) = 0.;
1154  n(2) = 0.;
1155  /* for (auto f : mySurfaceMesh->incidentFaces(v))
1156  n += vectorArea(f);
1157  return n.normalized();
1158  */
1159  auto faces = mySurfaceMesh->incidentFaces(v);
1160  for (auto f : faces)
1161  n += vectorArea(f);
1162 
1163  if (fabs(n.norm() - 0.0) < 0.00001)
1164  {
1165  //On non-manifold edges touching the boundary, n may be null.
1166  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;
1167  n << Vector::Random(3);
1168  }
1169  n = n.normalized();
1170  return n;
1171  }
Vector vectorArea(const Face f) const
std::ostream & warning()
Trace trace
Definition: Common.h:153
const Faces & incidentFaces(Vertex v) const
Definition: SurfaceMesh.h:321

References DGtal::SurfaceMesh< TRealPoint, TRealVector >::incidentFaces(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh, DGtal::trace, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::vectorArea(), and DGtal::Trace::warning().

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus().

◆ connectionLaplacian()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::connectionLaplacian ( const Face f,
double  lambda = 1.0 
) const
inline

L∇ := -(afG∇tG∇+λP∇tP∇)

Returns
Connection/Vector laplacian at face f
Note
The sign convention for the divergence and the Laplacian operator is opposite to the one of [44]. This is to match the usual mathematical convention that the laplacian (and the Laplacian-Beltrami) has negative eigenvalues (and is the sum of second derivatives in the cartesian grid). It also follows the formal adjointness of exterior derivative and opposite of divergence as relation \( \langle \mathrm{d} u, v \rangle = - \langle u, \mathrm{div} v \rangle \). See also https://en.wikipedia.org/wiki/Laplace–Beltrami_operator

Definition at line 752 of file PolygonalCalculus.h.

753  {
754  if (checkCache(CON_L_,f))
755  return myGlobalCache[CON_L_][f];
758  DenseMatrix L = -(faceArea(f) * G.transpose() * G + lambda * P.transpose() * P);
759  setInCache(CON_L_,f,L);
760  return L;
761  }
DenseMatrix P(const Face f) const
double faceArea(const Face f) const
DenseMatrix covariantGradient_f(const Face &f) const
DenseMatrix covariantProjection_f(const Face &f) const

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::checkCache(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::CON_L_, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantGradient_f(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantProjection_f(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceArea(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCache, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::P(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::setInCache().

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::globalConnectionLaplace(), and TEST_CASE().

◆ covariantGradient()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantGradient ( const Face  f,
const Vector uf 
)
inline

Covarient gradient at a face a f of intrinsic vectors uf.

Parameters
uflist of all intrinsic vectors per vertex concatenated in a column vector
fthe face
Returns
the covariant gradient of the given vector field uf (expressed in corresponding vertex tangent frames), wrt face f
Note
Unlike the rest of the per face operators, the covariant operators need to be applied directly to the restriction of the vector field to the face,

Definition at line 698 of file PolygonalCalculus.h.

699  {
700  return Tf(f).transpose() * gradient(f) *
702  }
DenseMatrix transportAndFormatVectorField(const Face f, const Vector &uf)
DenseMatrix Tf(const Face &f) const
DenseMatrix gradient(const Face f) const

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::gradient(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Tf(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::transportAndFormatVectorField().

Referenced by TEST_CASE().

◆ covariantGradient_f()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantGradient_f ( const Face f) const
inline
Returns
Covariant Gradient Operator, returns the operator that acts on the concatenated vectors. When applied, gives the associated 2x2 matrix in the isomorphic vector form (a b c d)^t to be used in the dirichlet energy (vector laplacian) G∇_f. Used to define the connection Laplacian.

Definition at line 724 of file PolygonalCalculus.h.

725  {
726  return kroneckerWithI2(Tf(f).transpose() * gradient(f)) * blockConnection(f);
727  }
DenseMatrix blockConnection(const Face &f) const
DenseMatrix kroneckerWithI2(const DenseMatrix &M) const

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::blockConnection(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::gradient(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::kroneckerWithI2(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Tf().

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::connectionLaplacian().

◆ covariantProjection()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantProjection ( const Face  f,
const Vector uf 
)
inline

Compute the covariance projection at a face f of intrinsic vectors uf.

Parameters
uflist of all intrinsic vectors per vertex concatenated in a column vector
fthe face
Returns
the covariant projection of the given vector field uf ( restricted to face f and expressed in corresponding vertex tangent frames)
Note
Unlike the rest of the per face operators, the covariant operators need to be applied directly to the restriction of the vector field to the face

Definition at line 714 of file PolygonalCalculus.h.

715  {
716  return P(f) * D(f) * transportAndFormatVectorField(f,uf);
717  }
DenseMatrix D(const Face f) const

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::D(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::P(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::transportAndFormatVectorField().

Referenced by TEST_CASE().

◆ covariantProjection_f()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantProjection_f ( const Face f) const
inline
Returns
Projection Gradient Operator, returns the operator that acts on the concatenated vectors. When applied, gives the associated nfx2 matrix in the isomorphic vector form (a b c d ...)^t to be used in the dirichlet energy (vector laplacian) P∇_f. Used to define the connection Laplacian.

Definition at line 734 of file PolygonalCalculus.h.

735  {
736  return kroneckerWithI2(P(f) * D(f)) * blockConnection(f);
737  ;
738  }

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::blockConnection(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::D(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::kroneckerWithI2(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::P().

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::connectionLaplacian().

◆ curl()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::curl ( const Face  f) const
inline

Curl operator of a one-form (identity matrix).

Parameters
fthe face
Returns
a degree x degree matrix

Definition at line 518 of file PolygonalCalculus.h.

519  {
520  if (checkCache(CURL_,f))
521  return myGlobalCache[CURL_][f];
522 
523  DenseMatrix op = DenseMatrix::Identity(myFaceDegree[f],myFaceDegree[f]);
524 
525  setInCache(CURL_,f,op);
526  return op;
527  }

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::checkCache(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::CURL_, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myFaceDegree, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCache, and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::setInCache().

Referenced by TEST_CASE().

◆ D()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::D ( const Face  f) const
inline

◆ degree()

template<typename TRealPoint , typename TRealVector >
size_t DGtal::PolygonalCalculus< TRealPoint, TRealVector >::degree ( const Face  f) const
inline
Returns
the degree of the face f (number of vertices)
Parameters
fthe face

Definition at line 1024 of file PolygonalCalculus.h.

1025  {
1026  return myFaceDegree[f];
1027  }

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myFaceDegree.

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::blockConnection(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::globalConnectionLaplace().

◆ disableInternalGlobalCache()

template<typename TRealPoint , typename TRealVector >
void DGtal::PolygonalCalculus< TRealPoint, TRealVector >::disableInternalGlobalCache ( )
inline

Disable the internal global cache for operators. This method will also clean up the

Definition at line 982 of file PolygonalCalculus.h.

983  {
984  myGlobalCacheEnabled = false;
985  myGlobalCache.clear();
986  }

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCache, and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCacheEnabled.

◆ divergence()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::divergence ( const Face  f) const
inline

Divergence operator of a one-form.

Parameters
fthe face
Returns
a degree x degree matrix
Note
The sign convention for the divergence and the Laplacian operator is opposite to the one of [44] . This is to match the usual mathematical convention that the Laplacian (and the Laplacian-Beltrami) has negative eigenvalues (and is the sum of second derivatives in the cartesian grid). It also follows the formal adjointness of exterior derivative and opposite of divergence as relation \( \langle \mathrm{d} u, v \rangle = - \langle u, \mathrm{div} v \rangle \). See also https://en.wikipedia.org/wiki/Laplace–Beltrami_operator

Definition at line 504 of file PolygonalCalculus.h.

505  {
506  if (checkCache(DIVERGENCE_,f))
507  return myGlobalCache[DIVERGENCE_][f];
508 
509  DenseMatrix op = -1.0 * D(f).transpose() * M(f);
510  setInCache(DIVERGENCE_,f,op);
511 
512  return op;
513  }
DenseMatrix M(const Face f, const double lambda=1.0) const

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::checkCache(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::D(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::DIVERGENCE_, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::M(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCache, and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::setInCache().

◆ doubledGlobalLumpedMassMatrix()

template<typename TRealPoint , typename TRealVector >
SparseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::doubledGlobalLumpedMassMatrix ( ) const
inline

Compute and returns the global lumped mass matrix tensorized with Id_2 (used for connection laplacian) (diagonal matrix with Max's weights for each vertex). M(2*i,2*i) = ∑_{adjface f} faceArea(f)/degree(f) ; M(2*i+1,2*i+1) = M(2*i,2*i)

Returns
the global lumped mass matrix.

Definition at line 902 of file PolygonalCalculus.h.

903  {
904  auto nv = mySurfaceMesh->nbVertices();
905  SparseMatrix M(2 * nv, 2 * nv);
906  std::vector<Triplet> triplets;
907  for (typename MySurfaceMesh::Index v = 0; v < mySurfaceMesh->nbVertices(); ++v)
908  {
909  auto faces = mySurfaceMesh->incidentFaces(v);
910  auto varea = 0.0;
911  for (auto f : faces)
912  varea += faceArea(f) / (double)myFaceDegree[f];
913  triplets.emplace_back(Triplet(2 * v, 2 * v, varea));
914  triplets.emplace_back(Triplet(2 * v + 1, 2 * v + 1, varea));
915  }
916  M.setFromTriplets(triplets.begin(), triplets.end());
917  return M;
918  }
LinAlg::Triplet Triplet
Type of sparse matrix triplet.
std::size_t Index
The type used for numbering vertices and faces.
Definition: SurfaceMesh.h:105
Size nbVertices() const
Definition: SurfaceMesh.h:288
EigenLinearAlgebraBackend::SparseMatrix SparseMatrix

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceArea(), DGtal::SurfaceMesh< TRealPoint, TRealVector >::incidentFaces(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::M(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myFaceDegree, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh, and DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbVertices().

◆ E()

◆ enableInternalGlobalCache()

template<typename TRealPoint , typename TRealVector >
void DGtal::PolygonalCalculus< TRealPoint, TRealVector >::enableInternalGlobalCache ( )
inline

Enable the internal global cache for operators.

Definition at line 975 of file PolygonalCalculus.h.

976  {
977  myGlobalCacheEnabled = true;
978  }

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCacheEnabled.

◆ faceArea()

◆ faceDegree()

template<typename TRealPoint , typename TRealVector >
size_t DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceDegree ( Face  f) const
inline

Helper to retrieve the degree of the face from the cache.

Parameters
fthe face
Returns
the number of vertices of the face.

Definition at line 1005 of file PolygonalCalculus.h.

1006  {
1007  return myFaceDegree[f];
1008  }

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myFaceDegree.

Referenced by TEST_CASE().

◆ faceNormal()

template<typename TRealPoint , typename TRealVector >
Vector DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceNormal ( const Face  f) const
inline

◆ faceNormalAsDGtalVector()

template<typename TRealPoint , typename TRealVector >
Real3dVector DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceNormalAsDGtalVector ( const Face  f) const
inline

Corrected normal vector of a face.

Parameters
fthe face
Returns
a vector (DGtal RealVector/RealPoint)

Definition at line 358 of file PolygonalCalculus.h.

359  {
360  Vector v = faceNormal(f);
361  return {v(0),v(1),v(2)};
362  }
Vector faceNormal(const Face f) const

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceNormal().

Referenced by TEST_CASE().

◆ flat()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::flat ( const Face  f) const
inline

◆ form0()

template<typename TRealPoint , typename TRealVector >
Form DGtal::PolygonalCalculus< TRealPoint, TRealVector >::form0 ( ) const
inline
Returns
a 0-form initialized to zero

Definition at line 770 of file PolygonalCalculus.h.

771  {
772  return Form::Zero( nbVertices() );
773  }

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::nbVertices().

Referenced by TEST_CASE().

◆ getOperatorCacheMatrix()

template<typename TRealPoint , typename TRealVector >
std::vector<DenseMatrix> DGtal::PolygonalCalculus< TRealPoint, TRealVector >::getOperatorCacheMatrix ( const std::function< DenseMatrix(Face)> &  perFaceOperator) const
inline

Generic method to compute all the per face DenseMatrices and store them in an indexed container.

Usage example:

auto opM = [&](const PolygonalCalculus<Mesh>::Face f){ return calculus.M(f);};
auto cacheM = boxCalculus.getOperatorCacheMatrix(opM);
...
//Now you have access to the cached values and mixed them with un-cached ones
Face f = ...;
auto res = cacheM[f] * calculus.D(f) * phi;
...
MySurfaceMesh::Face Face
Face type.
Parameters
perFaceOperatorthe per face operator
Returns
an indexed container of all DenseMatrix operators (indexed per Face).

Definition at line 941 of file PolygonalCalculus.h.

942  {
943  std::vector<DenseMatrix> cache;
944  for( typename MySurfaceMesh::Index f=0; f < mySurfaceMesh->nbFaces(); ++f)
945  cache.push_back(perFaceOperator(f));
946  return cache;
947  }
Size nbFaces() const
Definition: SurfaceMesh.h:296

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh, and DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbFaces().

Referenced by TEST_CASE().

◆ getOperatorCacheVector()

template<typename TRealPoint , typename TRealVector >
std::vector<Vector> DGtal::PolygonalCalculus< TRealPoint, TRealVector >::getOperatorCacheVector ( const std::function< Vector(Face)> &  perFaceVectorOperator) const
inline

Generic method to compute all the per face Vector and store them in an indexed container.

Usage example:

auto opCentroid = [&](const PolygonalCalculus<Mesh>::Face f){ return calculus.centroid(f);};
auto cacheCentroid = boxCalculus.getOperatorCacheVector(opCentroid);
...
//Now you have access to the cached values and mixed them with un-cached ones
Face f = ...;
auto res = calculus.P(f) * cacheCentroid[f] ;
...
Parameters
perFaceVectorOperatorthe per face operator
Returns
an indexed container of all Vector quantities (indexed per Face).

Definition at line 965 of file PolygonalCalculus.h.

966  {
967  std::vector<Vector> cache;
968  for( typename MySurfaceMesh::Index f=0; f < mySurfaceMesh->nbFaces(); ++f)
969  cache.push_back(perFaceVectorOperator(f));
970  return cache;
971  }

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh, and DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbFaces().

Referenced by TEST_CASE().

◆ getSurfaceMeshPtr()

template<typename TRealPoint , typename TRealVector >
const MySurfaceMesh* DGtal::PolygonalCalculus< TRealPoint, TRealVector >::getSurfaceMeshPtr ( ) const
inline

◆ globalConnectionLaplace()

template<typename TRealPoint , typename TRealVector >
SparseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::globalConnectionLaplace ( const double  lambda = 1.0) const
inline

Computes the global Connection-Laplace-Beltrami operator by accumulating the per face operators.

Parameters
lambdathe regualrization parameter for the local Connection-Laplace-Beltrami operators
Returns
a sparse 2*nbVertices x 2*nbVertices matrix
Note
The sign convention for the divergence is opposite to the one of [44]. This is also true for the Laplacian operator. This is to match the usual mathematical convention that the Laplacian (and the Laplacian-Beltrami) has negative eigenvalues (and is the sum of second derivatives in the cartesian grid). It also follows the formal adjointness of exterior derivative and opposite of divergence as relation \( \langle \mathrm{d} u, v \rangle = - \langle u, \mathrm{div} v \rangle \). See also https://en.wikipedia.org/wiki/Laplace–Beltrami_operator

Definition at line 870 of file PolygonalCalculus.h.

871  {
872  auto nv = mySurfaceMesh->nbVertices();
873  SparseMatrix lapGlobal(2 * nv, 2 * nv);
874  SparseMatrix local(2 * nv, 2 * nv);
875  std::vector<Triplet> triplets;
876  for (typename MySurfaceMesh::Index f = 0; f < mySurfaceMesh->nbFaces(); f++)
877  {
878  auto nf = degree(f);
879  DenseMatrix Lap = connectionLaplacian(f,lambda);
880  const auto vertices = mySurfaceMesh->incidentVertices(f);
881  for (auto i = 0u; i < nf; ++i)
882  for (auto j = 0u; j < nf; ++j)
883  for (short k1 = 0; k1 < 2; k1++)
884  for (short k2 = 0; k2 < 2; k2++)
885  {
886  auto v = Lap(2 * i + k1, 2 * j + k2);
887  if (v != 0.0)
888  triplets.emplace_back(Triplet(2 * vertices[i] + k1,
889  2 * vertices[j] + k2, v));
890  }
891  }
892  lapGlobal.setFromTriplets(triplets.begin(), triplets.end());
893  return lapGlobal;
894  }
DenseMatrix connectionLaplacian(const Face &f, double lambda=1.0) const
std::pair< typename graph_traits< DGtal::DigitalSurface< TDigitalSurfaceContainer > >::vertex_iterator, typename graph_traits< DGtal::DigitalSurface< TDigitalSurfaceContainer > >::vertex_iterator > vertices(const DGtal::DigitalSurface< TDigitalSurfaceContainer > &digSurf)
const Vertices & incidentVertices(Face f) const
Definition: SurfaceMesh.h:315

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::connectionLaplacian(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::degree(), DGtal::SurfaceMesh< TRealPoint, TRealVector >::incidentVertices(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh, DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbFaces(), and DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbVertices().

◆ globalInverseLumpedMassMatrix()

template<typename TRealPoint , typename TRealVector >
SparseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::globalInverseLumpedMassMatrix ( ) const
inline

Compute and returns the inverse of the global lumped mass matrix (diagonal matrix with Max's weights for each vertex).

Returns
the inverse of the global lumped mass matrix.

Definition at line 844 of file PolygonalCalculus.h.

845  {
847  for ( int k = 0; k < iM0.outerSize(); ++k )
848  for ( typename SparseMatrix::InnerIterator it( iM0, k ); it; ++it )
849  it.valueRef() = 1.0 / it.value();
850  return iM0;
851  }
SparseMatrix globalLumpedMassMatrix() const

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::globalLumpedMassMatrix().

◆ globalLaplaceBeltrami()

template<typename TRealPoint , typename TRealVector >
SparseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::globalLaplaceBeltrami ( const double  lambda = 1.0) const
inline

Computes the global Laplace-Beltrami operator by assembling the per face operators.

Parameters
lambdathe regularization parameter for the local Laplace-Beltrami operators
Returns
a sparse nbVertices x nbVertices matrix
Note
The sign convention for the divergence is opposite to the one of [44]. This is also true for the Laplacian operator. This is to match the usual mathematical convention that the Laplacian (and the Laplacian-Beltrami) has negative eigenvalues (and is the sum of second derivatives in the cartesian grid). It also follows the formal adjointness of exterior derivative and opposite of divergence as relation \( \langle \mathrm{d} u, v \rangle = - \langle u, \mathrm{div} v \rangle \). See also https://en.wikipedia.org/wiki/Laplace–Beltrami_operator

Definition at line 797 of file PolygonalCalculus.h.

798  {
800  std::vector<Triplet> triplets;
801  for( typename MySurfaceMesh::Index f = 0; f < mySurfaceMesh->nbFaces(); ++f )
802  {
803  auto nf = myFaceDegree[f];
804  DenseMatrix Lap = this->laplaceBeltrami(f,lambda);
805  const auto vertices = mySurfaceMesh->incidentVertices(f);
806  for(auto i=0u; i < nf; ++i)
807  for(auto j=0u; j < nf; ++j)
808  {
809  auto v = Lap(i,j);
810  if (v!= 0.0)
811  triplets.emplace_back( Triplet( (SparseMatrix::StorageIndex)vertices[ i ], (SparseMatrix::StorageIndex)vertices[ j ],
812  Lap( i, j ) ) );
813  }
814  }
815  lapGlobal.setFromTriplets(triplets.begin(), triplets.end());
816  return lapGlobal;
817  }
DenseMatrix laplaceBeltrami(const Face f, const double lambda=1.0) const

References DGtal::SurfaceMesh< TRealPoint, TRealVector >::incidentVertices(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::laplaceBeltrami(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myFaceDegree, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh, DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbFaces(), and DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbVertices().

Referenced by TEST_CASE().

◆ globalLumpedMassMatrix()

template<typename TRealPoint , typename TRealVector >
SparseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::globalLumpedMassMatrix ( ) const
inline

Compute and returns the global lumped mass matrix (diagonal matrix with Max's weights for each vertex). M(i,i) = ∑_{adjface f} faceArea(f)/degree(f) ;

Returns
the global lumped mass matrix.

Definition at line 824 of file PolygonalCalculus.h.

825  {
827  std::vector<Triplet> triplets;
828  for ( typename MySurfaceMesh::Index v = 0; v < mySurfaceMesh->nbVertices(); ++v )
829  {
830  auto faces = mySurfaceMesh->incidentFaces(v);
831  auto varea = 0.0;
832  for(auto f: faces)
833  varea += faceArea(f) /(double)myFaceDegree[f];
834  triplets.emplace_back(Triplet(v,v,varea));
835  }
836  M.setFromTriplets(triplets.begin(),triplets.end());
837  return M;
838  }

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceArea(), DGtal::SurfaceMesh< TRealPoint, TRealVector >::incidentFaces(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::M(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myFaceDegree, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh, and DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbVertices().

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::globalInverseLumpedMassMatrix(), and TEST_CASE().

◆ gradient()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::gradient ( const Face  f) const
inline

◆ identity0()

template<typename TRealPoint , typename TRealVector >
SparseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::identity0 ( ) const
inline
Returns
the identity linear operator for 0-forms

Definition at line 775 of file PolygonalCalculus.h.

776  {
777  SparseMatrix Id0( nbVertices(), nbVertices() );
778  Id0.setIdentity();
779  return Id0;
780  }

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::nbVertices().

◆ init()

template<typename TRealPoint , typename TRealVector >
void DGtal::PolygonalCalculus< TRealPoint, TRealVector >::init ( )
inline

Update the internal cache structures (e.g. degree of each face).

Definition at line 997 of file PolygonalCalculus.h.

998  {
1000  }
void updateFaceDegree()
Update the face degree cache.

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::updateFaceDegree().

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus().

◆ isValid()

template<typename TRealPoint , typename TRealVector >
bool DGtal::PolygonalCalculus< TRealPoint, TRealVector >::isValid ( ) const
inline

Checks the validity/consistency of the object.

Returns
'true' if the object is valid, 'false' otherwise.

Definition at line 1053 of file PolygonalCalculus.h.

1054  {
1055  return true;
1056  }

Referenced by TEST_CASE().

◆ kroneckerWithI2()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::kroneckerWithI2 ( const DenseMatrix M) const
inlineprotected
Returns
the tensor-kronecker product of M with 2x2 identity matrix

Definition at line 1197 of file PolygonalCalculus.h.

1198  {
1199  size_t h = M.rows();
1200  size_t w = M.cols();
1201  DenseMatrix MK = DenseMatrix::Zero(h * 2,w * 2);
1202  for (size_t j = 0; j < h; j++)
1203  for (size_t i = 0; i < w; i++)
1204  {
1205  MK(2 * j, 2 * i) = M(j, i);
1206  MK(2 * j + 1, 2 * i + 1) = M(j, i);
1207  }
1208  return MK;
1209  }

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::M().

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantGradient_f(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantProjection_f().

◆ laplaceBeltrami()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::laplaceBeltrami ( const Face  f,
const double  lambda = 1.0 
) const
inline

(weak) Laplace-Beltrami operator for the face.

Parameters
fthe face
lambdathe regularization parameter
Returns
a degree x degree matrix
Note
The sign convention for the divergence and the Laplacian operator is opposite to the one of [44] . This is to match the usual mathematical convention that the Laplacian (and the Laplacian-Beltrami) has negative eigenvalues (and is the sum of second derivatives in the cartesian grid). It also follows the formal adjointness of exterior derivative and opposite of divergence as relation \( \langle \mathrm{d} u, v \rangle = - \langle u, \mathrm{div} v \rangle \). See also https://en.wikipedia.org/wiki/Laplace–Beltrami_operator

Definition at line 544 of file PolygonalCalculus.h.

545  {
546  if (checkCache(L_,f))
547  return myGlobalCache[L_][f];
548 
549  DenseMatrix Df = D(f);
550  // Laplacian is a negative operator.
551  DenseMatrix op = -1.0 * Df.transpose() * M(f,lambda) * Df;
552 
553  setInCache(L_, f, op);
554  return op;
555  }

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::checkCache(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::D(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::L_, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::M(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCache, and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::setInCache().

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::globalLaplaceBeltrami(), and TEST_CASE().

◆ M()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::M ( const Face  f,
const double  lambda = 1.0 
) const
inline

◆ n_v()

template<typename TRealPoint , typename TRealVector >
Eigen::Vector3d DGtal::PolygonalCalculus< TRealPoint, TRealVector >::n_v ( const Vertex v) const
inlineprotected
Returns
the normal vector at vertex v, if no normal vertex embedder is set, the normal will be computed

Definition at line 1175 of file PolygonalCalculus.h.

1176  {
1177  return toVec3(myVertexNormalEmbedder(v));
1178  }
static Eigen::Vector3d toVec3(const Real3dPoint &x)
toVec3 convert Real3dPoint to Eigen::Vector3d

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myVertexNormalEmbedder, and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::toVec3().

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Qvf(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::shapeOperator(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Tv().

◆ nbFaces()

template<typename TRealPoint , typename TRealVector >
size_t DGtal::PolygonalCalculus< TRealPoint, TRealVector >::nbFaces ( ) const
inline
Returns
the number of faces of the underlying surface mesh.

Definition at line 1017 of file PolygonalCalculus.h.

1018  {
1019  return mySurfaceMesh->nbFaces();
1020  }

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh, and DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbFaces().

◆ nbVertices()

template<typename TRealPoint , typename TRealVector >
size_t DGtal::PolygonalCalculus< TRealPoint, TRealVector >::nbVertices ( ) const
inline

◆ operator=() [1/2]

template<typename TRealPoint , typename TRealVector >
PolygonalCalculus& DGtal::PolygonalCalculus< TRealPoint, TRealVector >::operator= ( const PolygonalCalculus< TRealPoint, TRealVector > &  other)
delete

Deleted copy assignment operator.

Parameters
otherthe object to copy.
Returns
a reference on 'this'.

◆ operator=() [2/2]

template<typename TRealPoint , typename TRealVector >
PolygonalCalculus& DGtal::PolygonalCalculus< TRealPoint, TRealVector >::operator= ( PolygonalCalculus< TRealPoint, TRealVector > &&  other)
delete

Deleted move assignment operator.

Parameters
otherthe object to move.
Returns
a reference on 'this'.

◆ P()

◆ project()

template<typename TRealPoint , typename TRealVector >
static Vector DGtal::PolygonalCalculus< TRealPoint, TRealVector >::project ( const Vector u,
const Vector n 
)
inlinestaticprotected

Project u on the orthgonal of n

Parameters
uvector to project
nvector to build orthogonal space from
Returns
projected vector

Definition at line 1109 of file PolygonalCalculus.h.

1110  {
1111  return u - (u.dot(n) / n.squaredNorm()) * n;
1112  }

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Tf(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Tv().

◆ Qvf()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Qvf ( const Vertex v,
const Face f 
) const
inline

https://math.stackexchange.com/questions/180418/calculate-rotation-matrix-to-align-vector-a-to-vector-b-in-3d

Returns
3x3 Rotation matrix to align n_v to n_f

Definition at line 630 of file PolygonalCalculus.h.

631  {
632  Eigen::Vector3d nf = faceNormal(f);
633  Eigen::Vector3d nv = n_v(v);
634  double c = nv.dot(nf);
635  ASSERT(std::abs( c + 1.0) > 0.0001);
636  //Special case for opposite nv and nf vectors.
637  if (std::abs( c + 1.0) < 0.00001)
638  return -Eigen::Matrix3d::Identity();
639 
640  auto vv = nv.cross(nf);
641  DenseMatrix skew = bracket(vv);
642  return Eigen::Matrix3d::Identity() + skew +
643  1.0 / (1.0 + c) * skew * skew;
644  }
Eigen::Vector3d n_v(const Vertex &v) const

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::bracket(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceNormal(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::n_v().

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Rvf().

◆ Rvf()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Rvf ( const Vertex v,
const Face f 
) const
inline
Returns
Levi-Civita connection from vertex v tangent space to face f tangent space (2x2 rotation matrix)

Definition at line 648 of file PolygonalCalculus.h.

649  {
650  return Tf(f).transpose() * Qvf(v,f) * Tv(v);
651  }
DenseMatrix Tv(const Vertex &v) const
DenseMatrix Qvf(const Vertex &v, const Face &f) const

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Qvf(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Tf(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Tv().

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::blockConnection(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::transportAndFormatVectorField().

◆ selfDisplay()

template<typename TRealPoint , typename TRealVector >
void DGtal::PolygonalCalculus< TRealPoint, TRealVector >::selfDisplay ( std::ostream &  out) const
inline

Writes/Displays the object on an output stream.

Parameters
outthe output stream where the object is written.

Definition at line 1039 of file PolygonalCalculus.h.

1040  {
1041  out << "[PolygonalCalculus]: ";
1043  out<< "internal cache enabled, ";
1044  else
1045  out<<"internal cache disabled, ";
1046  out <<"SurfaceMesh="<<*mySurfaceMesh;
1047  }

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCacheEnabled, and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh.

◆ setEmbedder()

template<typename TRealPoint , typename TRealVector >
void DGtal::PolygonalCalculus< TRealPoint, TRealVector >::setEmbedder ( const std::function< Real3dPoint(Face, Vertex)> &  externalFunctor)
inline

Update the embedding function.

Parameters
externalFunctora new embedding functor (Face,Vertex)->RealPoint.

Definition at line 222 of file PolygonalCalculus.h.

223  {
224  myEmbedder = externalFunctor;
225  }

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myEmbedder.

◆ setInCache()

template<typename TRealPoint , typename TRealVector >
void DGtal::PolygonalCalculus< TRealPoint, TRealVector >::setInCache ( OPERATOR  key,
const Face  f,
const DenseMatrix ope 
) const
inlineprotected

◆ shapeOperator()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::shapeOperator ( const Face  f) const
inline

Shape operator on the face f (2x2 operator).

Returns
the shape operator at face f

Definition at line 655 of file PolygonalCalculus.h.

656  {
657  DenseMatrix N(myFaceDegree[f],3);
658  uint cpt = 0;
660  {
661  N.block(cpt,0,3,1) = n_v(v).transpose();
662  cpt++;
663  }
664  DenseMatrix GN = gradient(f) * N, Tf = T_f(f);
665 
666  return 0.5 * Tf.transpose() * (GN + GN.transpose()) * Tf;
667  }

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::gradient(), DGtal::SurfaceMesh< TRealPoint, TRealVector >::incidentVertices(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myFaceDegree, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::n_v(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Tf().

◆ sharp()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::sharp ( const Face  f) const
inline

◆ Tf()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Tf ( const Face f) const
inline
Returns
3x2 matrix defining the tangent space at face f, with basis vectors in columns

Definition at line 586 of file PolygonalCalculus.h.

587  {
588  Eigen::Vector3d nf = faceNormal(f);
589  ASSERT(std::abs(nf.norm() - 1.0) < 0.001);
590  const auto & N = getSurfaceMeshPtr()->incidentVertices(f);
591  auto v1 = *(N.begin());
592  auto v2 = *(N.begin() + 1);
593  Real3dPoint tangentVector =
595  Eigen::Vector3d w = toVec3(tangentVector);
596  Eigen::Vector3d uu = project(w,nf).normalized();
597  Eigen::Vector3d vv = nf.cross(uu);
598 
599  DenseMatrix tanB(3,2);
600  tanB.col(0) = uu;
601  tanB.col(1) = vv;
602  return tanB;
603  }
MySurfaceMesh::RealPoint Real3dPoint
Position type.
static Vector project(const Vector &u, const Vector &n)

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceNormal(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::getSurfaceMeshPtr(), DGtal::SurfaceMesh< TRealPoint, TRealVector >::incidentVertices(), DGtal::SurfaceMesh< TRealPoint, TRealVector >::position(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::project(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::toVec3().

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantGradient(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantGradient_f(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Rvf(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::shapeOperator().

◆ toExtrinsicVector()

template<typename TRealPoint , typename TRealVector >
Vector DGtal::PolygonalCalculus< TRealPoint, TRealVector >::toExtrinsicVector ( const Vertex  v,
const Vector I 
) const
inline

toExtrinsicVector

Parameters
vthe vertex
Ithe intrinsic vector at Tv
Returns
3D extrinsic vector from intrinsic 2D vector I expressed from tangent frame at vertex v

Definition at line 610 of file PolygonalCalculus.h.

611  {
612  DenseMatrix T = Tv(v);
613  return T.col(0) * I(0) + T.col(1) * I(1);
614  }

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Tv().

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::toExtrinsicVectors().

◆ toExtrinsicVectors()

template<typename TRealPoint , typename TRealVector >
std::vector<Vector> DGtal::PolygonalCalculus< TRealPoint, TRealVector >::toExtrinsicVectors ( const std::vector< Vector > &  I) const
inline
Parameters
Iset of intrinsic vectors, vectors indices must be the same as their associated vertex
Returns
converts a set of intrinsic vectors to their extrinsic equivalent, expressed in correponding tangent frame

Definition at line 620 of file PolygonalCalculus.h.

621  {
622  std::vector<Vector> ext(mySurfaceMesh->nbVertices());
623  for (auto v = 0; v < mySurfaceMesh->nbVertices(); v++)
624  ext[v] = toExtrinsicVector(v,I[v]);
625  return ext;
626  }
Vector toExtrinsicVector(const Vertex v, const Vector &I) const
toExtrinsicVector

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh, DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbVertices(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::toExtrinsicVector().

◆ toReal3dVector()

template<typename TRealPoint , typename TRealVector >
static Real3dVector DGtal::PolygonalCalculus< TRealPoint, TRealVector >::toReal3dVector ( const Eigen::Vector3d &  x)
inlinestaticprotected

toReal3dVector converts Eigen::Vector3d to Real3dVector.

Conversion routines.

Parameters
xthe vector
Returns
the same vector in DGtal type

Definition at line 1138 of file PolygonalCalculus.h.

1139  {
1140  return { x(0), x(1), x(2)};
1141  }

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus().

◆ toVec3()

template<typename TRealPoint , typename TRealVector >
static Eigen::Vector3d DGtal::PolygonalCalculus< TRealPoint, TRealVector >::toVec3 ( const Real3dPoint x)
inlinestaticprotected

toVec3 convert Real3dPoint to Eigen::Vector3d

Parameters
xthe vector
Returns
the same vector in eigen type

Definition at line 1129 of file PolygonalCalculus.h.

1130  {
1131  return Eigen::Vector3d(x(0),x(1),x(2));
1132  }

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::n_v(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Tf(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Tv().

◆ toVector()

template<typename TRealPoint , typename TRealVector >
static Vector DGtal::PolygonalCalculus< TRealPoint, TRealVector >::toVector ( const Eigen::Vector3d &  x)
inlinestaticprotected

toVector convert Real3dPoint to Eigen::VectorXd

Conversion routines.

Parameters
xthe vector
Returns
the same vector in eigen type

Definition at line 1118 of file PolygonalCalculus.h.

1119  {
1120  Vector X(3);
1121  for (int i = 0; i < 3; i++)
1122  X(i) = x(i);
1123  return X;
1124  }

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::X().

◆ transportAndFormatVectorField()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::transportAndFormatVectorField ( const Face  f,
const Vector uf 
)
inline
Returns
to fit [44] paper's notations, this function maps all the per vertex vectors (expressed in the (2*nf) vector form) into the nfx2 matrix with transported vectors (to face f) in each row.
Note
Unlike the rest of the per face operators, the covariant operators need to be applied directly to the restriction of the vector field to the face,

Definition at line 676 of file PolygonalCalculus.h.

677  {
678  DenseMatrix uf_nabla(myFaceDegree[f], 2);
679  size_t cpt = 0;
680  for (auto v : mySurfaceMesh->incidentVertices(f))
681  {
682  uf_nabla.block(cpt,0,1,2) =
683  (Rvf(v,f) * uf.block(2 * cpt,0,2,1)).transpose();
684  ++cpt;
685  }
686  return uf_nabla;
687  }

References DGtal::SurfaceMesh< TRealPoint, TRealVector >::incidentVertices(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myFaceDegree, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh, and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Rvf().

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantGradient(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantProjection().

◆ Tv()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Tv ( const Vertex v) const
inline
Returns
3x2 matrix defining the tangent space at vertex v, with basis vectors in columns

Definition at line 566 of file PolygonalCalculus.h.

567  {
568  Eigen::Vector3d nv = n_v(v);
569  ASSERT(std::abs(nv.norm() - 1.0) < 0.001);
570  const auto & N = getSurfaceMeshPtr()->neighborVertices(v);
571  auto neighbor = *N.begin();
572  Real3dPoint tangentVector = getSurfaceMeshPtr()->position(v) -
573  getSurfaceMeshPtr()->position(neighbor);
574  Eigen::Vector3d w = toVec3(tangentVector);
575  Eigen::Vector3d uu = project(w,nv).normalized();
576  Eigen::Vector3d vv = nv.cross(uu);
577 
578  DenseMatrix tanB(3,2);
579  tanB.col(0) = uu;
580  tanB.col(1) = vv;
581  return tanB;
582  }
const Vertices & neighborVertices(Vertex v) const
Definition: SurfaceMesh.h:331

References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::getSurfaceMeshPtr(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::n_v(), DGtal::SurfaceMesh< TRealPoint, TRealVector >::neighborVertices(), DGtal::SurfaceMesh< TRealPoint, TRealVector >::position(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::project(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::toVec3().

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Rvf(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::toExtrinsicVector().

◆ updateFaceDegree()

template<typename TRealPoint , typename TRealVector >
void DGtal::PolygonalCalculus< TRealPoint, TRealVector >::updateFaceDegree ( )
inlineprotected

◆ vectorArea()

template<typename TRealPoint , typename TRealVector >
Vector DGtal::PolygonalCalculus< TRealPoint, TRealVector >::vectorArea ( const Face  f) const
inline

Polygonal (corrected) vector area.

Parameters
fthe face
Returns
a vector oriented in the (corrected) normal direction and with length equal to the (corrected) area of the face f.

Definition at line 316 of file PolygonalCalculus.h.

317  {
318  Real3dPoint af(0.0,0.0,0.0);
319  const auto vertices = mySurfaceMesh->incidentVertices(f);
320  auto it = vertices.cbegin();
321  auto itnext = vertices.cbegin();
322  ++itnext;
323  while (it != vertices.cend())
324  {
325  auto xi = myEmbedder(f,*it);
326  auto xip = myEmbedder(f,*itnext);
327  af += xi.crossProduct(xip);
328  ++it;
329  ++itnext;
330  if (itnext == vertices.cend())
331  itnext = vertices.cbegin();
332  }
333  Eigen::Vector3d output = {af[0],af[1],af[2]};
334  return 0.5*output;
335  }

References DGtal::SurfaceMesh< TRealPoint, TRealVector >::incidentVertices(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myEmbedder, and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh.

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::computeVertexNormal(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceArea(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceNormal(), and TEST_CASE().

◆ X()

template<typename TRealPoint , typename TRealVector >
DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::X ( const Face  f) const
inline

Field Documentation

◆ dimension

template<typename TRealPoint , typename TRealVector >
const Dimension DGtal::PolygonalCalculus< TRealPoint, TRealVector >::dimension = TRealPoint::dimension
static

Concept checking.

Definition at line 76 of file PolygonalCalculus.h.

◆ myEmbedder

template<typename TRealPoint , typename TRealVector >
std::function<Real3dPoint(Face, Vertex)> DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myEmbedder
private

◆ myFaceDegree

◆ myGlobalCache

◆ myGlobalCacheEnabled

◆ mySurfaceMesh

template<typename TRealPoint , typename TRealVector >
const MySurfaceMesh* DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh
private

◆ myVertexNormalEmbedder

template<typename TRealPoint , typename TRealVector >
std::function<Real3dVector(Vertex)> DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myVertexNormalEmbedder
private

Embedding function (vertex)->R^3 for the vertex normal.

Definition at line 1226 of file PolygonalCalculus.h.

Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::n_v(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus().


The documentation for this class was generated from the following file: