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

Implements differential operators on polygonal surfaces from [45]. 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 [45].

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 [45]. 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 128 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 158 of file PolygonalCalculus.h.

◆ Face

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

Face type.

Definition at line 145 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 156 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 152 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 141 of file PolygonalCalculus.h.

◆ Real3dPoint

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

Position type.

Definition at line 147 of file PolygonalCalculus.h.

◆ Real3dVector

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

Real vector type.

Definition at line 149 of file PolygonalCalculus.h.

◆ Self

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

Self type.

Definition at line 138 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 165 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 160 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 162 of file PolygonalCalculus.h.

◆ Vector

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

Type of Vector.

Definition at line 154 of file PolygonalCalculus.h.

◆ Vertex

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

Vertex type.

Definition at line 143 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 174 of file PolygonalCalculus.h.

175  :
176  mySurfaceMesh(&surf), myGlobalCacheEnabled(globalInternalCacheEnabled)
177  {
178  myEmbedder =[&](Face f,Vertex v){ (void)f; return mySurfaceMesh->position(v);};
180  init();
181  };
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 189 of file PolygonalCalculus.h.

191  :
192  mySurfaceMesh(&surf), myEmbedder(embedder), myGlobalCacheEnabled(globalInternalCacheEnabled)
193  {
195  init();
196  };

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 204 of file PolygonalCalculus.h.

206  :
207  mySurfaceMesh(&surf), myVertexNormalEmbedder(embedder),
208  myGlobalCacheEnabled(globalInternalCacheEnabled)
209  {
210  myEmbedder = [&](Face f,Vertex v){(void)f; return mySurfaceMesh->position(v); };
211  init();
212  };

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 224 of file PolygonalCalculus.h.

227  :
228  mySurfaceMesh(&surf), myEmbedder(pos_embedder),
229  myVertexNormalEmbedder(normal_embedder),
230  myGlobalCacheEnabled(globalInternalCacheEnabled)
231  {
232  init();
233  };

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 353 of file PolygonalCalculus.h.

354  {
355  if (checkCache(A_,f))
356  return myGlobalCache[A_][f];
357 
358  const auto nf = myFaceDegree[f];
359  DenseMatrix a = DenseMatrix::Zero(nf ,nf);
360  for(auto i=0u; i < nf; ++i)
361  {
362  a(i, (i+1)%nf) = 0.5;
363  a(i,i) = 0.5;
364  }
365 
366  setInCache(A_,f,a);
367  return a;
368  }
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(), initQuantities(), 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 1240 of file PolygonalCalculus.h.

1241  {
1242  auto nf = degree(f);
1243  DenseMatrix RU_fO = DenseMatrix::Zero(nf * 2,nf * 2);
1244  size_t cpt = 0;
1245  for (auto v : getSurfaceMeshPtr()->incidentVertices(f))
1246  {
1247  auto Rv = Rvf(v,f);
1248  RU_fO.block<2,2>(2 * cpt,2 * cpt) = Rv;
1249  ++cpt;
1250  }
1251  return RU_fO;
1252  }
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 436 of file PolygonalCalculus.h.

437  {
438  DenseMatrix brack(3,3);
439  brack << 0.0 , -n(2), n(1),
440  n(2), 0.0 , -n(0),
441  -n(1) , n(0),0.0 ;
442  return brack;
443  }

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 486 of file PolygonalCalculus.h.

487  {
488  const auto nf = myFaceDegree[f];
489  return 1.0/(double)nf * X(f).transpose() * Vector::Ones(nf);
490  }

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 494 of file PolygonalCalculus.h.

495  {
496  const Vector c = centroid(f);
497  return {c(0),c(1),c(2)};
498  }
Vector centroid(const Face f) const
DigitalPlane::Point Vector

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

Referenced by initQuantities(), and 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 1207 of file PolygonalCalculus.h.

1208  {
1209  Vector n(3);
1210  n(0) = 0.;
1211  n(1) = 0.;
1212  n(2) = 0.;
1213  /* for (auto f : mySurfaceMesh->incidentFaces(v))
1214  n += vectorArea(f);
1215  return n.normalized();
1216  */
1217  auto faces = mySurfaceMesh->incidentFaces(v);
1218  for (auto f : faces)
1219  n += vectorArea(f);
1220 
1221  if (fabs(n.norm() - 0.0) < 0.00001)
1222  {
1223  //On non-manifold edges touching the boundary, n may be null.
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);
1226  }
1227  n = n.normalized();
1228  return n;
1229  }
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 [45]. 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 810 of file PolygonalCalculus.h.

811  {
812  if (checkCache(CON_L_,f))
813  return myGlobalCache[CON_L_][f];
816  DenseMatrix L = -(faceArea(f) * G.transpose() * G + lambda * P.transpose() * P);
817  setInCache(CON_L_,f,L);
818  return L;
819  }
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 756 of file PolygonalCalculus.h.

757  {
758  return Tf(f).transpose() * gradient(f) *
760  }
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 782 of file PolygonalCalculus.h.

783  {
784  return kroneckerWithI2(Tf(f).transpose() * gradient(f)) * blockConnection(f);
785  }
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 772 of file PolygonalCalculus.h.

773  {
774  return P(f) * D(f) * transportAndFormatVectorField(f,uf);
775  }
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 792 of file PolygonalCalculus.h.

793  {
794  return kroneckerWithI2(P(f) * D(f)) * blockConnection(f);
795  ;
796  }

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 576 of file PolygonalCalculus.h.

577  {
578  if (checkCache(CURL_,f))
579  return myGlobalCache[CURL_][f];
580 
581  DenseMatrix op = DenseMatrix::Identity(myFaceDegree[f],myFaceDegree[f]);
582 
583  setInCache(CURL_,f,op);
584  return op;
585  }

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 1082 of file PolygonalCalculus.h.

1083  {
1084  return myFaceDegree[f];
1085  }

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 1040 of file PolygonalCalculus.h.

1041  {
1042  myGlobalCacheEnabled = false;
1043  myGlobalCache.clear();
1044  }

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 [45] . 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 562 of file PolygonalCalculus.h.

563  {
564  if (checkCache(DIVERGENCE_,f))
565  return myGlobalCache[DIVERGENCE_][f];
566 
567  DenseMatrix op = -1.0 * D(f).transpose() * M(f);
568  setInCache(DIVERGENCE_,f,op);
569 
570  return op;
571  }
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 960 of file PolygonalCalculus.h.

961  {
962  auto nv = mySurfaceMesh->nbVertices();
963  SparseMatrix M(2 * nv, 2 * nv);
964  std::vector<Triplet> triplets;
965  for (typename MySurfaceMesh::Index v = 0; v < mySurfaceMesh->nbVertices(); ++v)
966  {
967  auto faces = mySurfaceMesh->incidentFaces(v);
968  auto varea = 0.0;
969  for (auto f : faces)
970  varea += faceArea(f) / (double)myFaceDegree[f];
971  triplets.emplace_back(Triplet(2 * v, 2 * v, varea));
972  triplets.emplace_back(Triplet(2 * v + 1, 2 * v + 1, varea));
973  }
974  M.setFromTriplets(triplets.begin(), triplets.end());
975  return M;
976  }
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 1033 of file PolygonalCalculus.h.

1034  {
1035  myGlobalCacheEnabled = true;
1036  }

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 1063 of file PolygonalCalculus.h.

1064  {
1065  return myFaceDegree[f];
1066  }

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 416 of file PolygonalCalculus.h.

417  {
418  Vector v = faceNormal(f);
419  return {v(0),v(1),v(2)};
420  }
Vector faceNormal(const Face f) const

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

Referenced by initQuantities(), initQuantitiesCached(), and 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 828 of file PolygonalCalculus.h.

829  {
830  return Form::Zero( nbVertices() );
831  }

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

Referenced by computeLaplace(), and 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.
PolygonalCalculus< SH3::RealPoint, SH3::RealVector >::Vector phi(const Face f)
PolyCalculus * calculus
Parameters
perFaceOperatorthe per face operator
Returns
an indexed container of all DenseMatrix operators (indexed per Face).

Definition at line 999 of file PolygonalCalculus.h.

1000  {
1001  std::vector<DenseMatrix> cache;
1002  for( typename MySurfaceMesh::Index f=0; f < mySurfaceMesh->nbFaces(); ++f)
1003  cache.push_back(perFaceOperator(f));
1004  return cache;
1005  }
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 1023 of file PolygonalCalculus.h.

1024  {
1025  std::vector<Vector> cache;
1026  for( typename MySurfaceMesh::Index f=0; f < mySurfaceMesh->nbFaces(); ++f)
1027  cache.push_back(perFaceVectorOperator(f));
1028  return cache;
1029  }

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 [45]. 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 928 of file PolygonalCalculus.h.

929  {
930  auto nv = mySurfaceMesh->nbVertices();
931  SparseMatrix lapGlobal(2 * nv, 2 * nv);
932  SparseMatrix local(2 * nv, 2 * nv);
933  std::vector<Triplet> triplets;
934  for (typename MySurfaceMesh::Index f = 0; f < mySurfaceMesh->nbFaces(); f++)
935  {
936  auto nf = degree(f);
937  DenseMatrix Lap = connectionLaplacian(f,lambda);
938  const auto vertices = mySurfaceMesh->incidentVertices(f);
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++)
943  {
944  auto v = Lap(2 * i + k1, 2 * j + k2);
945  if (v != 0.0)
946  triplets.emplace_back(Triplet(2 * vertices[i] + k1,
947  2 * vertices[j] + k2, v));
948  }
949  }
950  lapGlobal.setFromTriplets(triplets.begin(), triplets.end());
951  return lapGlobal;
952  }
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 902 of file PolygonalCalculus.h.

903  {
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();
908  return iM0;
909  }
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 [45]. 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 855 of file PolygonalCalculus.h.

856  {
858  std::vector<Triplet> triplets;
859  for( typename MySurfaceMesh::Index f = 0; f < mySurfaceMesh->nbFaces(); ++f )
860  {
861  auto nf = myFaceDegree[f];
862  DenseMatrix Lap = this->laplaceBeltrami(f,lambda);
863  const auto vertices = mySurfaceMesh->incidentVertices(f);
864  for(auto i=0u; i < nf; ++i)
865  for(auto j=0u; j < nf; ++j)
866  {
867  auto v = Lap(i,j);
868  if (v!= 0.0)
869  triplets.emplace_back( Triplet( (SparseMatrix::StorageIndex)vertices[ i ], (SparseMatrix::StorageIndex)vertices[ j ],
870  Lap( i, j ) ) );
871  }
872  }
873  lapGlobal.setFromTriplets(triplets.begin(), triplets.end());
874  return lapGlobal;
875  }
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 computeLaplace(), and 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 882 of file PolygonalCalculus.h.

883  {
885  std::vector<Triplet> triplets;
886  for ( typename MySurfaceMesh::Index v = 0; v < mySurfaceMesh->nbVertices(); ++v )
887  {
888  auto faces = mySurfaceMesh->incidentFaces(v);
889  auto varea = 0.0;
890  for(auto f: faces)
891  varea += faceArea(f) /(double)myFaceDegree[f];
892  triplets.emplace_back(Triplet(v,v,varea));
893  }
894  M.setFromTriplets(triplets.begin(),triplets.end());
895  return M;
896  }

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 833 of file PolygonalCalculus.h.

834  {
835  SparseMatrix Id0( nbVertices(), nbVertices() );
836  Id0.setIdentity();
837  return Id0;
838  }

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 1055 of file PolygonalCalculus.h.

1056  {
1057  updateFaceDegree();
1058  }
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 1111 of file PolygonalCalculus.h.

1112  {
1113  return true;
1114  }

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 1255 of file PolygonalCalculus.h.

1256  {
1257  size_t h = M.rows();
1258  size_t w = M.cols();
1259  DenseMatrix MK = DenseMatrix::Zero(h * 2,w * 2);
1260  for (size_t j = 0; j < h; j++)
1261  for (size_t i = 0; i < w; i++)
1262  {
1263  MK(2 * j, 2 * i) = M(j, i);
1264  MK(2 * j + 1, 2 * i + 1) = M(j, i);
1265  }
1266  return MK;
1267  }

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 [45] . 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 602 of file PolygonalCalculus.h.

603  {
604  if (checkCache(L_,f))
605  return myGlobalCache[L_][f];
606 
607  DenseMatrix Df = D(f);
608  // Laplacian is a negative operator.
609  DenseMatrix op = -1.0 * Df.transpose() * M(f,lambda) * Df;
610 
611  setInCache(L_, f, op);
612  return op;
613  }

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 1233 of file PolygonalCalculus.h.

1234  {
1235  return toVec3(myVertexNormalEmbedder(v));
1236  }
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 1075 of file PolygonalCalculus.h.

1076  {
1077  return mySurfaceMesh->nbFaces();
1078  }

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 1167 of file PolygonalCalculus.h.

1168  {
1169  return u - (u.dot(n) / n.squaredNorm()) * n;
1170  }

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 688 of file PolygonalCalculus.h.

689  {
690  Eigen::Vector3d nf = faceNormal(f);
691  Eigen::Vector3d nv = n_v(v);
692  double c = nv.dot(nf);
693  ASSERT(std::abs( c + 1.0) > 0.0001);
694  //Special case for opposite nv and nf vectors.
695  if (std::abs( c + 1.0) < 0.00001)
696  return -Eigen::Matrix3d::Identity();
697 
698  auto vv = nv.cross(nf);
699  DenseMatrix skew = bracket(vv);
700  return Eigen::Matrix3d::Identity() + skew +
701  1.0 / (1.0 + c) * skew * skew;
702  }
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 706 of file PolygonalCalculus.h.

707  {
708  return Tf(f).transpose() * Qvf(v,f) * Tv(v);
709  }
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 1097 of file PolygonalCalculus.h.

1098  {
1099  out << "[PolygonalCalculus]: ";
1101  out<< "internal cache enabled, ";
1102  else
1103  out<<"internal cache disabled, ";
1104  out <<"SurfaceMesh="<<*mySurfaceMesh;
1105  }

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 280 of file PolygonalCalculus.h.

281  {
282  myEmbedder = externalFunctor;
283  }

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

Referenced by precompute().

◆ 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 713 of file PolygonalCalculus.h.

714  {
715  DenseMatrix N(myFaceDegree[f],3);
716  uint cpt = 0;
718  {
719  N.block(cpt,0,3,1) = n_v(v).transpose();
720  cpt++;
721  }
722  DenseMatrix GN = gradient(f) * N, Tf = T_f(f);
723 
724  return 0.5 * Tf.transpose() * (GN + GN.transpose()) * Tf;
725  }

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 644 of file PolygonalCalculus.h.

645  {
646  Eigen::Vector3d nf = faceNormal(f);
647  ASSERT(std::abs(nf.norm() - 1.0) < 0.001);
648  const auto & N = getSurfaceMeshPtr()->incidentVertices(f);
649  auto v1 = *(N.begin());
650  auto v2 = *(N.begin() + 1);
651  Real3dPoint tangentVector =
653  Eigen::Vector3d w = toVec3(tangentVector);
654  Eigen::Vector3d uu = project(w,nf).normalized();
655  Eigen::Vector3d vv = nf.cross(uu);
656 
657  DenseMatrix tanB(3,2);
658  tanB.col(0) = uu;
659  tanB.col(1) = vv;
660  return tanB;
661  }
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 668 of file PolygonalCalculus.h.

669  {
670  DenseMatrix T = Tv(v);
671  return T.col(0) * I(0) + T.col(1) * I(1);
672  }

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 678 of file PolygonalCalculus.h.

679  {
680  std::vector<Vector> ext(mySurfaceMesh->nbVertices());
681  for (auto v = 0; v < mySurfaceMesh->nbVertices(); v++)
682  ext[v] = toExtrinsicVector(v,I[v]);
683  return ext;
684  }
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 1196 of file PolygonalCalculus.h.

1197  {
1198  return { x(0), x(1), x(2)};
1199  }

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 1187 of file PolygonalCalculus.h.

1188  {
1189  return Eigen::Vector3d(x(0),x(1),x(2));
1190  }

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 1176 of file PolygonalCalculus.h.

1177  {
1178  Vector X(3);
1179  for (int i = 0; i < 3; i++)
1180  X(i) = x(i);
1181  return X;
1182  }

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 [45] 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 734 of file PolygonalCalculus.h.

735  {
736  DenseMatrix uf_nabla(myFaceDegree[f], 2);
737  size_t cpt = 0;
738  for (auto v : mySurfaceMesh->incidentVertices(f))
739  {
740  uf_nabla.block(cpt,0,1,2) =
741  (Rvf(v,f) * uf.block(2 * cpt,0,2,1)).transpose();
742  ++cpt;
743  }
744  return uf_nabla;
745  }

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 624 of file PolygonalCalculus.h.

625  {
626  Eigen::Vector3d nv = n_v(v);
627  ASSERT(std::abs(nv.norm() - 1.0) < 0.001);
628  const auto & N = getSurfaceMeshPtr()->neighborVertices(v);
629  auto neighbor = *N.begin();
630  Real3dPoint tangentVector = getSurfaceMeshPtr()->position(v) -
631  getSurfaceMeshPtr()->position(neighbor);
632  Eigen::Vector3d w = toVec3(tangentVector);
633  Eigen::Vector3d uu = project(w,nv).normalized();
634  Eigen::Vector3d vv = nv.cross(uu);
635 
636  DenseMatrix tanB(3,2);
637  tanB.col(0) = uu;
638  tanB.col(1) = vv;
639  return tanB;
640  }
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 374 of file PolygonalCalculus.h.

375  {
376  Real3dPoint af(0.0,0.0,0.0);
377  const auto vertices = mySurfaceMesh->incidentVertices(f);
378  auto it = vertices.cbegin();
379  auto itnext = vertices.cbegin();
380  ++itnext;
381  while (it != vertices.cend())
382  {
383  auto xi = myEmbedder(f,*it);
384  auto xip = myEmbedder(f,*itnext);
385  af += xi.crossProduct(xip);
386  ++it;
387  ++itnext;
388  if (itnext == vertices.cend())
389  itnext = vertices.cbegin();
390  }
391  Eigen::Vector3d output = {af[0],af[1],af[2]};
392  return 0.5*output;
393  }

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(), initQuantities(), initQuantitiesCached(), 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 134 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 1284 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: