v0.15.0
Loading...
Searching...
No Matches
Classes
Loops

Manages complexities for integrating over finite elements and dofs. More...

Collaboration diagram for Loops:

Classes

struct  MoFEM::PetscData
 Base data structure for PETSc-related contexts. More...
 
struct  MoFEM::KspMethod
 Data structure for KSP (linear solver) context. More...
 
struct  MoFEM::SnesMethod
 Data structure for SNES (nonlinear solver) context. More...
 
struct  MoFEM::TSMethod
 Data structure for TS (time stepping) context. More...
 
struct  MoFEM::TaoMethod
 Data structure for TAO (optimization) context. More...
 
struct  MoFEM::BasicMethod
 Data structure to exchange data between MoFEM and User Loop Methods. More...
 
struct  MoFEM::FEMethod
 Structure for user loop methods on finite elements. More...
 
struct  MoFEM::EntityMethod
 Data structure for user loop methods on entities. More...
 
struct  MoFEM::DofMethod
 Data structure for user loop methods on degrees of freedom (DOFs) More...
 

Making loops on elements and entities

virtual MoFEMErrorCode MoFEM::CoreInterface::problem_basic_method_postProcess (const Problem *problem_ptr, BasicMethod &method, int verb=DEFAULT_VERBOSITY)=0
 Set data for BasicMethod.
 
virtual MoFEMErrorCode MoFEM::CoreInterface::problem_basic_method_postProcess (const std::string &problem_name, BasicMethod &method, int verb=DEFAULT_VERBOSITY)=0
 Set data for BasicMethod.
 
virtual MoFEMErrorCode MoFEM::CoreInterface::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.
 
virtual MoFEMErrorCode MoFEM::CoreInterface::loop_finite_elements (const Problem *problem_ptr, const std::string &fe_name, FEMethod &method, int lower_rank, int upper_rank, 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 on partitions from upper to lower rank.
 
virtual MoFEMErrorCode MoFEM::CoreInterface::loop_finite_elements (const std::string problem_name, const std::string &fe_name, FEMethod &method, int lower_rank, int upper_rank, 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 on partitions from upper to lower rank.
 
virtual MoFEMErrorCode MoFEM::CoreInterface::loop_dofs (const Problem *problem_ptr, const std::string &field_name, RowColData rc, DofMethod &method, int lower_rank, int upper_rank, int verb=DEFAULT_VERBOSITY)=0
 Make a loop over dofs.
 
virtual MoFEMErrorCode MoFEM::CoreInterface::loop_dofs (const std::string &problem_name, const std::string &field_name, RowColData rc, DofMethod &method, int lower_rank, int upper_rank, int verb=DEFAULT_VERBOSITY)=0
 Make a loop over dofs.
 
virtual MoFEMErrorCode MoFEM::CoreInterface::loop_dofs (const std::string &problem_name, const std::string &field_name, RowColData rc, DofMethod &method, int verb=DEFAULT_VERBOSITY)=0
 Make a loop over dofs.
 

Detailed Description

Manages complexities for integrating over finite elements and dofs.

Function Documentation

◆ loop_dofs() [1/3]

virtual MoFEMErrorCode MoFEM::CoreInterface::loop_dofs ( const Problem problem_ptr,
const std::string &  field_name,
RowColData  rc,
DofMethod method,
int  lower_rank,
int  upper_rank,
int  verb = DEFAULT_VERBOSITY 
)
pure virtual

#include <src/interfaces/Interface.hpp>

Make a loop over dofs.

Parameters
problem_ptrPointer to problem structure
field_nameName of the field containing DOFs to loop over
rcRow/column data specification
methodUser-defined method to execute on each DOF
lower_rankLower processor rank range (inclusive)
upper_rankUpper processor rank range (inclusive)
verbVerbosity level (default: DEFAULT_VERBOSITY)
Returns
MoFEMErrorCode Error code indicating success or failure

Implemented in MoFEM::CoreTmp< 0 >.

Examples
mofem/atom_tests/forces_and_sources_testing_edge_element.cpp, mofem/atom_tests/forces_and_sources_testing_flat_prism_element.cpp, mofem/atom_tests/hcurl_curl_operator.cpp, mofem/atom_tests/hdiv_divergence_operator.cpp, mofem/atom_tests/simple_interface.cpp, mofem/atom_tests/tensor_divergence_operator.cpp, mofem/tools/field_to_vertices.cpp, mofem/tutorials/adv-0/plastic.cpp, mofem/tutorials/adv-4/dynamic_first_order_con_law.cpp, mofem/tutorials/cor-10/navier_stokes.cpp, mofem/tutorials/cor-7/elasticity_mixed_formulation.cpp, mofem/tutorials/max-0/src/MagneticElement.hpp, mofem/tutorials/scl-12/electrostatics.cpp, mofem/tutorials/vec-2/src/NonlinearElasticExample.hpp, mofem/tutorials/vec-3/nonlinear_dynamic_elastic.cpp, mofem/tutorials/vec-4/approx_sphere.cpp, mofem/tutorials/vec-7/adjoint.cpp, mofem/users_modules/adolc-plasticity/adolc_plasticity.cpp, mofem/users_modules/basic_finite_elements/atom_tests/testing_jacobian_of_hook_element.cpp, mofem/users_modules/basic_finite_elements/atom_tests/testing_jacobian_of_hook_scaled_with_density_element.cpp, mofem/users_modules/basic_finite_elements/elasticity/elasticity.cpp, mofem/users_modules/basic_finite_elements/nonlinear_elasticity/nonlinear_dynamics.cpp, mofem/users_modules/bone_remodelling/src/impl/Remodeling.cpp, and mofem/users_modules/eshelbian_plasticity/src/impl/EshelbianPlasticity.cpp.

◆ loop_dofs() [2/3]

virtual MoFEMErrorCode MoFEM::CoreInterface::loop_dofs ( const std::string &  problem_name,
const std::string &  field_name,
RowColData  rc,
DofMethod method,
int  lower_rank,
int  upper_rank,
int  verb = DEFAULT_VERBOSITY 
)
pure virtual

#include <src/interfaces/Interface.hpp>

Make a loop over dofs.

Parameters
problem_nameName of the problem containing DOFs
field_nameName of the field containing DOFs to loop over
rcRow/column data specification
methodUser-defined method to execute on each DOF
lower_rankLower processor rank range (inclusive)
upper_rankUpper processor rank range (inclusive)
verbVerbosity level (default: DEFAULT_VERBOSITY)
Returns
MoFEMErrorCode Error code indicating success or failure

Implemented in MoFEM::CoreTmp< 0 >.

◆ loop_dofs() [3/3]

virtual MoFEMErrorCode MoFEM::CoreInterface::loop_dofs ( const std::string &  problem_name,
const std::string &  field_name,
RowColData  rc,
DofMethod method,
int  verb = DEFAULT_VERBOSITY 
)
pure virtual

#include <src/interfaces/Interface.hpp>

Make a loop over dofs.

Parameters
problem_nameName of the problem containing DOFs
field_nameName of the field containing DOFs to loop over
rcRow/column data specification
methodUser-defined method to execute on each DOF
verbVerbosity level (default: DEFAULT_VERBOSITY)
Returns
MoFEMErrorCode Error code indicating success or failure

Implemented in MoFEM::CoreTmp< 0 >.

◆ loop_finite_elements() [1/3]

virtual MoFEMErrorCode MoFEM::CoreInterface::loop_finite_elements ( const Problem problem_ptr,
const std::string &  fe_name,
FEMethod method,
int  lower_rank,
int  upper_rank,
boost::shared_ptr< NumeredEntFiniteElement_multiIndex fe_ptr = nullptr,
MoFEMTypes  bh = MF_EXIST,
CacheTupleWeakPtr  cache_ptr = CacheTupleSharedPtr(),
int  verb = DEFAULT_VERBOSITY 
)
pure virtual

#include <src/interfaces/Interface.hpp>

Make a loop over finite elements on partitions from upper to lower rank.

This function is like a Swiss knife; it can be used for post-processing or matrix and vector assembly. It iterates over the given finite element for a given problem. The particular methods executed on each element are provided by a class derived from Interface::FEMethod. At the beginning of each loop the user- defined function (method) preProcess() is called; for each element operator() is executed; at the end, the loop finalizes with the user-defined function (method) postProcess().

Note
If fe_ptr is given, it is expected that the multi-index is a subset of problem multi-index. If this is not the case behavior of the code is undetermined.

For more details please look to examples.

Interface::FEMethod

Parameters
problem_ptrpointer to problem consisting set of elements
fe_namename of element in problem
methodclass derived from Interface::FEMethod
lower_ranklower rank of process owning the finite element
upper_rankupper rank of process owning the finite element
fe_ptrpointer to finite elements multi-index
bhif bH = MF_EXIST, throws error if fe_name does not exist
cache_datacache data vector
cache_rowcache row vector
cache_colcache row vector
verbverbosity level
Returns
error code

Implemented in MoFEM::CoreTmp< 0 >.

◆ loop_finite_elements() [2/3]

virtual MoFEMErrorCode MoFEM::CoreInterface::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 
)
pure virtual

#include <src/interfaces/Interface.hpp>

Make a loop over finite elements.

This function is like a Swiss knife; it can be used for post-processing or matrix and vector assembly. It iterates over the given finite element for a given problem. The particular methods executed on each element are provided by a class derived from Interface::FEMethod. At the beginning of each loop the user- defined function (method) preProcess() is called; for each element operator() is executed; at the end, the loop finalizes with the user-defined function (method) postProcess().

Methods are executed only for local elements on a given processor.

For more details please look at examples.

Note
If fe_ptr is given, it is expected that the multi-index is a subset of problem multi-index. If this is not the case behavior of the code is undetermined.
Parameters
problem_nameproblem consisting set of elements
fe_namename of element in problem
methodclass derived form Interface::FEMethod
fe_ptrpointer to finite elements multi-index
bhif bH = MF_EXIST, throws error if fe_name does not exist
cache_tuple_ptrcache
verbverbosity level
Returns
error code

Implemented in MoFEM::CoreTmp< 0 >.

Examples
mofem/atom_tests/build_large_problem.cpp, mofem/atom_tests/continuity_check_on_contact_prism_side_ele.cpp, mofem/atom_tests/continuity_check_on_skeleton_3d.cpp, mofem/atom_tests/dm_partitioned_no_field.cpp, mofem/atom_tests/forces_and_sources_testing_edge_element.cpp, mofem/atom_tests/forces_and_sources_testing_flat_prism_element.cpp, mofem/atom_tests/forces_and_sources_testing_users_base.cpp, mofem/atom_tests/hcurl_divergence_operator_2d.cpp, mofem/atom_tests/prism_elements_from_surface.cpp, mofem/atom_tests/prism_polynomial_approximation.cpp, mofem/atom_tests/quad_polynomial_approximation.cpp, mofem/users_modules/basic_finite_elements/nonlinear_elasticity/nonlinear_dynamics.cpp, and mofem/users_modules/eshelbian_plasticity/src/impl/EshelbianPlasticity.cpp.

◆ loop_finite_elements() [3/3]

virtual MoFEMErrorCode MoFEM::CoreInterface::loop_finite_elements ( const std::string  problem_name,
const std::string &  fe_name,
FEMethod method,
int  lower_rank,
int  upper_rank,
boost::shared_ptr< NumeredEntFiniteElement_multiIndex fe_ptr = nullptr,
MoFEMTypes  bh = MF_EXIST,
CacheTupleWeakPtr  cache_ptr = CacheTupleSharedPtr(),
int  verb = DEFAULT_VERBOSITY 
)
pure virtual

#include <src/interfaces/Interface.hpp>

Make a loop over finite elements on partitions from upper to lower rank.

This function is like a Swiss knife; it can be used for post-processing or matrix and vector assembly. It iterates over the given finite element for a given problem. The particular methods executed on each element are provided by a class derived from Interface::FEMethod. At the beginning of each loop the user- defined function (method) preProcess() is called; for each element operator() is executed; at the end, the loop finalizes with the user-defined function (method) postProcess().

Note
If fe_ptr is given, it is expected that the multi-index is a subset of problem multi-index. If this is not the case behavior of the code is undetermined.

For more details please look to examples.

Parameters
problem_namepointer to problem consisting set of elements
fe_namename of element in problem
methodclass derived from Interface::FEMethod
lower_ranklower rank of process owning the finite element
upper_rankupper rank of process owning the finite element
fe_ptrpointer to finite elements multi-index
bhif bH = MF_EXIST, throws error if fe_name does not exist
cache_datacache data vector
cache_rowcache row vector
cache_colcache row vector
verbverbosity level
Returns
error code

Implemented in MoFEM::CoreTmp< 0 >.

◆ problem_basic_method_postProcess() [1/2]

virtual MoFEMErrorCode MoFEM::CoreInterface::problem_basic_method_postProcess ( const Problem problem_ptr,
BasicMethod method,
int  verb = DEFAULT_VERBOSITY 
)
pure virtual

#include <src/interfaces/Interface.hpp>

Set data for BasicMethod.

This function set data about problem, adjacencies and other multi-indices in database. This function can be used a special case when user need to do some pre- and post-processing before matrix or vector is initiated, or to assemble matrix for group of FEMethods. Is used by classes SnesCtx and TsCtx. Look for more details there.

FIXME: Here we need example

Parameters
pointerto problem data structure
methoduser method derived from BasicMethod

Implemented in MoFEM::CoreTmp< 0 >.

◆ problem_basic_method_postProcess() [2/2]

virtual MoFEMErrorCode MoFEM::CoreInterface::problem_basic_method_postProcess ( const std::string &  problem_name,
BasicMethod method,
int  verb = DEFAULT_VERBOSITY 
)
pure virtual

#include <src/interfaces/Interface.hpp>

Set data for BasicMethod.

This function set data about problem, adjacencies and other multi-indices in database. This function can be used a special case when user need to do some pre- and post-processing before matrix or vector is initiated, or to assemble matrix for group of FEMethods. Is used by classes SnesCtx and TsCtx. Look for more details there.

FIXME: Here we need example

Parameters
problem_namename of the problem
methoduser method derived from BasicMethod

Implemented in MoFEM::CoreTmp< 0 >.