DGtal  1.5.beta
dgtalCalculus-bunny.cpp
Go to the documentation of this file.
1 
26 #include <iostream>
27 #include <string>
28 #include <DGtal/base/Common.h>
29 #include <DGtal/helpers/StdDefs.h>
30 #include <DGtal/helpers/Shortcuts.h>
31 #include <DGtal/helpers/ShortcutsGeometry.h>
32 #include <DGtal/shapes/SurfaceMesh.h>
33 #include <DGtal/geometry/surfaces/DigitalSurfaceRegularization.h>
34 #include <DGtal/dec/PolygonalCalculus.h>
35 #include <DGtal/math/linalg/DirichletConditions.h>
36 
37 #include <polyscope/polyscope.h>
38 #include <polyscope/surface_mesh.h>
39 #include <polyscope/point_cloud.h>
40 #include <Eigen/Dense>
41 #include <Eigen/Sparse>
42 
43 #include "ConfigExamples.h"
44 
45 using namespace DGtal;
46 using namespace Z3i;
47 
48 // Using standard 3D digital space.
51 // The following typedefs are useful
55 
56 //Polyscope global
57 polyscope::SurfaceMesh *psMesh;
59 float scale = 0.1;
61 
62 //Restriction of a scalar function to vertices
63 double phiVertex(const Vertex v)
64 {
65  return cos(scale*(surfmesh.position(v)[0]))*sin(scale*surfmesh.position(v)[1]);
66 }
67 
68 //Restriction of a scalar function to vertices
70 {
71  auto vertices = surfmesh.incidentVertices(f);
72  auto nf = vertices.size();
73  Eigen::VectorXd ph(nf);
74  size_t cpt=0;
75  for(auto v: vertices)
76  {
77  ph(cpt) = phiVertex(v);
78  ++cpt;
79  }
80  return ph;
81 }
82 
83 void initPhi()
84 {
85  phiEigen.resize(surfmesh.nbVertices());
86  for(auto i = 0; i < surfmesh.nbVertices(); ++i)
87  phiEigen(i) = phiVertex(i);
88  psMesh->addVertexScalarQuantity("Phi", phiEigen);
89 }
90 
92 {
93  trace.beginBlock("Basic quantities");
95 
98  std::vector<PolygonalCalculus<SH3::RealPoint,SH3::RealVector>::Real3dPoint> normals;
99  std::vector<PolygonalCalculus<SH3::RealPoint,SH3::RealVector>::Real3dPoint> vectorArea;
100  std::vector<PolygonalCalculus<SH3::RealPoint,SH3::RealVector>::Real3dPoint> centroids;
101 
102  std::vector<double> faceArea;
103 
104  for(auto f=0; f < surfmesh.nbFaces(); ++f)
105  {
107  gradients.push_back( grad );
108 
110  cogradients.push_back( cograd );
111 
112  normals.push_back(calculus.faceNormalAsDGtalVector(f));
113 
115  vectorArea.push_back({vA(0) , vA(1), vA(2)});
116 
117  faceArea.push_back( calculus.faceArea(f));
118  }
119  trace.endBlock();
120 
121  psMesh->addFaceVectorQuantity("Gradients", gradients);
122  psMesh->addFaceVectorQuantity("co-Gradients", cogradients);
123  psMesh->addFaceVectorQuantity("Normals", normals);
124  psMesh->addFaceScalarQuantity("Face area", faceArea);
125  psMesh->addFaceVectorQuantity("Vector area", vectorArea);
126 }
127 
128 
130 {
131  trace.beginBlock("Basic quantities (cached)");
133 
136  std::vector<PolygonalCalculus<SH3::RealPoint,SH3::RealVector>::Real3dPoint> normals;
137  std::vector<PolygonalCalculus<SH3::RealPoint,SH3::RealVector>::Real3dPoint> vectorArea;
138  std::vector<PolygonalCalculus<SH3::RealPoint,SH3::RealVector>::Real3dPoint> centroids;
139 
140  std::vector<double> faceArea;
141 
142  for(auto f=0; f < surfmesh.nbFaces(); ++f)
143  {
145  gradients.push_back( grad );
146 
148  cogradients.push_back( cograd );
149 
150  normals.push_back(calculus.faceNormalAsDGtalVector(f));
151 
153  vectorArea.push_back({vA(0) , vA(1), vA(2)});
154 
155  faceArea.push_back( calculus.faceArea(f));
156  }
157  trace.endBlock();
158 
159  psMesh->addFaceVectorQuantity("Gradients", gradients);
160  psMesh->addFaceVectorQuantity("co-Gradients", cogradients);
161  psMesh->addFaceVectorQuantity("Normals", normals);
162  psMesh->addFaceScalarQuantity("Face area", faceArea);
163  psMesh->addFaceVectorQuantity("Vector area", vectorArea);
164 }
165 
166 
168 {
170  trace.beginBlock("Operator construction...");
172  trace.endBlock();
173 
174  const auto nbv = surfmesh.nbVertices();
176 
177  //Setting some random sources
179  DC::IntegerVector p = DC::nullBoundaryVector( g );
180  for(auto cpt=0; cpt< 10;++cpt)
181  {
182  int idx = rand() % nbv;
183  g( idx ) = rand() % 100;
184  p( idx ) = 1.0;
185  }
186 
187  //Solve Δu=0 with g as boundary conditions
189 
190  trace.beginBlock("Prefactorization...");
191  DC::SparseMatrix L_dirichlet = DC::dirichletOperator( L, p );
192  solver.compute( L_dirichlet );
193  ASSERT(solver.info()==Eigen::Success);
194  trace.endBlock();
195 
196  trace.beginBlock("Solve...");
197  PolygonalCalculus<SH3::RealPoint,SH3::RealVector>::Vector g_dirichlet = DC::dirichletVector( L, g, p, g );
198  PolygonalCalculus<SH3::RealPoint,SH3::RealVector>::Vector x_dirichlet = solver.solve( g_dirichlet );
199  PolygonalCalculus<SH3::RealPoint,SH3::RealVector>::Vector u = DC::dirichletSolution( x_dirichlet, p, g );
200  ASSERT(solver.info()==Eigen::Success);
201  trace.endBlock();
202 
203  psMesh->addVertexScalarQuantity("g", g);
204  psMesh->addVertexScalarQuantity("u", u);
205 }
206 
208 {
209  ImGui::SliderFloat("Phi scale", &scale, 0., 1.);
210  if (ImGui::Button("Phi and basic operators"))
211  {
212  initPhi();
213  initQuantities();
214  }
215  if (ImGui::Button("Phi and basic operators (cached)"))
216  {
217  initPhi();
219  }
220  if(ImGui::Button("Compute Laplace problem"))
221  computeLaplace();
222 }
223 
224 int main()
225 {
226  auto params = SH3::defaultParameters() | SHG3::defaultParameters() | SHG3::parametersGeometryEstimation();
227  params("surfaceComponents", "All");
228 
229  std::string filename = examplesPath + std::string("/samples/bunny-32.vol");
230  auto binary_image = SH3::makeBinaryImage(filename, params );
231  auto K = SH3::getKSpace( binary_image, params );
232  auto surface = SH3::makeDigitalSurface( binary_image, K, params );
233  SH3::Cell2Index c2i;
234  auto primalSurface = SH3::makePrimalSurfaceMesh(surface);
235 
236  //Need to convert the faces
237  std::vector<std::vector<SH3::SurfaceMesh::Vertex>> faces;
238  std::vector<RealPoint> positions;
239 
240  for(auto face= 0 ; face < primalSurface->nbFaces(); ++face)
241  faces.push_back(primalSurface->incidentVertices( face ));
242 
243  //Recasting to vector of vertices
244  positions = primalSurface->positions();
245 
246  surfmesh = SurfMesh(positions.begin(),
247  positions.end(),
248  faces.begin(),
249  faces.end());
250 
251  std::cout<<"number of non-manifold Edges = " << surfmesh.computeNonManifoldEdges().size()<<std::endl;
252  psMesh = polyscope::registerSurfaceMesh("digital surface", positions, faces);
253 
254  // Initialize polyscope
255  polyscope::init();
256 
257  // Set the callback function
258  polyscope::state::userCallback = myCallback;
259  polyscope::show();
260  return EXIT_SUCCESS;
261 }
Aim: A helper class to solve a system with Dirichlet boundary conditions.
Implements differential operators on polygonal surfaces from .
LinAlg::SparseMatrix SparseMatrix
Type of sparse matrix.
Real3dVector faceNormalAsDGtalVector(const Face f) const
Vector vectorArea(const Face f) const
SparseMatrix globalLaplaceBeltrami(const double lambda=1.0) const
LinAlg::SolverSimplicialLDLT Solver
Type of a sparse matrix solver.
double faceArea(const Face f) const
DenseMatrix coGradient(const Face f) const
DenseMatrix gradient(const Face f) const
LinAlg::DenseVector Vector
Type of Vector.
Aim: This class is used to simplify shape and surface creation. With it, you can create new shapes an...
Aim: This class is used to simplify shape and surface creation. With it, you can create new shapes an...
Definition: Shortcuts.h:105
std::map< Cell, IdxVertex > Cell2Index
Definition: Shortcuts.h:189
void beginBlock(const std::string &keyword="")
double endBlock()
SurfaceMesh< RealPoint, RealVector > SurfMesh
PolygonalCalculus< SH3::RealPoint, SH3::RealVector >::Vector phiEigen
double phiVertex(const Vertex v)
SurfMesh::Face Face
void computeLaplace()
void initPhi()
void initQuantitiesCached()
void myCallback()
polyscope::SurfaceMesh * psMesh
SurfMesh::Vertex Vertex
PolygonalCalculus< SH3::RealPoint, SH3::RealVector >::Vector phi(const Face f)
Shortcuts< Z3i::KSpace > SH3
SurfMesh surfmesh
ShortcutsGeometry< Z3i::KSpace > SHG3
void initQuantities()
int main()
PolyCalculus * calculus
CountedPtr< SH3::DigitalSurface > surface
CountedPtr< SH3::BinaryImage > binary_image
DigitalPlane::Point Vector
DGtal is the top-level namespace which contains all DGtal functions and types.
Trace trace
Definition: Common.h:153
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
const Vertices & incidentVertices(Face f) const
Definition: SurfaceMesh.h:315
Size nbFaces() const
Definition: SurfaceMesh.h:296
Size nbVertices() const
Definition: SurfaceMesh.h:288
Edges computeNonManifoldEdges() const
K init(Point(0, 0, 0), Point(512, 512, 512), true)
KSpace K
EigenLinearAlgebraBackend::SparseMatrix SparseMatrix
TriMesh::Face Face
TriMesh::Vertex Vertex