DGtal 1.4.0
|
Implements differential operators on polygonal surfaces from [45]. More...
#include <DGtal/dec/PolygonalCalculus.h>
Public Types | |
typedef PolygonalCalculus< TRealPoint, TRealVector > | Self |
Self type. | |
typedef SurfaceMesh< TRealPoint, TRealVector > | MySurfaceMesh |
Type of SurfaceMesh. | |
typedef MySurfaceMesh::Vertex | Vertex |
Vertex type. | |
typedef MySurfaceMesh::Face | Face |
Face type. | |
typedef MySurfaceMesh::RealPoint | Real3dPoint |
Position type. | |
typedef MySurfaceMesh::RealVector | Real3dVector |
Real vector type. | |
typedef EigenLinearAlgebraBackend | LinAlg |
Linear Algebra Backend from Eigen. | |
typedef LinAlg::DenseVector | Vector |
Type of Vector. | |
typedef Vector | Form |
Global 0-form, 1-form, 2-form are Vector. | |
typedef LinAlg::DenseMatrix | DenseMatrix |
Type of dense matrix. | |
typedef LinAlg::SparseMatrix | SparseMatrix |
Type of sparse matrix. | |
typedef LinAlg::Triplet | Triplet |
Type of sparse matrix triplet. | |
typedef LinAlg::SolverSimplicialLDLT | Solver |
Type of a sparse matrix solver. | |
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 | |
PolygonalCalculus & | operator= (const PolygonalCalculus &other)=delete |
PolygonalCalculus & | operator= (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 | |
std::vector< Vector > | toExtrinsicVectors (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< DenseMatrix > | getOperatorCacheMatrix (const std::function< DenseMatrix(Face)> &perFaceOperator) const |
std::vector< Vector > | getOperatorCacheVector (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 MySurfaceMesh * | getSurfaceMeshPtr () const |
void | selfDisplay (std::ostream &out) const |
bool | isValid () const |
Static Public Attributes | |
static const Dimension | dimension = TRealPoint::dimension |
Concept checking. | |
Private Attributes | |
const MySurfaceMesh * | mySurfaceMesh |
Underlying SurfaceMesh. | |
std::function< Real3dPoint(Face, Vertex)> | myEmbedder |
Embedding function (face,vertex)->R^3 for the vertex position wrt. the face. | |
std::function< Real3dVector(Vertex)> | myVertexNormalEmbedder |
Embedding function (vertex)->R^3 for the vertex normal. | |
std::vector< size_t > | myFaceDegree |
Cache containing the face degree. | |
bool | myGlobalCacheEnabled |
Global cache. | |
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. | |
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 | |
static Eigen::Vector3d | toVec3 (const Real3dPoint &x) |
toVec3 convert Real3dPoint to Eigen::Vector3d | |
static Real3dVector | toReal3dVector (const Eigen::Vector3d &x) |
toReal3dVector converts Eigen::Vector3d to Real3dVector. | |
Implements differential operators on polygonal surfaces from [45].
Description of template class 'PolygonalCalculus'
See Discrete differential calculus on polygonal surfaces for details.
TRealPoint | a model of points \(\mathbb{R}^3\) (e.g. PointVector). |
TRealVector | a model of vectors in \(\mathbb{R}^3\) (e.g. PointVector). |
Definition at line 128 of file PolygonalCalculus.h.
typedef LinAlg::DenseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::DenseMatrix |
Type of dense matrix.
Definition at line 158 of file PolygonalCalculus.h.
typedef MySurfaceMesh::Face DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Face |
Face type.
Definition at line 145 of file PolygonalCalculus.h.
typedef Vector DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Form |
Global 0-form, 1-form, 2-form are Vector.
Definition at line 156 of file PolygonalCalculus.h.
typedef EigenLinearAlgebraBackend DGtal::PolygonalCalculus< TRealPoint, TRealVector >::LinAlg |
Linear Algebra Backend from Eigen.
Definition at line 152 of file PolygonalCalculus.h.
typedef SurfaceMesh<TRealPoint, TRealVector> DGtal::PolygonalCalculus< TRealPoint, TRealVector >::MySurfaceMesh |
Type of SurfaceMesh.
Definition at line 141 of file PolygonalCalculus.h.
typedef MySurfaceMesh::RealPoint DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Real3dPoint |
Position type.
Definition at line 147 of file PolygonalCalculus.h.
typedef MySurfaceMesh::RealVector DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Real3dVector |
Real vector type.
Definition at line 149 of file PolygonalCalculus.h.
typedef PolygonalCalculus<TRealPoint, TRealVector> DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Self |
Self type.
Definition at line 138 of file PolygonalCalculus.h.
typedef LinAlg::SolverSimplicialLDLT DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Solver |
Type of a sparse matrix solver.
Definition at line 165 of file PolygonalCalculus.h.
typedef LinAlg::SparseMatrix DGtal::PolygonalCalculus< TRealPoint, TRealVector >::SparseMatrix |
Type of sparse matrix.
Definition at line 160 of file PolygonalCalculus.h.
typedef LinAlg::Triplet DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Triplet |
Type of sparse matrix triplet.
Definition at line 162 of file PolygonalCalculus.h.
typedef LinAlg::DenseVector DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Vector |
Type of Vector.
Definition at line 154 of file PolygonalCalculus.h.
typedef MySurfaceMesh::Vertex DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Vertex |
Vertex type.
Definition at line 143 of file PolygonalCalculus.h.
|
protected |
Enum for operators in the internal cache strategy.
Enumerator | |
---|---|
X_ | |
D_ | |
E_ | |
A_ | |
COGRAD_ | |
GRAD_ | |
FLAT_ | |
B_ | |
SHARP_ | |
P_ | |
M_ | |
DIVERGENCE_ | |
CURL_ | |
L_ | |
CON_L_ |
Definition at line 1126 of file PolygonalCalculus.h.
|
inline |
Create a Polygonal DEC structure from a surface mesh (surf) using an default identity embedder.
surf | an instance of SurfaceMesh |
globalInternalCacheEnabled | enable the internal cache for all operators (default: false) |
Definition at line 174 of file PolygonalCalculus.h.
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().
|
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.
surf | an instance of SurfaceMesh |
embedder | an embedder for the vertex position |
globalInternalCacheEnabled | if true, the global operator cache is enabled |
Definition at line 189 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::computeVertexNormal(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::init(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myVertexNormalEmbedder, and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::toReal3dVector().
|
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.
surf | an instance of SurfaceMesh |
embedder | an embedder for the vertex position |
globalInternalCacheEnabled | if true, the global operator cache is enabled |
Definition at line 204 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::init(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myEmbedder, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh, and DGtal::SurfaceMesh< TRealPoint, TRealVector >::position().
|
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.
surf | an instance of SurfaceMesh |
pos_embedder | an embedder for the position |
normal_embedder | an embedder for the position |
globalInternalCacheEnabled |
Definition at line 224 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::init().
|
delete |
Deleted default constructor.
|
default |
Destructor (default).
|
delete |
Deleted copy constructor.
other | the object to clone. |
|
delete |
Deleted move constructor.
other | the object to move. |
|
inline |
Average operator to average, per edge, its vertex values.
f | the face |
Definition at line 353 of file PolygonalCalculus.h.
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().
|
inline |
Edge mid-point operator of the face.
f | the face |
Definition at line 475 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::A(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::B_, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::checkCache(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCache, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::setInCache(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::X().
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::sharp().
|
inlineprotected |
Definition at line 1240 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::degree(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::getSurfaceMeshPtr(), DGtal::SurfaceMesh< TRealPoint, TRealVector >::incidentVertices(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Rvf().
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantGradient_f(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantProjection_f().
DGtal::PolygonalCalculus< TRealPoint, TRealVector >::BOOST_STATIC_ASSERT | ( | (dimension==3) | ) |
|
inline |
Return [n] as the 3x3 operator such that [n]q = n x q
n | a vector |
Definition at line 436 of file PolygonalCalculus.h.
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::gradient(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Qvf(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::sharp().
|
inline |
f | the face |
Definition at line 486 of file PolygonalCalculus.h.
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().
|
inline |
f | the face |
Definition at line 494 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::centroid().
Referenced by initQuantities(), and TEST_CASE().
|
inlineprotected |
Check internal cache if enabled.
key | the operator name |
f | the face |
Definition at line 1144 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCache, and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCacheEnabled.
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::A(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::B(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::coGradient(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::connectionLaplacian(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::curl(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::D(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::divergence(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::E(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::flat(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::gradient(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::laplaceBeltrami(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::M(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::P(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::sharp(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::X().
|
inline |
co-Gradient operator of the face
f | the face |
Definition at line 425 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::A(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::checkCache(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::COGRAD_, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::E(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCache, and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::setInCache().
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::gradient(), initQuantities(), initQuantitiesCached(), and TEST_CASE().
|
inlineprotected |
Compute the (normalized) normal vector at a Vertex by averaging the adjacent face normal vectors.
v | the vertex to compute the normal from |
Definition at line 1207 of file PolygonalCalculus.h.
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(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus().
|
inline |
L∇ := -(afG∇tG∇+λP∇tP∇)
Definition at line 810 of file PolygonalCalculus.h.
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().
|
inline |
Covarient gradient at a face a f of intrinsic vectors uf.
uf | list of all intrinsic vectors per vertex concatenated in a column vector |
f | the face |
Definition at line 756 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::gradient(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Tf(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::transportAndFormatVectorField().
Referenced by TEST_CASE().
|
inline |
Definition at line 782 of file PolygonalCalculus.h.
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().
|
inline |
Compute the covariance projection at a face f of intrinsic vectors uf.
uf | list of all intrinsic vectors per vertex concatenated in a column vector |
f | the face |
Definition at line 772 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::D(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::P(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::transportAndFormatVectorField().
Referenced by TEST_CASE().
|
inline |
Definition at line 792 of file PolygonalCalculus.h.
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().
|
inline |
Curl operator of a one-form (identity matrix).
f | the face |
Definition at line 576 of file PolygonalCalculus.h.
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().
|
inline |
Derivative operator (d_0) of a face.
f | the face |
Definition at line 319 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::checkCache(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::D_, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myFaceDegree, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCache, and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::setInCache().
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantProjection(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantProjection_f(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::divergence(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::E(), initQuantities(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::laplaceBeltrami(), and TEST_CASE().
|
inline |
f | the face |
Definition at line 1082 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myFaceDegree.
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::blockConnection(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::globalConnectionLaplace().
|
inline |
Disable the internal global cache for operators. This method will also clean up the
Definition at line 1040 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCache, and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCacheEnabled.
|
inline |
Divergence operator of a one-form.
f | the face |
Definition at line 562 of file PolygonalCalculus.h.
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().
|
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)
Definition at line 960 of file PolygonalCalculus.h.
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().
|
inline |
Edge vector operator per face.
f | the face |
Definition at line 339 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::checkCache(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::D(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::E_, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCache, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::setInCache(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::X().
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::coGradient(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::flat().
|
inline |
Enable the internal global cache for operators.
Definition at line 1033 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCacheEnabled.
|
inline |
Area of a face from the vector area.
f | the face |
Definition at line 398 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::vectorArea().
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::connectionLaplacian(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::doubledGlobalLumpedMassMatrix(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::globalLumpedMassMatrix(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::gradient(), initQuantities(), initQuantitiesCached(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::M(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::sharp(), and TEST_CASE().
|
inline |
Helper to retrieve the degree of the face from the cache.
f | the face |
Definition at line 1063 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myFaceDegree.
Referenced by TEST_CASE().
|
inline |
Corrected normal vector of a face.
f | the face |
Definition at line 406 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::vectorArea().
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceNormalAsDGtalVector(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::flat(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::gradient(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Qvf(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::sharp(), TEST_CASE(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Tf().
|
inline |
Corrected normal vector of a face.
f | the face |
Definition at line 416 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceNormal().
Referenced by initQuantities(), initQuantitiesCached(), and TEST_CASE().
|
inline |
Flat operator for the face.
f | the face |
Definition at line 462 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::checkCache(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::E(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceNormal(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::FLAT_, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCache, and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::setInCache().
Referenced by initQuantities(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::P(), and TEST_CASE().
|
inline |
Definition at line 828 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::nbVertices().
Referenced by computeLaplace(), and TEST_CASE().
|
inline |
Generic method to compute all the per face DenseMatrices and store them in an indexed container.
Usage example:
perFaceOperator | the per face operator |
Definition at line 999 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh, and DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbFaces().
Referenced by TEST_CASE().
|
inline |
Generic method to compute all the per face Vector and store them in an indexed container.
Usage example:
perFaceVectorOperator | the per face operator |
Definition at line 1023 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh, and DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbFaces().
Referenced by TEST_CASE().
|
inline |
Definition at line 1088 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh.
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::blockConnection(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Tf(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Tv().
|
inline |
Computes the global Connection-Laplace-Beltrami operator by accumulating the per face operators.
lambda | the regualrization parameter for the local Connection-Laplace-Beltrami operators |
Definition at line 928 of file PolygonalCalculus.h.
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().
|
inline |
Compute and returns the inverse of the global lumped mass matrix (diagonal matrix with Max's weights for each vertex).
Definition at line 902 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::globalLumpedMassMatrix().
|
inline |
Computes the global Laplace-Beltrami operator by assembling the per face operators.
lambda | the regularization parameter for the local Laplace-Beltrami operators |
Definition at line 855 of file PolygonalCalculus.h.
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(), TEST_CASE(), and TEST_CASE().
|
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) ;
Definition at line 882 of file PolygonalCalculus.h.
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().
|
inline |
Gradient operator of the face.
f | the face |
Definition at line 448 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::bracket(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::checkCache(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::coGradient(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceArea(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceNormal(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::GRAD_, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCache, and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::setInCache().
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantGradient(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantGradient_f(), initQuantities(), initQuantitiesCached(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::shapeOperator(), and TEST_CASE().
|
inline |
Definition at line 833 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::nbVertices().
|
inline |
Update the internal cache structures (e.g. degree of each face).
Definition at line 1055 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::updateFaceDegree().
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus().
|
inline |
Checks the validity/consistency of the object.
Definition at line 1111 of file PolygonalCalculus.h.
Referenced by TEST_CASE().
|
inlineprotected |
Definition at line 1255 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::M().
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantGradient_f(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantProjection_f().
|
inline |
(weak) Laplace-Beltrami operator for the face.
f | the face |
lambda | the regularization parameter |
Definition at line 602 of file PolygonalCalculus.h.
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().
|
inline |
Inner product on 1-forms associated with the face
f | the face |
lambda | the regularization parameter |
Definition at line 535 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::checkCache(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceArea(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::M_, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCache, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::P(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::setInCache(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::sharp().
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::divergence(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::doubledGlobalLumpedMassMatrix(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::globalLumpedMassMatrix(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::kroneckerWithI2(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::laplaceBeltrami().
|
inlineprotected |
Definition at line 1233 of file PolygonalCalculus.h.
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().
|
inline |
Definition at line 1075 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh, and DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbFaces().
|
inline |
Definition at line 1069 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh, and DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbVertices().
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::form0(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::identity0(), and TEST_CASE().
|
delete |
Deleted copy assignment operator.
other | the object to copy. |
|
delete |
Deleted move assignment operator.
other | the object to move. |
|
inline |
Projection operator for the face.
f | the face |
Definition at line 519 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::checkCache(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::flat(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myFaceDegree, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCache, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::P_, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::setInCache(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::sharp().
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::connectionLaplacian(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantProjection(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::covariantProjection_f(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::M(), and TEST_CASE().
|
inlinestaticprotected |
Project u on the orthgonal of n
u | vector to project |
n | vector to build orthogonal space from |
Definition at line 1167 of file PolygonalCalculus.h.
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Tf(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Tv().
|
inline |
https://math.stackexchange.com/questions/180418/calculate-rotation-matrix-to-align-vector-a-to-vector-b-in-3d@return 3x3 Rotation matrix to align n_v to n_f
Definition at line 688 of file PolygonalCalculus.h.
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().
|
inline |
Definition at line 706 of file PolygonalCalculus.h.
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().
|
inline |
Writes/Displays the object on an output stream.
out | the output stream where the object is written. |
Definition at line 1097 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCacheEnabled, and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh.
|
inline |
Update the embedding function.
externalFunctor | a new embedding functor (Face,Vertex)->RealPoint. |
Definition at line 280 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myEmbedder.
Referenced by precompute().
|
inlineprotected |
Set an operator in the internal cache.
key | the operator name |
f | the face |
ope | the operator to store |
Definition at line 1156 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCache, and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCacheEnabled.
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::A(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::B(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::coGradient(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::connectionLaplacian(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::curl(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::D(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::divergence(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::E(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::flat(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::gradient(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::laplaceBeltrami(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::M(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::P(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::sharp(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::X().
|
inline |
Shape operator on the face f (2x2 operator).
Definition at line 713 of file PolygonalCalculus.h.
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().
|
inline |
Sharp operator for the face.
f | the face |
Definition at line 503 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::B(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::bracket(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::centroid(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::checkCache(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceArea(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceNormal(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myFaceDegree, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCache, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::setInCache(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::SHARP_.
Referenced by initQuantities(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::M(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::P(), TEST_CASE(), and test_linear_structure().
|
inline |
Definition at line 644 of file PolygonalCalculus.h.
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().
|
inline |
toExtrinsicVector
v | the vertex |
I | the intrinsic vector at Tv |
Definition at line 668 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Tv().
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::toExtrinsicVectors().
|
inline |
I | set of intrinsic vectors, vectors indices must be the same as their associated vertex |
Definition at line 678 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh, DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbVertices(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::toExtrinsicVector().
|
inlinestaticprotected |
toReal3dVector converts Eigen::Vector3d to Real3dVector.
Conversion routines.
x | the vector |
Definition at line 1196 of file PolygonalCalculus.h.
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus().
|
inlinestaticprotected |
toVec3 convert Real3dPoint to Eigen::Vector3d
x | the vector |
Definition at line 1187 of file PolygonalCalculus.h.
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::n_v(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Tf(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::Tv().
|
inlinestaticprotected |
toVector convert Real3dPoint to Eigen::VectorXd
Conversion routines.
x | the vector |
Definition at line 1176 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::X().
|
inline |
Definition at line 734 of file PolygonalCalculus.h.
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().
|
inline |
Definition at line 624 of file PolygonalCalculus.h.
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().
|
inlineprotected |
Update the face degree cache.
Definition at line 1129 of file PolygonalCalculus.h.
References DGtal::SurfaceMesh< TRealPoint, TRealVector >::incidentVertices(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myFaceDegree, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh, and DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbFaces().
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::init().
|
inline |
Polygonal (corrected) vector area.
f | the face |
Definition at line 374 of file PolygonalCalculus.h.
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().
|
inline |
Return the vertex position matrix degree x 3 of the face.
f | a face |
Definition at line 294 of file PolygonalCalculus.h.
References DGtal::PolygonalCalculus< TRealPoint, TRealVector >::checkCache(), DGtal::SurfaceMesh< TRealPoint, TRealVector >::incidentVertices(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myEmbedder, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myFaceDegree, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::myGlobalCache, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::mySurfaceMesh, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::setInCache(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::X_.
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::B(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::centroid(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::E(), TEST_CASE(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::toVector().
|
static |
Concept checking.
Definition at line 134 of file PolygonalCalculus.h.
|
private |
Embedding function (face,vertex)->R^3 for the vertex position wrt. the face.
Definition at line 1281 of file PolygonalCalculus.h.
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::setEmbedder(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::vectorArea(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::X().
|
private |
Cache containing the face degree.
Definition at line 1287 of file PolygonalCalculus.h.
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::A(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::centroid(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::curl(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::D(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::degree(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::doubledGlobalLumpedMassMatrix(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceDegree(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::globalLaplaceBeltrami(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::globalLumpedMassMatrix(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::P(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::shapeOperator(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::sharp(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::transportAndFormatVectorField(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::updateFaceDegree(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::X().
|
mutableprivate |
Definition at line 1291 of file PolygonalCalculus.h.
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::A(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::B(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::checkCache(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::coGradient(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::connectionLaplacian(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::curl(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::D(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::disableInternalGlobalCache(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::divergence(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::E(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::flat(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::gradient(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::laplaceBeltrami(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::M(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::P(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::setInCache(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::sharp(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::X().
|
private |
Global cache.
Definition at line 1290 of file PolygonalCalculus.h.
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::checkCache(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::disableInternalGlobalCache(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::enableInternalGlobalCache(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::selfDisplay(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::setInCache().
|
private |
Underlying SurfaceMesh.
Definition at line 1278 of file PolygonalCalculus.h.
Referenced by DGtal::PolygonalCalculus< TRealPoint, TRealVector >::computeVertexNormal(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::doubledGlobalLumpedMassMatrix(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::getOperatorCacheMatrix(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::getOperatorCacheVector(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::getSurfaceMeshPtr(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::globalConnectionLaplace(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::globalLaplaceBeltrami(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::globalLumpedMassMatrix(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::nbFaces(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::nbVertices(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::selfDisplay(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::shapeOperator(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::toExtrinsicVectors(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::transportAndFormatVectorField(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::updateFaceDegree(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::vectorArea(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::X().
|
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(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus(), and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::PolygonalCalculus().