DGtal 1.4.0
Loading...
Searching...
No Matches
dgtalCalculus-poisson.cpp File Reference
#include <iostream>
#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 dependency graph for dgtalCalculus-poisson.cpp:

Go to the source code of this file.

Typedefs

typedef Shortcuts< Z3i::KSpaceSH3
 
typedef ShortcutsGeometry< Z3i::KSpaceSHG3
 
typedef SurfaceMesh< RealPoint, RealVectorSurfMesh
 
typedef std::size_t Index
 

Functions

void computeLaplace ()
 
void myCallback ()
 
int main ()
 

Variables

polyscope::SurfaceMesh * psMesh
 
SurfMesh surfmesh
 
float scale = 0.1
 

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-poisson.cpp.

Typedef Documentation

◆ Index

typedef std::size_t Index

Definition at line 52 of file dgtalCalculus-poisson.cpp.

◆ SH3

Definition at line 48 of file dgtalCalculus-poisson.cpp.

◆ SHG3

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

◆ SurfMesh

Definition at line 51 of file dgtalCalculus-poisson.cpp.

Function Documentation

◆ computeLaplace()

void computeLaplace ( )

[PolyDEC-init]

[PolyDEC-init]

Definition at line 58 of file dgtalCalculus-poisson.cpp.

59{
63 PolyDEC calculus(surfmesh);
64 PolyDEC::SparseMatrix L = calculus.globalLaplaceBeltrami();
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}
Aim: A helper class to solve a system with Dirichlet boundary conditions.
Implements differential operators on polygonal surfaces from degoes2020discrete.
SparseMatrix globalLaplaceBeltrami(const double lambda=1.0) const
PolyCalculus * calculus
polyscope::SurfaceMesh * psMesh
SurfMesh surfmesh
Edges computeManifoldBoundaryEdges() const
VertexPair edgeVertices(Edge e) const
RealPoint & position(Vertex v)

References calculus, DGtal::SurfaceMesh< TRealPoint, TRealVector >::computeManifoldBoundaryEdges(), DGtal::SurfaceMesh< TRealPoint, TRealVector >::edgeVertices(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::form0(), DGtal::PolygonalCalculus< TRealPoint, TRealVector >::globalLaplaceBeltrami(), DGtal::SurfaceMesh< TRealPoint, TRealVector >::position(), psMesh, and surfmesh.

Referenced by myCallback().

◆ main()

int main ( void )

Definition at line 104 of file dgtalCalculus-poisson.cpp.

105{
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 );
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}
static Parameters parametersGeometryEstimation()
static Parameters defaultParameters()
static KSpace getKSpace(const Point &low, const Point &up, Parameters params=parametersKSpace())
Definition Shortcuts.h:332
static CountedPtr< DigitizedImplicitShape3D > makeDigitizedImplicitShape3D(CountedPtr< ImplicitShape3D > shape, Parameters params=parametersDigitizedImplicitShape3D())
Definition Shortcuts.h:523
static CountedPtr< DigitalSurface > makeDigitalSurface(CountedPtr< TPointPredicate > bimage, const KSpace &K, const Parameters &params=parametersDigitalSurface())
Definition Shortcuts.h:1209
static Parameters defaultParameters()
Definition Shortcuts.h:203
static CountedPtr< SurfaceMesh > makePrimalSurfaceMesh(Cell2Index &c2i, CountedPtr< ::DGtal::DigitalSurface< TContainer > > aSurface)
Definition Shortcuts.h:2372
static CountedPtr< BinaryImage > makeBinaryImage(Domain shapeDomain)
Definition Shortcuts.h:561
static CountedPtr< ImplicitShape3D > makeImplicitShape3D(const Parameters &params=parametersImplicitShape3D())
Definition Shortcuts.h:282
CountedPtr< SH3::DigitalSurface > surface
CountedPtr< SH3::BinaryImage > binary_image
SurfaceMesh< RealPoint, RealVector > SurfMesh
void myCallback()
KSpace K

References binary_image, DGtal::Shortcuts< TKSpace >::defaultParameters(), DGtal::ShortcutsGeometry< TKSpace >::defaultParameters(), DGtal::Shortcuts< TKSpace >::getKSpace(), K, DGtal::Shortcuts< TKSpace >::makeBinaryImage(), DGtal::Shortcuts< TKSpace >::makeDigitalSurface(), DGtal::Shortcuts< TKSpace >::makeDigitizedImplicitShape3D(), DGtal::Shortcuts< TKSpace >::makeImplicitShape3D(), DGtal::Shortcuts< TKSpace >::makePrimalSurfaceMesh(), myCallback(), DGtal::ShortcutsGeometry< TKSpace >::parametersGeometryEstimation(), psMesh, surface, and surfmesh.

◆ myCallback()

void myCallback ( )

Definition at line 97 of file dgtalCalculus-poisson.cpp.

98{
99 ImGui::SliderFloat("Phi scale", &scale, 0., 10.);
100 if(ImGui::Button("Compute Laplace problem"))
102}
void computeLaplace()

References computeLaplace().

Referenced by main().

Variable Documentation

◆ psMesh

polyscope::SurfaceMesh* psMesh

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

Referenced by computeLaplace(), and main().

◆ scale

float scale = 0.1

Definition at line 56 of file dgtalCalculus-poisson.cpp.

◆ surfmesh

SurfMesh surfmesh

Definition at line 55 of file dgtalCalculus-poisson.cpp.

Referenced by computeLaplace(), and main().