34#include "DGtal/base/Common.h"
35#include "ConfigTest.h"
36#include "DGtalCatch.h"
37#include "DGtal/helpers/StdDefs.h"
39#include "DGtal/dec/PolygonalCalculus.h"
40#include "DGtal/shapes/SurfaceMesh.h"
41#include "DGtal/shapes/MeshHelpers.h"
42#include "DGtal/helpers/Shortcuts.h"
43#include "DGtal/math/linalg/DirichletConditions.h"
62 std::vector<RealPoint> positions = {
RealPoint( 0, 0, 0 ) ,
72 std::vector<Mesh::Vertices> faces = { { 1, 0, 2, 3 },
79 Mesh box(positions.cbegin(), positions.cend(),
80 faces.cbegin(), faces.cend());
84 SECTION(
"Construction and basic operators")
90 auto x = boxCalculus.
X(f);
91 auto d = boxCalculus.
D(f);
92 auto a = boxCalculus.
A(f);
118 for(
auto ff=0; ff<6; ++ff)
121 box.computeFaceNormalsFromPositions();
122 for(
auto ff=0; ff<6; ++ff)
125 auto n = box.faceNormal(ff);
137 auto d = boxCalculus.
D(f);
141 phi << 1.0, 3.0, 2.0, 6.0;
142 expected << 2,-1,4,-5;
148 SECTION(
"Structural propertes")
153 phi << 1.0, 3.0, 2.0, 6.0;
158 auto cogphi = coG*
phi;
161 REQUIRE( gphi.dot(cogphi) == 0.0);
171 auto P = boxCalculus.
P(f);
181 auto curl = boxCalculus.
curl(f);
186 SECTION(
"Local Laplace-Beltrami")
189 auto nf = box.incidentVertices(f).size();
193 phi << 1.0, 1.0, 1.0, 1.0;
199 SECTION(
"Local Connection-Laplace-Beltrami")
202 auto nf = box.incidentVertices(f).size();
206 phi << 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0;
211 double det = L.determinant()+1.;
213 REQUIRE( lphi[2] == Approx(-3.683));
219 auto nf = box.incidentVertices(f).size();
222 phi << 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0;
229 REQUIRE( CP.rows() == (Eigen::Index)faces[f].size());
232 REQUIRE( CG(0,0) == Approx(0.707106));
233 REQUIRE( CP(0,0) == Approx(1.224744));
235 SECTION(
"Check lumped mass matrix")
240 a += M.coeffRef(v,v);
244 fa += box.faceArea(f);
260 trace.
info()<< boxCalculusCached <<std::endl;
264 for(
auto i=0; i < 1000 ; ++i)
270 for(
auto i=0; i < 1000 ; ++i)
274 REQUIRE(L.norm() == Approx(LC.norm()));
280TEST_CASE(
"Testing PolygonalCalculus and DirichletConditions" )
291 params(
"polynomial",
"0.1*y*y -0.1*x*x - 2.0*z" )(
"gridstep", 2.0 );
299 std::vector<std::vector< Index > > faces;
300 std::vector<RealPoint> positions = primalSurface->positions();
302 faces.push_back(primalSurface->incidentVertices( face ));
305 faces.begin(), faces.end() );
315 DC::IntegerVector b = DC::IntegerVector::Zero( g.rows() );
317 SECTION(
"Solve Poisson problem with boundary Dirichlet conditions")
319 for (
double scale = 0.1; scale < 2.0; scale *= 2.0 )
321 std::cout <<
"scale=" << scale << std::endl;
328 for(
auto &e: boundaryEdges)
331 auto v1 = adjVertices.first;
332 auto v2 = adjVertices.second;
339 PolyDEC::Solver solver;
340 PolyDEC::SparseMatrix L_dirichlet = DC::dirichletOperator( L, b );
341 solver.compute( L_dirichlet );
342 REQUIRE( solver.info() == Eigen::Success );
343 PolyDEC::Vector g_dirichlet = DC::dirichletVector( L, g, b, g );
344 PolyDEC::Vector x_dirichlet = solver.solve( g_dirichlet );
345 REQUIRE( solver.info() == Eigen::Success );
346 PolyDEC::Vector u = DC::dirichletSolution( x_dirichlet, b, g );
347 double min_phi = 0.0;
348 double max_phi = 0.0;
351 double min_i_u = 0.0;
352 double max_i_u = 0.0;
355 min_phi = std::min( min_phi,
phi( v ) );
356 max_phi = std::max( max_phi,
phi( v ) );
357 min_u = std::min( min_u , u ( v ) );
358 max_u = std::max( max_u , u ( v ) );
361 min_i_u = std::min( min_i_u, u ( v ) );
362 max_i_u = std::max( max_i_u, u ( v ) );
Aim: A helper class to solve a system with Dirichlet boundary conditions.
Aim: This class is defined to represent a surface mesh through a set of vertices and faces....
VertexStorage::size_type Index
Aim: Implements basic operations that will be used in Point and Vector classes.
Implements differential operators on polygonal surfaces from degoes2020discrete.
DenseMatrix D(const Face f) const
DenseMatrix laplaceBeltrami(const Face f, const double lambda=1.0) const
LinAlg::SparseMatrix SparseMatrix
Type of sparse matrix.
DenseMatrix X(const Face f) const
size_t nbVertices() const
std::vector< DenseMatrix > getOperatorCacheMatrix(const std::function< DenseMatrix(Face)> &perFaceOperator) const
Vector faceNormal(const Face f) const
Vector centroid(const Face f) const
DenseMatrix curl(const Face f) const
Real3dVector faceNormalAsDGtalVector(const Face f) const
Real3dPoint centroidAsDGtalPoint(const Face f) const
DenseMatrix covariantProjection(const Face f, const Vector &uf)
Vector vectorArea(const Face f) const
size_t faceDegree(Face f) const
DenseMatrix P(const Face f) const
SparseMatrix globalLaplaceBeltrami(const double lambda=1.0) const
DenseMatrix sharp(const Face f) const
double faceArea(const Face f) const
std::vector< Vector > getOperatorCacheVector(const std::function< Vector(Face)> &perFaceVectorOperator) const
DenseMatrix coGradient(const Face f) const
DenseMatrix connectionLaplacian(const Face &f, double lambda=1.0) const
SparseMatrix globalLumpedMassMatrix() const
DenseMatrix gradient(const Face f) const
MySurfaceMesh::Face Face
Face type.
LinAlg::DenseVector Vector
Type of Vector.
DenseMatrix A(const Face f) const
DenseMatrix covariantGradient(const Face f, const Vector &uf)
DenseMatrix flat(const Face f) const
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< DigitizedImplicitShape3D > makeDigitizedImplicitShape3D(CountedPtr< ImplicitShape3D > shape, Parameters params=parametersDigitizedImplicitShape3D())
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)
static CountedPtr< ImplicitShape3D > makeImplicitShape3D(const Parameters ¶ms=parametersImplicitShape3D())
void beginBlock(const std::string &keyword="")
PolygonalCalculus< SH3::RealPoint, SH3::RealVector >::Vector phi(const Face f)
CountedPtr< SH3::DigitalSurface > surface
CountedPtr< SH3::BinaryImage > binary_image
Space::RealPoint RealPoint
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::size_t Index
The type used for numbering vertices and faces.
Edges computeManifoldBoundaryEdges() const
VertexPair edgeVertices(Edge e) const
RealPoint & position(Vertex v)
TEST_CASE("Testing PolygonalCalculus")
RealPoint vecToRealPoint(PolygonalCalculus< RealPoint, RealPoint >::Vector &v)
SECTION("Testing constant forward iterators")
REQUIRE(domain.isInside(aPoint))