30#include <DGtal/base/Common.h>
31#include <DGtal/helpers/StdDefs.h>
32#include <DGtal/helpers/Shortcuts.h>
33#include <DGtal/helpers/ShortcutsGeometry.h>
34#include <DGtal/shapes/SurfaceMesh.h>
35#include <DGtal/dec/PolygonalCalculus.h>
36#include "DGtal/math/linalg/DirichletConditions.h"
38#include <polyscope/polyscope.h>
39#include <polyscope/surface_mesh.h>
41#include "ConfigExamples.h"
63typedef std::vector<Vertex> chain;
69polyscope::SurfaceMesh *
psMesh;
70polyscope::SurfaceMesh *psParam;
77std::vector<std::vector<size_t>> faces;
78std::vector<RealPoint> positions;
87std::pair<Vector,Vector> FixBoundaryParametrization(
const std::vector<Vertex>& boundary)
89 auto nb = boundary.size();
91 Vector u = Vector::Zero(n),v = Vector::Zero(n);
92 double totalBoundaryLength = 0;
93 for (
Vertex i = 0;i<nb;i++)
96 double partialSum = 0;
97 for (
Vertex i = 0;i<nb;i++)
99 double th = 2*M_PI*partialSum/totalBoundaryLength;
100 auto vi = boundary[i];
101 auto vj = boundary[(i+1)%nb];
102 u(vi) = std::cos(th);
103 v(vj) = std::sin(th);
114void VisualizeParametrizationOnCircle(
const DenseMatrix& UV)
117 std::vector<RealPoint> pos(n);
120 for (
size_t v = 0;v<n;v++)
124 if (p.norm() > scale)
129 for (
size_t v = 0;v<n;v++)
130 pos[v] =
RealPoint{scale*UV(v,0),scale*UV(v,1),0.} + avg;
131 polyscope::registerSurfaceMesh(
"On circle parametrization", pos, faces)->setEnabled(
false);
146 std::cout<<
"Nb boundaries = "<< chains.size() << std::endl;
149 auto B = *std::max_element(chains.begin(),chains.end(),[] (
const chain& A,
const chain& B) {return A.size() < B.size();});
151 IntegerVector boundary = IntegerVector::Zero(n);
155 std::pair<Vector,Vector> uv_b = FixBoundaryParametrization(B);
160 Vector Z = Vector::Zero(n);
170 Vector rslt_u_d = solver.solve(b_u);
171 Vector rslt_v_d = solver.solve(b_v);
183int main(
int argc,
char **argv)
187 params(
"surfaceComponents",
"All");
190 std::string inputFilename(examplesPath +
"samples/bunny-64.vol" );
203 for(
size_t face= 0 ; face < primalSurface->nbFaces(); ++face)
204 faces.push_back(primalSurface->incidentVertices( face ));
207 positions = primalSurface->positions();
222 VisualizeParametrizationOnCircle(UV);
224 psMesh = polyscope::registerSurfaceMesh(
"Digital Surface", positions, faces);
225 psMesh->addVertexParameterizationQuantity(
"Harmonic parametrization",UV)->setEnabled(
true);
228 polyscope::view::upDir = polyscope::view::UpDir::XUp;
Aim: A helper class to solve a system with Dirichlet boundary conditions.
LinearAlgebraBackend::IntegerVector IntegerVector
static DenseVector dirichletVector(const SparseMatrix &A, const DenseVector &b, const IntegerVector &p, const DenseVector &u)
static SparseMatrix dirichletOperator(const SparseMatrix &A, const IntegerVector &p)
static DenseVector dirichletSolution(const DenseVector &xd, const IntegerVector &p, const DenseVector &u)
const Point & lowerBound() const
Return the lower bound for digital points in this space.
bool init(const Point &lower, const Point &upper, bool isClosed)
Specifies the upper and lower bounds for the maximal cells in this space.
const Point & upperBound() const
Return the upper bound for digital points in this space.
double norm(const NormType type=L_2) const
Implements differential operators on polygonal surfaces from degoes2020discrete.
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::DenseMatrix DenseMatrix
Type of dense matrix.
LinAlg::Triplet Triplet
Type of sparse matrix triplet.
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...
static Parameters parametersGeometryEstimation()
static Parameters defaultParameters()
Aim: This class is used to simplify shape and surface creation. With it, you can create new shapes an...
static KSpace getKSpace(const Point &low, const Point &up, Parameters params=parametersKSpace())
static CountedPtr< DigitalSurface > makeDigitalSurface(CountedPtr< TPointPredicate > bimage, const KSpace &K, const Parameters ¶ms=parametersDigitalSurface())
static Parameters defaultParameters()
static CountedPtr< SurfaceMesh > makePrimalSurfaceMesh(Cell2Index &c2i, CountedPtr< ::DGtal::DigitalSurface< TContainer > > aSurface)
static CountedPtr< BinaryImage > makeBinaryImage(Domain shapeDomain)
polyscope::SurfaceMesh * psMesh
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.
Aim: Represents an embedded mesh as faces and a list of vertices. Vertices may be shared among faces ...
std::vector< Vertex > Vertices
The type that defines a list/range of vertices (e.g. to define faces)
Scalar distance(const Vertex i, const Vertex j) const
Edges computeManifoldBoundaryEdges() const
RealPoint & position(Vertex v)
std::vector< Vertices > computeManifoldBoundaryChains() const
EigenLinearAlgebraBackend::SparseMatrix SparseMatrix
EigenLinearAlgebraBackend::DenseMatrix DenseMatrix