v0.13.2
Loading...
Searching...
No Matches
Classes | Public Member Functions | Public Attributes | List of all members
AnalyticalDirichletBC Struct Reference

Analytical Dirichlet boundary conditions. More...

#include <users_modules/basic_finite_elements/src/AnalyticalDirichlet.hpp>

Collaboration diagram for AnalyticalDirichletBC:
[legend]

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
 

Detailed Description

Analytical Dirichlet boundary conditions.

Definition at line 18 of file AnalyticalDirichlet.hpp.

Constructor & Destructor Documentation

◆ AnalyticalDirichletBC()

AnalyticalDirichletBC::AnalyticalDirichletBC ( MoFEM::Interface m_field)

Definition at line 188 of file AnalyticalDirichlet.cpp.

189 : approxField(m_field){};

Member Function Documentation

◆ destroyProblem()

MoFEMErrorCode AnalyticalDirichletBC::destroyProblem ( )

Destroy problem.

Destroy matrices and vectors used to solve boundary problem, i.e. finding values of DOFs which approximate analytical solution.

Returns
error code

Definition at line 268 of file AnalyticalDirichlet.cpp.

268 {
270 CHKERR KSPDestroy(&kspSolver);
271 CHKERR MatDestroy(&A);
272 CHKERR VecDestroy(&F);
273 CHKERR VecDestroy(&D);
275}
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535

◆ setApproxOps()

template<typename FUNEVAL >
MoFEMErrorCode AnalyticalDirichletBC::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" 
)
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.

Parameters
m_fieldinterface
field_namefield name
function_evaluatoranalytical function to evaluate
field_numberfield index
nodals_positionsname of the field for ho-geometry description
Returns
error code

Definition at line 158 of file AnalyticalDirichlet.hpp.

161 {
164 if (m_field.check_field(nodals_positions))
166 false, false);
168 new ApproxField::OpLhs(field_name));
169 }
171 new ApproxField::OpRhs<FUNEVAL>(field_name, function_evaluator,
172 field_number));
174 }
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
virtual bool check_field(const std::string &name) const =0
check if field is in database
MoFEMErrorCode addHOOpsFace3D(const std::string field, E &e, bool hcurl, bool hdiv)
constexpr auto field_name
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.

◆ setFiniteElement()

MoFEMErrorCode AnalyticalDirichletBC::setFiniteElement ( MoFEM::Interface m_field,
string  fe,
string  field,
Range tris,
string  nodals_positions = "MESH_NODE_POSITIONS" 
)

set finite element

Parameters
m_fieldmofem interface
fefinite element name
fieldfield name
trisfaces where analytical boundary is given
nodals_positionsfield having higher order geometry description
Returns
error code

Definition at line 192 of file AnalyticalDirichlet.cpp.

194 {
196
197 CHKERR m_field.add_finite_element(fe, MF_ZERO);
201 if (m_field.check_field(nodals_positions)) {
202 CHKERR m_field.modify_finite_element_add_field_data(fe, nodals_positions);
203 }
204 CHKERR m_field.add_ents_to_finite_element_by_type(tris, MBTRI, fe);
206}
@ MF_ZERO
Definition: definitions.h:98
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string name_row)=0
set field col which finite element use
virtual MoFEMErrorCode add_ents_to_finite_element_by_type(const EntityHandle entities, const EntityType type, const std::string &name, const bool recursive=true)=0
add entities to finite element
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string name_row)=0
set field row which finite element use

◆ setUpProblem()

MoFEMErrorCode AnalyticalDirichletBC::setUpProblem ( MoFEM::Interface m_field,
string  problem 
)

set problem solver and create matrices and vectors

Parameters
m_fieldmofem interface
problemproblem name
Returns
error code

Definition at line 208 of file AnalyticalDirichlet.cpp.

209 {
211 CHKERR m_field.getInterface<VecManager>()->vecCreateGhost(problem, ROW, &F);
212 CHKERR m_field.getInterface<VecManager>()->vecCreateGhost(problem, COL, &D);
214 ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>(problem, &A);
215
216 CHKERR KSPCreate(PETSC_COMM_WORLD, &kspSolver);
217 CHKERR KSPSetOperators(kspSolver, A, A);
218 CHKERR KSPSetFromOptions(kspSolver);
220}
@ COL
Definition: definitions.h:123
@ ROW
Definition: definitions.h:123
Matrix manager is used to build and partition problems.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Vector manager is used to create vectors \mofem_vectors.
Definition: VecManager.hpp:23

◆ solveProblem() [1/2]

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.

Parameters
m_fieldmofem interface
problemproblem name
fefinite element name
bcDriblet boundary structure used to apply boundary conditions
Returns
[description]

Definition at line 257 of file AnalyticalDirichlet.cpp.

259 {
261 EntityHandle fe_meshset = m_field.get_finite_element_meshset("BC_FE");
262 Range bc_tris;
263 CHKERR m_field.get_moab().get_entities_by_type(fe_meshset, MBTRI, bc_tris);
264 return solveProblem(m_field, problem, fe, bc, bc_tris);
266}
virtual EntityHandle get_finite_element_meshset(const std::string name) const =0
MoFEMErrorCode solveProblem(MoFEM::Interface &m_field, string problem, string fe, DirichletBC &bc, Range &tris)
solve boundary problem
virtual moab::Interface & get_moab()=0

◆ solveProblem() [2/2]

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.

Parameters
m_fieldmofem interface
problemproblem name
fefinite element name
bcDriblet boundary structure used to apply boundary conditions
tristriangles on boundary
Returns
error code

Definition at line 222 of file AnalyticalDirichlet.cpp.

225 {
227
228 CHKERR VecZeroEntries(F);
229 CHKERR MatZeroEntries(A);
230
233 CHKERR m_field.loop_finite_elements(problem, fe,
235 CHKERR VecGhostUpdateBegin(F, ADD_VALUES, SCATTER_REVERSE);
236 CHKERR VecGhostUpdateEnd(F, ADD_VALUES, SCATTER_REVERSE);
237 CHKERR VecAssemblyBegin(F);
238 CHKERR VecAssemblyEnd(F);
239 CHKERR MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
240 CHKERR MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
241
242 CHKERR KSPSolve(kspSolver, F, D);
243 CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
244 CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
245
246 CHKERR m_field.getInterface<VecManager>()->setGlobalGhostVector(
247 problem, ROW, D, INSERT_VALUES, SCATTER_REVERSE);
248
249 bc.trisPtr = boost::shared_ptr<Range>(new Range(tris));
250 bc.mapZeroRows.clear();
251 bc.dofsIndices.clear();
252 bc.dofsValues.clear();
253
255}
virtual MoFEMErrorCode loop_finite_elements(const std::string problem_name, const std::string &fe_name, FEMethod &method, boost::shared_ptr< NumeredEntFiniteElement_multiIndex > fe_ptr=nullptr, MoFEMTypes bh=MF_EXIST, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr(), int verb=DEFAULT_VERBOSITY)=0
Make a loop over finite elements.
Vec & snes_f
residual
Mat & snes_B
preconditioner of jacobian matrix

Member Data Documentation

◆ A

Mat AnalyticalDirichletBC::A

Definition at line 219 of file AnalyticalDirichlet.hpp.

◆ approxField

ApproxField AnalyticalDirichletBC::approxField

Definition at line 137 of file AnalyticalDirichlet.hpp.

◆ D

Vec AnalyticalDirichletBC::D

Definition at line 220 of file AnalyticalDirichlet.hpp.

◆ F

Vec AnalyticalDirichletBC::F

Definition at line 220 of file AnalyticalDirichlet.hpp.

◆ kspSolver

KSP AnalyticalDirichletBC::kspSolver

Definition at line 221 of file AnalyticalDirichlet.hpp.


The documentation for this struct was generated from the following files: