LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Public Member Functions | Private Attributes | List of all members
projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR > Class Template Reference

Creates and solves the discretised Hodge Laplacian source problems. More...

#include </home/nico/bildung/SemVI/thesis/lehrfempp/projects/hldo_sphere/operators/hodge_laplacians_source_problems.h>

Public Member Functions

 HodgeLaplaciansSourceProblems ()
 Constructor initializes basic mesh (Octaeder with radius 1.0) More...
 
void Compute ()
 Computes the Galerkin LSE for all three source problems. More...
 
void Solve ()
 solves the linear systems build in Compute More...
 
void SetMesh (std::shared_ptr< const lf::mesh::Mesh > mesh_p)
 Sets the mesh and creates dof_handler. More...
 
void SetLoadFunctions (std::function< SCALAR(const Eigen::Matrix< double, 3, 1 > &)> &f0, std::function< Eigen::Matrix< SCALAR, 3, 1 >(const Eigen::Matrix< double, 3, 1 > &)> &f1, std::function< SCALAR(const Eigen::Matrix< double, 3, 1 > &)> &f2)
 Sets the load functions. More...
 
void SetK (double k)
 Sets k for the souce problems. More...
 
Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 > GetLoadVector (int index)
 returns the Loadvector More...
 
lf::assemble::COOMatrix< SCALAR > GetGalerkinMatrix (int index)
 returns the Galerkin Matrix More...
 
Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 > GetMuZero ()
 retunrs the basis expansion coefficiants for the solution of the zero form More...
 
std::tuple< Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 >, Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 > > GetMuOne ()
 retunrs the basis expansion coefficiants for the solution of the one form More...
 
std::tuple< Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 >, Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 > > GetMuTwo ()
 retunrs the basis expansion coefficiants for the solution of the two form More...
 

Private Attributes

std::shared_ptr< const lf::mesh::Meshmesh_p_
 
std::function< SCALAR(const Eigen::Matrix< double, 3, 1 > &)> f0_
 
std::function< Eigen::Matrix< SCALAR, 3, 1 >(const Eigen::Matrix< double, 3, 1 > &)> f1_
 
std::function< SCALAR(const Eigen::Matrix< double, 3, 1 > &)> f2_
 
std::vector< lf::assemble::COOMatrix< SCALAR > > coo_matrix_
 
std::vector< Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 > > phi_
 
std::vector< Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 > > mu_
 
double k_
 

Detailed Description

template<typename SCALAR>
class projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >

Creates and solves the discretised Hodge Laplacian source problems.

Details regarding the mathematical derivations can be found in the thesis Hodge-Laplacians and Dirac Operators on the Surface of the 3-Sphere section 4.4.

Note
Only triangular meshes are supported

Definition at line 34 of file hodge_laplacians_source_problems.h.

Constructor & Destructor Documentation

◆ HodgeLaplaciansSourceProblems()

template<typename SCALAR >
projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::HodgeLaplaciansSourceProblems ( )
inline

Member Function Documentation

◆ Compute()

template<typename SCALAR >
void projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::Compute ( )
inline

Computes the Galerkin LSE for all three source problems.

\[ - \Delta_l u + k^2 u = f_l \qquad l = 0, 1, 2 \]

The Galerkin matrix will be accessable with GetGalerkinMatrix(int index) The load vector will be accessable with GetLoadVector(int index)

Definition at line 95 of file hodge_laplacians_source_problems.h.

References lf::assemble::COOMatrix< SCALAR >::AddToEntry(), lf::assemble::COOMatrix< SCALAR >::cols(), projects::hldo_sphere::operators::WhitneyOneHodgeLaplace< SCALAR >::Compute(), projects::hldo_sphere::operators::WhitneyTwoHodgeLaplace< SCALAR >::Compute(), projects::hldo_sphere::operators::WhitneyZeroHodgeLaplace< SCALAR >::Compute(), projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::coo_matrix_, projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::f0_, projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::f1_, projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::f2_, projects::hldo_sphere::operators::WhitneyOneHodgeLaplace< SCALAR >::GetGalerkinMatrix(), projects::hldo_sphere::operators::WhitneyTwoHodgeLaplace< SCALAR >::GetGalerkinMatrix(), projects::hldo_sphere::operators::WhitneyZeroHodgeLaplace< SCALAR >::GetGalerkinMatrix(), projects::hldo_sphere::operators::WhitneyOneHodgeLaplace< SCALAR >::GetLoadVector(), projects::hldo_sphere::operators::WhitneyTwoHodgeLaplace< SCALAR >::GetLoadVector(), projects::hldo_sphere::operators::WhitneyZeroHodgeLaplace< SCALAR >::GetLoadVector(), projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::k_, lf::base::RefEl::kPoint(), lf::base::RefEl::kSegment(), lf::base::RefEl::kTria(), projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::mesh_p_, lf::assemble::DofHandler::NumDofs(), projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::phi_, lf::assemble::COOMatrix< SCALAR >::rows(), projects::hldo_sphere::operators::WhitneyOneHodgeLaplace< SCALAR >::SetLoadFunction(), projects::hldo_sphere::operators::WhitneyTwoHodgeLaplace< SCALAR >::SetLoadFunction(), projects::hldo_sphere::operators::WhitneyZeroHodgeLaplace< SCALAR >::SetLoadFunction(), projects::hldo_sphere::operators::WhitneyOneHodgeLaplace< SCALAR >::SetMesh(), projects::hldo_sphere::operators::WhitneyTwoHodgeLaplace< SCALAR >::SetMesh(), projects::hldo_sphere::operators::WhitneyZeroHodgeLaplace< SCALAR >::SetMesh(), lf::assemble::COOMatrix< SCALAR >::setZero(), and lf::assemble::COOMatrix< SCALAR >::triplets().

◆ GetGalerkinMatrix()

template<typename SCALAR >
lf::assemble::COOMatrix< SCALAR > projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::GetGalerkinMatrix ( int  index)
inline

returns the Galerkin Matrix

This is the Matrix of the LSE with index index

Note
The Galerkin matrix must be computed with Compute() before calling this funciton

Definition at line 357 of file hodge_laplacians_source_problems.h.

References projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::coo_matrix_.

◆ GetLoadVector()

template<typename SCALAR >
Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 > projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::GetLoadVector ( int  index)
inline

returns the Loadvector

This is the righthandside of LSE with index index

Note
The loadvector must be computed with Compute before calling this function

Definition at line 344 of file hodge_laplacians_source_problems.h.

References projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::phi_.

◆ GetMuOne()

template<typename SCALAR >
std::tuple< Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 >, Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 > > projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::GetMuOne ( )
inline

retunrs the basis expansion coefficiants for the solution of the one form

The first vector of the return value are the coefficiants with respect to whitney one form basis function. They represent an approximation of

\[ u \]

The second vector of the return value are the coefficiants with respect to to the barycentric basis functions. They represent the approximation for

\[ p := \text{div}_{\Gamma}(u) \]

on the mesh.

Returns
basis expansion coefficiants of one form
Note
the methods Compute() and Solve() must be called before the result is available

Definition at line 400 of file hodge_laplacians_source_problems.h.

References projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::mesh_p_, and projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::mu_.

◆ GetMuTwo()

template<typename SCALAR >
std::tuple< Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 >, Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 > > projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::GetMuTwo ( )
inline

retunrs the basis expansion coefficiants for the solution of the two form

The first vector of the return value are the coefficiants with respect to rotated (90 degree colckwise) whitney one form basis function. They represent an approximation of

\[ j = \textbf{grad}_{\Gamma}(u) \]

The second vector of the return value are the coefficiants with respect to to the piecewise constant basis functions. They represent the approximation for

\[ u \]

on the mesh.

Returns
basis expansion coefficiants of two form
Note
the methods Compute() and Solve() must be called before the result is available

Definition at line 433 of file hodge_laplacians_source_problems.h.

References projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::mesh_p_, and projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::mu_.

◆ GetMuZero()

template<typename SCALAR >
Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 > projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::GetMuZero ( )
inline

retunrs the basis expansion coefficiants for the solution of the zero form

The basis expansion coefficiants are with respect to the barycentric basis functions on the mesh.

Returns
basis expansion coefficiants of zero form
Note
the methods Compute() and Solve() must be called before the result is available

Definition at line 374 of file hodge_laplacians_source_problems.h.

References projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::mu_.

◆ SetK()

template<typename SCALAR >
void projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::SetK ( double  k)
inline

Sets k for the souce problems.

Parameters
kreal valued parameter scaling the mass matrix

Definition at line 333 of file hodge_laplacians_source_problems.h.

References projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::k_.

◆ SetLoadFunctions()

template<typename SCALAR >
void projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::SetLoadFunctions ( std::function< SCALAR(const Eigen::Matrix< double, 3, 1 > &)> &  f0,
std::function< Eigen::Matrix< SCALAR, 3, 1 >(const Eigen::Matrix< double, 3, 1 > &)> &  f1,
std::function< SCALAR(const Eigen::Matrix< double, 3, 1 > &)> &  f2 
)
inline

◆ SetMesh()

template<typename SCALAR >
void projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::SetMesh ( std::shared_ptr< const lf::mesh::Mesh mesh_p)
inline

Sets the mesh and creates dof_handler.

Parameters
mesh_ppointer to the mesh

requries all cells in the mesh are triangles requries mesh global dimension to be 3

Definition at line 296 of file hodge_laplacians_source_problems.h.

References lf::base::RefEl::kTria(), projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::mesh_p_, and lf::base::RefEl::RefEl().

◆ Solve()

template<typename SCALAR >
void projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::Solve ( )
inline

solves the linear systems build in Compute

Uses the Matrices stored in coo_matrix_ and the righthandside vectors stored in phi_ to compute the basis expansion coefficiants mu_ approximating

Note
the method Compute() has to be called prior to calling Solve()

\[ - \Delta_l u + k^2 u = f_l \qquad l = 0, 1, 2 \]

Definition at line 276 of file hodge_laplacians_source_problems.h.

References projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::coo_matrix_, projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::mu_, and projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::phi_.

Member Data Documentation

◆ coo_matrix_

template<typename SCALAR >
std::vector<lf::assemble::COOMatrix<SCALAR> > projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::coo_matrix_
private

◆ f0_

template<typename SCALAR >
std::function<SCALAR(const Eigen::Matrix<double, 3, 1> &)> projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::f0_
private

◆ f1_

template<typename SCALAR >
std::function<Eigen::Matrix<SCALAR, 3, 1>( const Eigen::Matrix<double, 3, 1> &)> projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::f1_
private

◆ f2_

template<typename SCALAR >
std::function<SCALAR(const Eigen::Matrix<double, 3, 1> &)> projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::f2_
private

◆ k_

template<typename SCALAR >
double projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::k_
private

◆ mesh_p_

template<typename SCALAR >
std::shared_ptr<const lf::mesh::Mesh> projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::mesh_p_
private

◆ mu_

template<typename SCALAR >
std::vector<Eigen::Matrix<SCALAR, Eigen::Dynamic, 1> > projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::mu_
private

◆ phi_

template<typename SCALAR >
std::vector<Eigen::Matrix<SCALAR, Eigen::Dynamic, 1> > projects::hldo_sphere::operators::HodgeLaplaciansSourceProblems< SCALAR >::phi_
private

The documentation for this class was generated from the following file: