31 #if defined(VectorsInHeat_RECURSES)
32 #error Recursive header files inclusion detected in VectorsInHeat.h
35 #define VectorsInHeat_RECURSES
37 #if !defined VectorsInHeat_h
39 #define VectorsInHeat_h
44 #include "DGtal/base/Common.h"
45 #include "DGtal/base/ConstAlias.h"
46 #include "DGtal/math/linalg/DirichletConditions.h"
61 template <
typename TPolygonalCalculus>
134 void init(
double dt,
double lambda = 1.0,
135 bool boundary_with_mixed_solution =
false )
159 if ( ! boundary_with_mixed_solution )
return;
162 const auto edges =
surfmesh->computeManifoldBoundaryEdges();
163 for (
auto e : edges )
165 const auto vtcs =
surfmesh->edgeVertices( e );
184 ASSERT_MSG(aV < myCalculus->nbVertices(),
"Vertex is not in the surface mesh vertex range");
186 v = v.normalized()*ev.norm();
207 FATAL_ERROR_MSG(
myIsInit,
"init() method must be called first");
217 FATAL_ERROR_MSG(
myIsInit,
"init() method must be called first");
227 FATAL_ERROR_MSG(
myIsInit,
"init() method must be called first");
239 FATAL_ERROR_MSG(
myIsInit,
"init() method must be called first");
254 Vector heatDiffusionDirichlet
256 scalarHeatDiffusion = 0.5 * ( scalarHeatDiffusion + heatDiffusionDirichlet );
259 std::vector<Vector> result(
surfmesh->nbVertices());
263 Y(0) = vectorHeatDiffusion(2*v);
264 Y(1) = vectorHeatDiffusion(2*v+1);
265 Y = Y.normalized()*(scalarHeatDiffusion(v)/diracHeatDiffusion(v));
266 result[v] =
myCalculus->toExtrinsicVector(v,Y);
320 #undef VectorsInHeat_RECURSES
Aim: This class encapsulates its parameter class so that to indicate to the user that the object/poin...
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)
LinAlg::SparseMatrix SparseMatrix
Type of sparse matrix.
MySurfaceMesh::Vertex Vertex
Vertex type.
LinAlg::SolverSimplicialLDLT Solver
Type of a sparse matrix solver.
LinAlg::DenseMatrix DenseMatrix
Type of dense matrix.
LinAlg::DenseVector Vector
Type of Vector.
This class implements on polygonal surfaces (using Discrete differential calculus on polygonal surfa...
Vector extrinsicVectorSourceAtVertex(const Vertex aV)
extrinsicVectorSourceAtVertex get extrinsic source at vertex
Vector intrinsicVectorSourceAtVertex(const Vertex aV)
intrinsicVectorSourceAtVertex get intrinsic source at vertex
PolygonalCalculus::Vector Vector
DirichletConditions< LinAlgBackend > Conditions
void addSource(const Vertex aV, const Vector &ev)
PolygonalCalculus::LinAlg LinAlgBackend
Vector myScalarSource
Source vectors.
PolygonalCalculus::Vertex Vertex
PolygonalCalculus::Solver Solver
PolygonalCalculus::DenseMatrix DenseMatrix
void init(double dt, double lambda=1.0, bool boundary_with_mixed_solution=false)
VectorsInHeat & operator=(const VectorsInHeat &other)=delete
SparseMatrix myScalarHeatOpe
The operators for heat diffusion.
std::vector< Vector > compute() const
PolygonalCalculus::SparseMatrix SparseMatrix
Vector vectorSource() const
VectorsInHeat(VectorsInHeat &&other)=delete
VectorsInHeat(const VectorsInHeat &other)=delete
Solver myVectorHeatSolver
IntegerVector myBoundary
The boundary characteristic vector.
bool myIsInit
Validitate flag.
Solver myScalarHeatSolver
Heat solvers.
TPolygonalCalculus PolygonalCalculus
Conditions::IntegerVector IntegerVector
VectorsInHeat(ConstAlias< PolygonalCalculus > calculus)
SparseMatrix myVectorHeatOpe
const PolygonalCalculus * myCalculus
The underlying PolygonalCalculus instance.
Solver myHeatDirichletSolver
Heat solver with Dirichlet boundary conditions.
DGtal is the top-level namespace which contains all DGtal functions and types.
void laplacian(Shape &shape, const Options &options, std::function< double(const RealPoint3D &)> input_function, std::function< double(const RealPoint3D &)> target_function, int argc, char **argv)
Aim: Provide linear algebra backend using Eigen dense and sparse matrix as well as dense vector....
std::size_t Index
The type used for numbering vertices and faces.