v0.14.0 |
Analytical Dirichlet boundary conditions. More...
#include <users_modules/basic_finite_elements/src/AnalyticalDirichlet.hpp>
Classes | |
struct | ApproxField |
finite element to approximate analytical solution on surface More... | |
struct | DirichletBC |
Structure used to enforce analytical boundary conditions. More... | |
Public Member Functions | |
AnalyticalDirichletBC (MoFEM::Interface &m_field) | |
template<typename FUNEVAL > | |
MoFEMErrorCode | setApproxOps (MoFEM::Interface &m_field, const std::string field_name, boost::shared_ptr< FUNEVAL > function_evaluator, const int field_number=0, const string nodals_positions="MESH_NODE_POSITIONS") |
Set operators used to calculate the rhs vector and the lhs matrix. More... | |
MoFEMErrorCode | setFiniteElement (MoFEM::Interface &m_field, string fe, string field, Range &tris, string nodals_positions="MESH_NODE_POSITIONS") |
set finite element More... | |
MoFEMErrorCode | setUpProblem (MoFEM::Interface &m_field, string problem) |
set problem solver and create matrices and vectors More... | |
MoFEMErrorCode | solveProblem (MoFEM::Interface &m_field, string problem, string fe, DirichletBC &bc, Range &tris) |
solve boundary problem More... | |
MoFEMErrorCode | solveProblem (MoFEM::Interface &m_field, string problem, string fe, DirichletBC &bc) |
solve boundary problem More... | |
MoFEMErrorCode | destroyProblem () |
Destroy problem. More... | |
Public Attributes | |
ApproxField | approxField |
Mat | A |
Vec | D |
Vec | F |
KSP | kspSolver |
Analytical Dirichlet boundary conditions.
Definition at line 18 of file AnalyticalDirichlet.hpp.
AnalyticalDirichletBC::AnalyticalDirichletBC | ( | MoFEM::Interface & | m_field | ) |
Definition at line 188 of file AnalyticalDirichlet.cpp.
MoFEMErrorCode AnalyticalDirichletBC::destroyProblem | ( | ) |
Destroy problem.
Destroy matrices and vectors used to solve boundary problem, i.e. finding values of DOFs which approximate analytical solution.
Definition at line 268 of file AnalyticalDirichlet.cpp.
|
inline |
Set operators used to calculate the rhs vector and the lhs matrix.
To enforce analytical function on boundary, first function has to be approximated by finite element base functions. This is done by solving system of linear equations, Following function set finite element operators to calculate the left hand side matrix and the left hand side matrix.
m_field | interface |
field_name | field name |
function_evaluator | analytical function to evaluate |
field_number | field index |
nodals_positions | name of the field for ho-geometry description |
Definition at line 158 of file AnalyticalDirichlet.hpp.
MoFEMErrorCode AnalyticalDirichletBC::setFiniteElement | ( | MoFEM::Interface & | m_field, |
string | fe, | ||
string | field, | ||
Range & | tris, | ||
string | nodals_positions = "MESH_NODE_POSITIONS" |
||
) |
set finite element
m_field | mofem interface |
fe | finite element name |
field | field name |
tris | faces where analytical boundary is given |
nodals_positions | field having higher order geometry description |
Definition at line 192 of file AnalyticalDirichlet.cpp.
MoFEMErrorCode AnalyticalDirichletBC::setUpProblem | ( | MoFEM::Interface & | m_field, |
string | problem | ||
) |
set problem solver and create matrices and vectors
m_field | mofem interface |
problem | problem name |
Definition at line 208 of file AnalyticalDirichlet.cpp.
MoFEMErrorCode AnalyticalDirichletBC::solveProblem | ( | MoFEM::Interface & | m_field, |
string | problem, | ||
string | fe, | ||
DirichletBC & | bc | ||
) |
solve boundary problem
This functions solve for DOFs on boundary where analytical solution is given.
m_field | mofem interface |
problem | problem name |
fe | finite element name |
bc | Driblet boundary structure used to apply boundary conditions |
Definition at line 257 of file AnalyticalDirichlet.cpp.
MoFEMErrorCode AnalyticalDirichletBC::solveProblem | ( | MoFEM::Interface & | m_field, |
string | problem, | ||
string | fe, | ||
DirichletBC & | bc, | ||
Range & | tris | ||
) |
solve boundary problem
This functions solve for DOFs on boundary where analytical solution is given. i.e. finding values of DOFs which approximate analytical solution.
m_field | mofem interface |
problem | problem name |
fe | finite element name |
bc | Driblet boundary structure used to apply boundary conditions |
tris | triangles on boundary |
Definition at line 222 of file AnalyticalDirichlet.cpp.
Mat AnalyticalDirichletBC::A |
Definition at line 219 of file AnalyticalDirichlet.hpp.
ApproxField AnalyticalDirichletBC::approxField |
Definition at line 137 of file AnalyticalDirichlet.hpp.
Vec AnalyticalDirichletBC::D |
Definition at line 220 of file AnalyticalDirichlet.hpp.
Vec AnalyticalDirichletBC::F |
Definition at line 220 of file AnalyticalDirichlet.hpp.
KSP AnalyticalDirichletBC::kspSolver |
Definition at line 221 of file AnalyticalDirichlet.hpp.