DGtal  1.5.beta
PolygonalCalculus.h
1 
17 #pragma once
18 
31 // Inclusions
32 #include <iostream>
33 #include <functional>
34 #include <vector>
35 #include <string>
36 #include <map>
37 #include <unordered_map>
38 #include "DGtal/base/ConstAlias.h"
39 #include "DGtal/base/Common.h"
40 #include "DGtal/shapes/SurfaceMesh.h"
41 #include "DGtal/math/linalg/EigenSupport.h"
43 
44 namespace DGtal
45 {
46 
47  namespace functors { template <typename TRealPoint, typename TRealVector>
62  {
66  typedef typename MySurfaceMesh::Vertex Vertex;
68  typedef typename MySurfaceMesh::Face Face;
71 
73 
77  EmbedderFromNormalVectors(ConstAlias<std::vector<Real3dPoint>> normals,
79  {
80  myNormals = &normals;
82  }
83 
89  Real3dPoint operator()(const Face &f,const Vertex &v)
90  {
91  const auto nn = (*myNormals)[f];
93  return p - nn.dot(p)*nn;
94  }
95 
97  const std::vector<Real3dPoint> *myNormals;
100  };
101  }
102 
103 
104 
106 // template class PolygonalCalculus
127 template <typename TRealPoint, typename TRealVector>
129 {
130  // ----------------------- Standard services ------------------------------
131 public:
132 
134  static const Dimension dimension = TRealPoint::dimension;
136 
139 
143  typedef typename MySurfaceMesh::Vertex Vertex;
145  typedef typename MySurfaceMesh::Face Face;
150 
156  typedef Vector Form;
163 
166 
169 
175  bool globalInternalCacheEnabled = false):
176  mySurfaceMesh(&surf), myGlobalCacheEnabled(globalInternalCacheEnabled)
177  {
178  myEmbedder =[&](Face f,Vertex v){ (void)f; return mySurfaceMesh->position(v);};
180  init();
181  };
182 
190  const std::function<Real3dPoint(Face,Vertex)> &embedder,
191  bool globalInternalCacheEnabled = false):
192  mySurfaceMesh(&surf), myEmbedder(embedder), myGlobalCacheEnabled(globalInternalCacheEnabled)
193  {
195  init();
196  };
197 
205  const std::function<Real3dVector(Vertex)> &embedder,
206  bool globalInternalCacheEnabled = false):
207  mySurfaceMesh(&surf), myVertexNormalEmbedder(embedder),
208  myGlobalCacheEnabled(globalInternalCacheEnabled)
209  {
210  myEmbedder = [&](Face f,Vertex v){(void)f; return mySurfaceMesh->position(v); };
211  init();
212  };
213 
225  const std::function<Real3dPoint(Face,Vertex)> &pos_embedder,
226  const std::function<Vector(Vertex)> &normal_embedder,
227  bool globalInternalCacheEnabled = false) :
228  mySurfaceMesh(&surf), myEmbedder(pos_embedder),
229  myVertexNormalEmbedder(normal_embedder),
230  myGlobalCacheEnabled(globalInternalCacheEnabled)
231  {
232  init();
233  };
234 
238  PolygonalCalculus() = delete;
239 
243  ~PolygonalCalculus() = default;
244 
249  PolygonalCalculus ( const PolygonalCalculus & other ) = delete;
250 
255  PolygonalCalculus ( PolygonalCalculus && other ) = delete;
256 
262  PolygonalCalculus & operator= ( const PolygonalCalculus & other ) = delete;
263 
270 
272 
273  // ----------------------- embedding services --------------------------
274  //MARK: Embedding services
277 
280  void setEmbedder(const std::function<Real3dPoint(Face,Vertex)> &externalFunctor)
281  {
282  myEmbedder = externalFunctor;
283  }
285 
286  // ----------------------- Per face operators --------------------------------------
287  //MARK: Per face operator on scalars
290 
294  DenseMatrix X(const Face f) const
295  {
296  if (checkCache(X_,f))
297  return myGlobalCache[X_][f];
298 
299  const auto vertices = mySurfaceMesh->incidentVertices(f);
300  const auto nf = myFaceDegree[f];
301  DenseMatrix Xt(nf,3);
302  size_t cpt=0;
303  for(auto v: vertices)
304  {
305  Xt(cpt,0) = myEmbedder(f,v)[0];
306  Xt(cpt,1) = myEmbedder(f,v)[1];
307  Xt(cpt,2) = myEmbedder(f,v)[2];
308  ++cpt;
309  }
310 
311  setInCache(X_,f,Xt);
312  return Xt;
313  }
314 
315 
319  DenseMatrix D(const Face f) const
320  {
321  if (checkCache(D_,f))
322  return myGlobalCache[D_][f];
323 
324  const auto nf = myFaceDegree[f];
325  DenseMatrix d = DenseMatrix::Zero(nf ,nf);
326  for(auto i=0u; i < nf; ++i)
327  {
328  d(i,i) = -1.;
329  d(i, (i+1)%nf) = 1.;
330  }
331 
332  setInCache(D_,f,d);
333  return d;
334  }
335 
339  DenseMatrix E(const Face f) const
340  {
341  if (checkCache(E_,f))
342  return myGlobalCache[E_][f];
343 
344  DenseMatrix op = D(f)*X(f);
345 
346  setInCache(E_,f,op);
347  return op;
348  }
349 
353  DenseMatrix A(const Face f) const
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  }
369 
370 
374  Vector vectorArea(const Face f) const
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  }
394 
398  double faceArea(const Face f) const
399  {
400  return vectorArea(f).norm();
401  }
402 
406  Vector faceNormal(const Face f) const
407  {
408  Vector v = vectorArea(f);
409  v.normalize();
410  return v;
411  }
412 
417  {
418  Vector v = faceNormal(f);
419  return {v(0),v(1),v(2)};
420  }
421 
425  DenseMatrix coGradient(const Face f) const
426  {
427  if (checkCache(COGRAD_,f))
428  return myGlobalCache[COGRAD_][f];
429  DenseMatrix op = E(f).transpose() * A(f);
430  setInCache(COGRAD_, f, op);
431  return op;
432  }
433 
436  DenseMatrix bracket(const Vector &n) const
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  }
444 
448  DenseMatrix gradient(const Face f) const
449  {
450  if (checkCache(GRAD_,f))
451  return myGlobalCache[GRAD_][f];
452 
453  DenseMatrix op = -1.0/faceArea(f) * bracket( faceNormal(f) ) * coGradient(f);
454 
455  setInCache(GRAD_,f,op);
456  return op;
457  }
458 
462  DenseMatrix flat(const Face f) const
463  {
464  if (checkCache(FLAT_,f))
465  return myGlobalCache[FLAT_][f];
466  DenseMatrix n = faceNormal(f);
467  DenseMatrix op = E(f)*( DenseMatrix::Identity(3,3) - n*n.transpose());
468  setInCache(FLAT_,f,op);
469  return op;
470  }
471 
475  DenseMatrix B(const Face f) const
476  {
477  if (checkCache(B_,f))
478  return myGlobalCache[B_][f];
479  DenseMatrix res = A(f) * X(f);
480  setInCache(B_,f,res);
481  return res;
482  }
483 
486  Vector centroid(const Face f) const
487  {
488  const auto nf = myFaceDegree[f];
489  return 1.0/(double)nf * X(f).transpose() * Vector::Ones(nf);
490  }
491 
495  {
496  const Vector c = centroid(f);
497  return {c(0),c(1),c(2)};
498  }
499 
503  DenseMatrix sharp(const Face f) const
504  {
505  if (checkCache(SHARP_,f))
506  return myGlobalCache[SHARP_][f];
507 
508  const auto nf = myFaceDegree[f];
509  DenseMatrix op = 1.0/faceArea(f) * bracket(faceNormal(f)) *
510  ( B(f).transpose() - centroid(f)* Vector::Ones(nf).transpose() );
511 
512  setInCache(SHARP_,f,op);
513  return op;
514  }
515 
519  DenseMatrix P(const Face f) const
520  {
521  if (checkCache(P_,f))
522  return myGlobalCache[P_][f];
523 
524  const auto nf = myFaceDegree[f];
525  DenseMatrix op = DenseMatrix::Identity(nf,nf) - flat(f)*sharp(f);
526 
527  setInCache(P_, f, op);
528  return op;
529  }
530 
535  DenseMatrix M(const Face f, const double lambda=1.0) const
536  {
537  if (checkCache(M_,f))
538  return myGlobalCache[M_][f];
539 
540  DenseMatrix Uf = sharp(f);
541  DenseMatrix Pf = P(f);
542  DenseMatrix op = faceArea(f) * Uf.transpose()*Uf + lambda * Pf.transpose()*Pf;
543 
544  setInCache(M_,f,op);
545  return op;
546  }
547 
562  DenseMatrix divergence(const Face f) const
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  }
572 
576  DenseMatrix curl(const Face f) const
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  }
586 
587 
602  DenseMatrix laplaceBeltrami(const Face f, const double lambda=1.0) const
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  }
615 
616  // ----------------------- Vector calculus----------------------------------
617  //MARK: Vector Field Calculus
620 
621 public:
624  DenseMatrix Tv(const Vertex & v) const
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  }
641 
644  DenseMatrix Tf(const Face & f) const
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  }
662 
668  Vector toExtrinsicVector(const Vertex v, const Vector & I) const
669  {
670  DenseMatrix T = Tv(v);
671  return T.col(0) * I(0) + T.col(1) * I(1);
672  }
673 
678  std::vector<Vector> toExtrinsicVectors(const std::vector<Vector> & I) const
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  }
685 
688  DenseMatrix Qvf(const Vertex & v, const Face & f) const
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  }
703 
706  DenseMatrix Rvf(const Vertex & v, const Face & f) const
707  {
708  return Tf(f).transpose() * Qvf(v,f) * Tv(v);
709  }
710 
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  }
726 
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  }
746 
757  {
758  return Tf(f).transpose() * gradient(f) *
760  }
761 
773  {
774  return P(f) * D(f) * transportAndFormatVectorField(f,uf);
775  }
776 
783  {
784  return kroneckerWithI2(Tf(f).transpose() * gradient(f)) * blockConnection(f);
785  }
786 
793  {
794  return kroneckerWithI2(P(f) * D(f)) * blockConnection(f);
795  ;
796  }
797 
810  DenseMatrix connectionLaplacian(const Face & f, double lambda = 1.0) const
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  }
821 
822  // ----------------------- Global operators --------------------------------------
823  //MARK: Global Operators
826 
828  Form form0() const
829  {
830  return Form::Zero( nbVertices() );
831  }
834  {
835  SparseMatrix Id0( nbVertices(), nbVertices() );
836  Id0.setIdentity();
837  return Id0;
838  }
839 
855  SparseMatrix globalLaplaceBeltrami(const double lambda=1.0) const
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  }
876 
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  }
897 
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  }
910 
928  SparseMatrix globalConnectionLaplace(const double lambda = 1.0) const
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  }
953 
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  }
978 
979  // ----------------------- Cache mechanism --------------------------------------
982 
999  std::vector<DenseMatrix> getOperatorCacheMatrix(const std::function<DenseMatrix(Face)> &perFaceOperator) const
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  }
1006 
1023  std::vector<Vector> getOperatorCacheVector(const std::function<Vector(Face)> &perFaceVectorOperator) const
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  }
1030 
1034  {
1035  myGlobalCacheEnabled = true;
1036  }
1037 
1041  {
1042  myGlobalCacheEnabled = false;
1043  myGlobalCache.clear();
1044  }
1045 
1047 
1048  // ----------------------- Common --------------------------------------
1049 public:
1052 
1055  void init()
1056  {
1057  updateFaceDegree();
1058  }
1059 
1063  size_t faceDegree(Face f) const
1064  {
1065  return myFaceDegree[f];
1066  }
1067 
1069  size_t nbVertices() const
1070  {
1071  return mySurfaceMesh->nbVertices();
1072  }
1073 
1075  size_t nbFaces() const
1076  {
1077  return mySurfaceMesh->nbFaces();
1078  }
1079 
1082  size_t degree(const Face f) const
1083  {
1084  return myFaceDegree[f];
1085  }
1086 
1089  {
1090  return mySurfaceMesh;
1091  }
1092 
1097  void selfDisplay ( std::ostream & out ) const
1098  {
1099  out << "[PolygonalCalculus]: ";
1101  out<< "internal cache enabled, ";
1102  else
1103  out<<"internal cache disabled, ";
1104  out <<"SurfaceMesh="<<*mySurfaceMesh;
1105  }
1106 
1111  bool isValid() const
1112  {
1113  return true;
1114  }
1115 
1117 
1118  // ------------------------- Protected Datas ------------------------------
1119  //MARK: Protected
1120 
1121 protected:
1124 
1127 
1130  {
1131  myFaceDegree.resize(mySurfaceMesh->nbFaces());
1132  for(typename MySurfaceMesh::Index f = 0; f < mySurfaceMesh->nbFaces(); ++f)
1133  {
1134  auto vertices = mySurfaceMesh->incidentVertices(f);
1135  auto nf = vertices.size();
1136  myFaceDegree[f] = nf;
1137  }
1138  }
1139 
1144  bool checkCache(OPERATOR key, const Face f) const
1145  {
1147  if (myGlobalCache[key].find(f) != myGlobalCache[key].end())
1148  return true;
1149  return false;
1150  }
1151 
1156  void setInCache(OPERATOR key, const Face f,
1157  const DenseMatrix &ope) const
1158  {
1160  myGlobalCache[key][f] = ope;
1161  }
1162 
1167  static Vector project(const Vector & u, const Vector & n)
1168  {
1169  return u - (u.dot(n) / n.squaredNorm()) * n;
1170  }
1171 
1176  static Vector toVector(const Eigen::Vector3d & x)
1177  {
1178  Vector X(3);
1179  for (int i = 0; i < 3; i++)
1180  X(i) = x(i);
1181  return X;
1182  }
1183 
1187  static Eigen::Vector3d toVec3(const Real3dPoint & x)
1188  {
1189  return Eigen::Vector3d(x(0),x(1),x(2));
1190  }
1191 
1196  static Real3dVector toReal3dVector(const Eigen::Vector3d & x)
1197  {
1198  return { x(0), x(1), x(2)};
1199  }
1200 
1201 
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  }
1230 
1233  Eigen::Vector3d n_v(const Vertex & v) const
1234  {
1235  return toVec3(myVertexNormalEmbedder(v));
1236  }
1237 
1238  //Covariant operators routines
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  }
1253 
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  }
1268 
1269 
1270 
1272 
1273  // ------------------------- Internals ------------------------------------
1274  //MARK: Internals
1275 private:
1276 
1279 
1281  std::function<Real3dPoint(Face, Vertex)> myEmbedder;
1282 
1285 
1287  std::vector<size_t> myFaceDegree;
1288 
1291  mutable std::array<std::unordered_map<Face,DenseMatrix>, 15> myGlobalCache;
1292 
1293 }; // end of class PolygonalCalculus
1294 
1301 template <typename TP, typename TV>
1302 std::ostream&
1303 operator<< ( std::ostream & out, const PolygonalCalculus<TP,TV> & object )
1304 {
1305  object.selfDisplay( out );
1306  return out;
1307 }
1308 
1309 } // namespace DGtal
Aim: This class encapsulates its parameter class so that to indicate to the user that the object/poin...
Definition: ConstAlias.h:187
Implements differential operators on polygonal surfaces from .
Eigen::Vector3d n_v(const Vertex &v) const
DenseMatrix Rvf(const Vertex &v, const Face &f) const
PolygonalCalculus(const ConstAlias< MySurfaceMesh > surf, const std::function< Real3dPoint(Face, Vertex)> &pos_embedder, const std::function< Vector(Vertex)> &normal_embedder, bool globalInternalCacheEnabled=false)
EigenLinearAlgebraBackend LinAlg
Linear Algebra Backend from Eigen.
DenseMatrix transportAndFormatVectorField(const Face f, const Vector &uf)
SurfaceMesh< TRealPoint, TRealVector > MySurfaceMesh
Type of SurfaceMesh.
std::function< Real3dPoint(Face, Vertex)> myEmbedder
Embedding function (face,vertex)->R^3 for the vertex position wrt. the face.
DenseMatrix D(const Face f) const
DenseMatrix laplaceBeltrami(const Face f, const double lambda=1.0) const
static Eigen::Vector3d toVec3(const Real3dPoint &x)
toVec3 convert Real3dPoint to Eigen::Vector3d
DenseMatrix Tv(const Vertex &v) const
LinAlg::SparseMatrix SparseMatrix
Type of sparse matrix.
const MySurfaceMesh * getSurfaceMeshPtr() const
OPERATOR
Enum for operators in the internal cache strategy.
std::vector< size_t > myFaceDegree
Cache containing the face degree.
DenseMatrix X(const Face f) const
void selfDisplay(std::ostream &out) const
DenseMatrix Tf(const Face &f) const
PolygonalCalculus(PolygonalCalculus &&other)=delete
Vector faceNormal(const Face f) const
Vector centroid(const Face f) const
const MySurfaceMesh * mySurfaceMesh
Underlying SurfaceMesh.
MySurfaceMesh::RealVector Real3dVector
Real vector type.
DenseMatrix curl(const Face f) const
static const Dimension dimension
Concept checking.
Real3dVector faceNormalAsDGtalVector(const Face f) const
PolygonalCalculus(const ConstAlias< MySurfaceMesh > surf, const std::function< Real3dPoint(Face, Vertex)> &embedder, bool globalInternalCacheEnabled=false)
Real3dPoint centroidAsDGtalPoint(const Face f) const
DenseMatrix covariantProjection(const Face f, const Vector &uf)
Vector vectorArea(const Face f) const
DenseMatrix E(const Face f) const
DenseMatrix blockConnection(const Face &f) const
MySurfaceMesh::Vertex Vertex
Vertex type.
SparseMatrix globalInverseLumpedMassMatrix() const
Vector computeVertexNormal(const Vertex &v) const
std::array< std::unordered_map< Face, DenseMatrix >, 15 > myGlobalCache
size_t faceDegree(Face f) const
Vector toExtrinsicVector(const Vertex v, const Vector &I) const
toExtrinsicVector
std::vector< Vector > toExtrinsicVectors(const std::vector< Vector > &I) const
static Real3dVector toReal3dVector(const Eigen::Vector3d &x)
toReal3dVector converts Eigen::Vector3d to Real3dVector.
std::vector< DenseMatrix > getOperatorCacheMatrix(const std::function< DenseMatrix(Face)> &perFaceOperator) const
SparseMatrix doubledGlobalLumpedMassMatrix() const
bool myGlobalCacheEnabled
Global cache.
DenseMatrix bracket(const Vector &n) const
PolygonalCalculus(const PolygonalCalculus &other)=delete
DenseMatrix kroneckerWithI2(const DenseMatrix &M) const
DenseMatrix P(const Face f) const
SparseMatrix globalLaplaceBeltrami(const double lambda=1.0) const
PolygonalCalculus(const ConstAlias< MySurfaceMesh > surf, const std::function< Real3dVector(Vertex)> &embedder, bool globalInternalCacheEnabled=false)
DenseMatrix Qvf(const Vertex &v, const Face &f) const
void setEmbedder(const std::function< Real3dPoint(Face, Vertex)> &externalFunctor)
LinAlg::SolverSimplicialLDLT Solver
Type of a sparse matrix solver.
DenseMatrix sharp(const Face f) const
bool checkCache(OPERATOR key, const Face f) const
double faceArea(const Face f) const
DenseMatrix covariantGradient_f(const Face &f) const
DenseMatrix coGradient(const Face f) const
DenseMatrix divergence(const Face f) const
DenseMatrix M(const Face f, const double lambda=1.0) const
DenseMatrix shapeOperator(const Face f) const
DenseMatrix connectionLaplacian(const Face &f, double lambda=1.0) const
size_t degree(const Face f) const
void updateFaceDegree()
Update the face degree cache.
SparseMatrix globalLumpedMassMatrix() const
DenseMatrix gradient(const Face f) const
MySurfaceMesh::Face Face
Face type.
PolygonalCalculus< TRealPoint, TRealVector > Self
Self type.
LinAlg::DenseMatrix DenseMatrix
Type of dense matrix.
MySurfaceMesh::RealPoint Real3dPoint
Position type.
static Vector toVector(const Eigen::Vector3d &x)
toVector convert Real3dPoint to Eigen::VectorXd
void setInCache(OPERATOR key, const Face f, const DenseMatrix &ope) const
std::vector< Vector > getOperatorCacheVector(const std::function< Vector(Face)> &perFaceVectorOperator) const
LinAlg::Triplet Triplet
Type of sparse matrix triplet.
std::function< Real3dVector(Vertex)> myVertexNormalEmbedder
Embedding function (vertex)->R^3 for the vertex normal.
PolygonalCalculus(const ConstAlias< MySurfaceMesh > surf, bool globalInternalCacheEnabled=false)
DenseMatrix B(const Face f) const
Vector Form
Global 0-form, 1-form, 2-form are Vector.
PolygonalCalculus & operator=(const PolygonalCalculus &other)=delete
DenseMatrix covariantProjection_f(const Face &f) const
LinAlg::DenseVector Vector
Type of Vector.
SparseMatrix globalConnectionLaplace(const double lambda=1.0) const
static Vector project(const Vector &u, const Vector &n)
DenseMatrix A(const Face f) const
DenseMatrix covariantGradient(const Face f, const Vector &uf)
BOOST_STATIC_ASSERT((dimension==3))
DenseMatrix flat(const Face f) const
SparseMatrix identity0() const
std::ostream & warning()
SurfMesh surfmesh
DGtal is the top-level namespace which contains all DGtal functions and types.
std::ostream & operator<<(std::ostream &out, const ATu0v1< TKSpace, TLinearAlgebra > &object)
DGtal::uint32_t Dimension
Definition: Common.h:136
Trace trace
Definition: Common.h:153
Aim: Provide linear algebra backend using Eigen dense and sparse matrix as well as dense vector....
Definition: EigenSupport.h:96
Eigen::SimplicialLDLT< SparseMatrix > SolverSimplicialLDLT
Definition: EigenSupport.h:106
Eigen::Triplet< double, SparseMatrix::StorageIndex > Triplet
Definition: EigenSupport.h:102
Eigen::SparseMatrix< DenseVector::Scalar, Eigen::ColMajor, DenseVector::Index > SparseMatrix
Definition: EigenSupport.h:101
Aim: Represents an embedded mesh as faces and a list of vertices. Vertices may be shared among faces ...
Definition: SurfaceMesh.h:92
RealPoint & position(Vertex v)
Definition: SurfaceMesh.h:647
TRealPoint RealPoint
Definition: SurfaceMesh.h:93
std::size_t Index
The type used for numbering vertices and faces.
Definition: SurfaceMesh.h:105
const Vertices & incidentVertices(Face f) const
Definition: SurfaceMesh.h:315
Size nbFaces() const
Definition: SurfaceMesh.h:296
const Faces & incidentFaces(Vertex v) const
Definition: SurfaceMesh.h:321
Size nbVertices() const
Definition: SurfaceMesh.h:288
const Vertices & neighborVertices(Vertex v) const
Definition: SurfaceMesh.h:331
TRealVector RealVector
Definition: SurfaceMesh.h:94
Functor that projects a face vertex of a surface mesh onto the tangent plane given by a per-face norm...
const MySurfaceMesh * mySurfaceMesh
Alias to the surface mesh.
const std::vector< Real3dPoint > * myNormals
Alias to the normal vectors.
SurfaceMesh< TRealPoint, TRealVector > MySurfaceMesh
Type of SurfaceMesh.
MySurfaceMesh::RealPoint Real3dPoint
Position type.
MySurfaceMesh::Vertex Vertex
Vertex type.
Real3dPoint operator()(const Face &f, const Vertex &v)
EmbedderFromNormalVectors(ConstAlias< std::vector< Real3dPoint >> normals, ConstAlias< MySurfaceMesh > surfmesh)