DGtal  1.5.beta
dgtalCalculus-poisson.cpp
Go to the documentation of this file.
1 
26 #include <iostream>
27 #include <DGtal/base/Common.h>
28 #include <DGtal/helpers/StdDefs.h>
29 #include <DGtal/helpers/Shortcuts.h>
30 #include <DGtal/helpers/ShortcutsGeometry.h>
31 #include <DGtal/shapes/SurfaceMesh.h>
32 #include <DGtal/geometry/surfaces/DigitalSurfaceRegularization.h>
33 #include <DGtal/dec/PolygonalCalculus.h>
34 #include <DGtal/math/linalg/DirichletConditions.h>
35 
36 #include <polyscope/polyscope.h>
37 #include <polyscope/surface_mesh.h>
38 #include <polyscope/point_cloud.h>
39 
40 
41 #include <Eigen/Dense>
42 #include <Eigen/Sparse>
43 
44 using namespace DGtal;
45 using namespace Z3i;
46 
47 // Using standard 3D digital space.
50 // The following typedefs are useful
52 typedef std::size_t Index;
53 //Polyscope global
54 polyscope::SurfaceMesh *psMesh;
56 float scale = 0.1;
57 
59 {
63  PolyDEC calculus(surfmesh);
65  PolyDEC::Form g = calculus.form0();
66  DC::IntegerVector b = DC::IntegerVector::Zero( g.rows() );
67 
68  //We set values on the boundary
69  auto boundaryEdges = surfmesh.computeManifoldBoundaryEdges();
70  std::cout<< "Number of boundary edges= "<<boundaryEdges.size()<<std::endl;
71 
72  auto pihVertex=[&](const SurfMesh::Vertex &v){return cos(scale*(surfmesh.position(v)[0]))*(scale*surfmesh.position(v)[1]);};
73 
74  for(auto &e: boundaryEdges)
75  {
76  auto adjVertices = surfmesh.edgeVertices(e);
77  g(adjVertices.first) = pihVertex(adjVertices.first);
78  g(adjVertices.second) = pihVertex(adjVertices.second);
79  b(adjVertices.first) = 1;
80  b(adjVertices.second) = 1;
81  }
82 
83  // Solve Δu=0 with g as boundary conditions
84  PolyDEC::Solver solver;
85  PolyDEC::SparseMatrix L_dirichlet = DC::dirichletOperator( L, b );
86  solver.compute( L_dirichlet );
87  ASSERT(solver.info()==Eigen::Success);
88  PolyDEC::Form g_dirichlet = DC::dirichletVector( L, g, b, g );
89  PolyDEC::Form x_dirichlet = solver.solve( g_dirichlet );
90  PolyDEC::Form u = DC::dirichletSolution( x_dirichlet, b, g );
92 
93  psMesh->addVertexScalarQuantity("g", g);
94  psMesh->addVertexScalarQuantity("u", u);
95 }
96 
97 void myCallback()
98 {
99  ImGui::SliderFloat("Phi scale", &scale, 0., 10.);
100  if(ImGui::Button("Compute Laplace problem"))
101  computeLaplace();
102 }
103 
104 int main()
105 {
106  auto params = SH3::defaultParameters() | SHG3::defaultParameters() | SHG3::parametersGeometryEstimation();
107 
108  auto h=.7 ; //gridstep
109  params( "polynomial", "0.1*y*y -0.1*x*x - 2.0*z" )( "gridstep", h );
110  auto implicit_shape = SH3::makeImplicitShape3D ( params );
111  auto digitized_shape = SH3::makeDigitizedImplicitShape3D( implicit_shape, params );
112  auto K = SH3::getKSpace( params );
113  auto binary_image = SH3::makeBinaryImage( digitized_shape, params );
114  auto surface = SH3::makeDigitalSurface( binary_image, K, params );
115  auto primalSurface = SH3::makePrimalSurfaceMesh(surface);
116 
117  //Need to convert the faces
118  std::vector<std::vector<SH3::SurfaceMesh::Vertex>> faces;
119  std::vector<RealPoint> positions = primalSurface->positions();
120  for(auto face= 0 ; face < primalSurface->nbFaces(); ++face)
121  faces.push_back(primalSurface->incidentVertices( face ));
122 
123  surfmesh = SurfMesh(positions.begin(),
124  positions.end(),
125  faces.begin(),
126  faces.end());
127  psMesh = polyscope::registerSurfaceMesh("digital surface", positions, faces);
128 
129  // Initialize polyscope
130  polyscope::init();
131 
132  // Set the callback function
133  polyscope::state::userCallback = myCallback;
134  polyscope::show();
135  return EXIT_SUCCESS;
136 }
Aim: A helper class to solve a system with Dirichlet boundary conditions.
Implements differential operators on polygonal surfaces from .
SparseMatrix globalLaplaceBeltrami(const double lambda=1.0) const
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
PolyCalculus * calculus
CountedPtr< SH3::DigitalSurface > surface
CountedPtr< SH3::BinaryImage > binary_image
SurfaceMesh< RealPoint, RealVector > SurfMesh
std::size_t Index
void computeLaplace()
void myCallback()
polyscope::SurfaceMesh * psMesh
Shortcuts< Z3i::KSpace > SH3
SurfMesh surfmesh
ShortcutsGeometry< Z3i::KSpace > SHG3
int main()
DGtal is the top-level namespace which contains all DGtal functions and types.
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
Edges computeManifoldBoundaryEdges() const
VertexPair edgeVertices(Edge e) const
Definition: SurfaceMesh.h:338
K init(Point(0, 0, 0), Point(512, 512, 512), true)
KSpace K
EigenLinearAlgebraBackend::SparseMatrix SparseMatrix