DGtal  1.5.beta
dgtalCalculus-bunny.cpp File Reference
#include <iostream>
#include <string>
#include <DGtal/base/Common.h>
#include <DGtal/helpers/StdDefs.h>
#include <DGtal/helpers/Shortcuts.h>
#include <DGtal/helpers/ShortcutsGeometry.h>
#include <DGtal/shapes/SurfaceMesh.h>
#include <DGtal/geometry/surfaces/DigitalSurfaceRegularization.h>
#include <DGtal/dec/PolygonalCalculus.h>
#include <DGtal/math/linalg/DirichletConditions.h>
#include <polyscope/polyscope.h>
#include <polyscope/surface_mesh.h>
#include <polyscope/point_cloud.h>
#include <Eigen/Dense>
#include <Eigen/Sparse>
#include "ConfigExamples.h"
Include dependency graph for dgtalCalculus-bunny.cpp:

Go to the source code of this file.

Typedefs

typedef Shortcuts< Z3i::KSpaceSH3
 
typedef ShortcutsGeometry< Z3i::KSpaceSHG3
 
typedef SurfaceMesh< RealPoint, RealVectorSurfMesh
 
typedef SurfMesh::Face Face
 
typedef SurfMesh::Vertex Vertex
 

Functions

double phiVertex (const Vertex v)
 
PolygonalCalculus< SH3::RealPoint, SH3::RealVector >::Vector phi (const Face f)
 
void initPhi ()
 
void initQuantities ()
 
void initQuantitiesCached ()
 
void computeLaplace ()
 
void myCallback ()
 
int main ()
 

Variables

polyscope::SurfaceMesh * psMesh
 
SurfMesh surfmesh
 
float scale = 0.1
 
PolygonalCalculus< SH3::RealPoint, SH3::RealVector >::Vector phiEigen
 

Detailed Description

This program is free software: you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/.

Author
David Coeurjolly (david.nosp@m..coe.nosp@m.urjol.nosp@m.ly@l.nosp@m.iris..nosp@m.cnrs.nosp@m..fr ) Laboratoire d'InfoRmatique en Image et Systemes d'information - LIRIS (CNRS, UMR 5205), CNRS, France
Date
2021/09/02

This file is part of the DGtal library.

Definition in file dgtalCalculus-bunny.cpp.

Typedef Documentation

◆ Face

Definition at line 53 of file dgtalCalculus-bunny.cpp.

◆ SH3

Definition at line 49 of file dgtalCalculus-bunny.cpp.

◆ SHG3

Definition at line 50 of file dgtalCalculus-bunny.cpp.

◆ SurfMesh

◆ Vertex

Definition at line 54 of file dgtalCalculus-bunny.cpp.

Function Documentation

◆ computeLaplace()

void computeLaplace ( )

Definition at line 167 of file dgtalCalculus-bunny.cpp.

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 }
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.
SparseMatrix globalLaplaceBeltrami(const double lambda=1.0) const
LinAlg::SolverSimplicialLDLT Solver
Type of a sparse matrix solver.
LinAlg::DenseVector Vector
Type of Vector.
void beginBlock(const std::string &keyword="")
double endBlock()
polyscope::SurfaceMesh * psMesh
SurfMesh surfmesh
PolyCalculus * calculus
Trace trace
Definition: Common.h:153
Size nbVertices() const
Definition: SurfaceMesh.h:288
EigenLinearAlgebraBackend::SparseMatrix SparseMatrix

References DGtal::Trace::beginBlock(), calculus, DGtal::Trace::endBlock(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::globalLaplaceBeltrami(), DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbVertices(), psMesh, surfmesh, and DGtal::trace.

Referenced by myCallback().

◆ initPhi()

void initPhi ( )

Definition at line 83 of file dgtalCalculus-bunny.cpp.

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 }
PolygonalCalculus< SH3::RealPoint, SH3::RealVector >::Vector phiEigen
double phiVertex(const Vertex v)

References DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbVertices(), phiEigen, phiVertex(), psMesh, and surfmesh.

Referenced by myCallback().

◆ initQuantities()

void initQuantities ( )

Definition at line 91 of file dgtalCalculus-bunny.cpp.

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 }
Real3dVector faceNormalAsDGtalVector(const Face f) const
Vector vectorArea(const Face f) const
double faceArea(const Face f) const
DenseMatrix coGradient(const Face f) const
DenseMatrix gradient(const Face f) const
PolygonalCalculus< SH3::RealPoint, SH3::RealVector >::Vector phi(const Face f)
DigitalPlane::Point Vector
Size nbFaces() const
Definition: SurfaceMesh.h:296

References DGtal::Trace::beginBlock(), calculus, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::coGradient(), DGtal::Trace::endBlock(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceArea(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceNormalAsDGtalVector(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::gradient(), DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbFaces(), phi(), psMesh, surfmesh, DGtal::trace, and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::vectorArea().

Referenced by myCallback().

◆ initQuantitiesCached()

void initQuantitiesCached ( )

Definition at line 129 of file dgtalCalculus-bunny.cpp.

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 }

References DGtal::Trace::beginBlock(), calculus, DGtal::PolygonalCalculus< TRealPoint, TRealVector >::coGradient(), DGtal::Trace::endBlock(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceArea(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::faceNormalAsDGtalVector(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::gradient(), DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbFaces(), phi(), psMesh, surfmesh, DGtal::trace, and DGtal::PolygonalCalculus< TRealPoint, TRealVector >::vectorArea().

Referenced by myCallback().

◆ main()

int main ( void  )

Definition at line 224 of file dgtalCalculus-bunny.cpp.

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 }
std::map< Cell, IdxVertex > Cell2Index
Definition: Shortcuts.h:189
SurfaceMesh< RealPoint, RealVector > SurfMesh
void myCallback()
CountedPtr< SH3::DigitalSurface > surface
CountedPtr< SH3::BinaryImage > binary_image
Edges computeNonManifoldEdges() const
K init(Point(0, 0, 0), Point(512, 512, 512), true)
KSpace K

References binary_image, DGtal::SurfaceMesh< TRealPoint, TRealVector >::computeNonManifoldEdges(), init(), K, myCallback(), psMesh, surface, and surfmesh.

◆ myCallback()

void myCallback ( )

Definition at line 207 of file dgtalCalculus-bunny.cpp.

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 }
void computeLaplace()
void initPhi()
void initQuantitiesCached()
void initQuantities()

References computeLaplace(), initPhi(), initQuantities(), and initQuantitiesCached().

Referenced by main().

◆ phi()

Examples
dec/exampleDECSurface.cpp, and geometry/volumes/checkFullConvexityTheorems.cpp.

Definition at line 69 of file dgtalCalculus-bunny.cpp.

70 {
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 }
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::SurfaceMesh< TRealPoint, TRealVector >::incidentVertices(), phiVertex(), and surfmesh.

Referenced by checkFullConvexityCharacterization(), checkProjectionFullConvexity(), initQuantities(), initQuantitiesCached(), DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::randomRotation(), and TEST_CASE().

◆ phiVertex()

double phiVertex ( const Vertex  v)

Definition at line 63 of file dgtalCalculus-bunny.cpp.

64 {
65  return cos(scale*(surfmesh.position(v)[0]))*sin(scale*surfmesh.position(v)[1]);
66 }
RealPoint & position(Vertex v)
Definition: SurfaceMesh.h:647

References DGtal::SurfaceMesh< TRealPoint, TRealVector >::position(), and surfmesh.

Referenced by initPhi(), and phi().

Variable Documentation

◆ phiEigen

Definition at line 60 of file dgtalCalculus-bunny.cpp.

Referenced by initPhi().

◆ psMesh

polyscope::SurfaceMesh* psMesh

◆ scale

float scale = 0.1

Definition at line 59 of file dgtalCalculus-bunny.cpp.

Referenced by testFFTScaling().

◆ surfmesh