v0.15.0
Loading...
Searching...
No Matches
Classes | Public Types | Public Member Functions | Private Member Functions | Private Attributes | List of all members
MoFEM::PipelineManager Struct Reference

PipelineManager interface. More...

#include "src/interfaces/PipelineManager.hpp"

Inheritance diagram for MoFEM::PipelineManager:
[legend]
Collaboration diagram for MoFEM::PipelineManager:
[legend]

Classes

struct  ElementsAndOpsByDim
 Template struct for dimension-specific finite element types. More...
 
struct  ElementsAndOpsByDim< 2 >
 Specialization for 2D problems. More...
 
struct  ElementsAndOpsByDim< 3 >
 Specialization for 3D problems. More...
 
struct  MeshsetFE
 

Public Types

enum  TSType { EX , IM , IM2 , IMEX }
 Enumeration of time solver types. More...
 
using UserDataOperator = MoFEM::ForcesAndSourcesCore::UserDataOperator
 
using RuleHookFun = MoFEM::ForcesAndSourcesCore::RuleHookFun
 
using VolEle = MoFEM::VolumeElementForcesAndSourcesCore
 
using FaceEle = MoFEM::FaceElementForcesAndSourcesCore
 
using EdgeEle = MoFEM::EdgeElementForcesAndSourcesCore
 

Public Member Functions

MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type-safe casting.
 
 PipelineManager (const MoFEM::Core &core)
 Construct PipelineManager.
 
boost::shared_ptr< FEMethod > & getDomainLhsFE ()
 Get domain left-hand side finite element.
 
boost::shared_ptr< FEMethod > & getDomainRhsFE ()
 Get domain right-hand side finite element.
 
boost::shared_ptr< FEMethod > & getBoundaryLhsFE ()
 Get boundary left-hand side finite element.
 
boost::shared_ptr< FEMethod > & getBoundaryRhsFE ()
 Get boundary right-hand side finite element.
 
boost::shared_ptr< FEMethod > & getSkeletonLhsFE ()
 Get skeleton left-hand side finite element.
 
boost::shared_ptr< FEMethod > & getSkeletonRhsFE ()
 Get skeleton right-hand side finite element.
 
boost::shared_ptr< FEMethod > & getDomainExplicitRhsFE ()
 Get domain explicit right-hand side finite element.
 
boost::shared_ptr< FEMethod > & getBoundaryExplicitRhsFE ()
 Get boundary explicit right-hand side finite element.
 
boost::shared_ptr< FEMethod > & getSkeletonExplicitRhsFE ()
 Get skeleton explicit right-hand side finite element.
 
boost::shared_ptr< FEMethod > & getMeshsetRhsFE ()
 Get meshset right-hand side finite element.
 
boost::shared_ptr< FEMethod > & getMeshsetLhsFE ()
 Get meshset left-hand side finite element.
 
boost::shared_ptr< FEMethod > & getMeshsetExplicitRhsFE ()
 Get meshset explicit right-hand side finite element.
 
template<typename T = ForcesAndSourcesCore, int DIM = -1>
auto getCastDomainLhsFE ()
 Get typed domain left-hand side finite element.
 
template<typename T = ForcesAndSourcesCore, int DIM = -1>
auto getCastDomainRhsFE ()
 Get typed domain right-hand side finite element.
 
template<typename T = ForcesAndSourcesCore, int DIM = -1>
auto getCastBoundaryLhsFE ()
 Get typed boundary left-hand side finite element.
 
template<typename T = ForcesAndSourcesCore, int DIM = -1>
auto getCastBoundaryRhsFE ()
 Get typed boundary right-hand side finite element.
 
template<typename T = ForcesAndSourcesCore, int DIM = -1>
auto getCastSkeletonLhsFE ()
 Get typed skeleton left-hand side finite element.
 
template<typename T = ForcesAndSourcesCore, int DIM = -1>
auto getCastSkeletonRhsFE ()
 Get typed skeleton right-hand side finite element.
 
template<typename T = ForcesAndSourcesCore, int DIM = -1>
auto getCastDomainExplicitRhsFE ()
 Get typed domain explicit right-hand side finite element.
 
template<typename T = ForcesAndSourcesCore, int DIM = -1>
auto getCastBoundaryExplicitRhsFE ()
 Get typed boundary explicit right-hand side finite element.
 
template<typename T = ForcesAndSourcesCore, int DIM = -1>
auto getCastSkeletonExplicitRhsFE ()
 Get typed skeleton explicit right-hand side finite element.
 
template<typename T = ForcesAndSourcesCore>
auto getCastMeshsetRhsFE ()
 Get typed meshset right-hand side finite element.
 
template<typename T = ForcesAndSourcesCore>
auto getCastMeshsetLhsFE ()
 Get typed meshset left-hand side finite element.
 
template<typename T = ForcesAndSourcesCore>
auto getCastMeshsetExplicitRhsFE ()
 Get typed meshset explicit right-hand side finite element.
 
template<int DIM = -1>
MoFEMErrorCode setDomainLhsIntegrationRule (RuleHookFun rule)
 Set integration rule for domain left-hand side finite element.
 
template<int DIM = -1>
MoFEMErrorCode setDomainRhsIntegrationRule (RuleHookFun rule)
 Set integration rule for domain right-hand side finite element.
 
template<int DIM = -1>
MoFEMErrorCode setBoundaryLhsIntegrationRule (RuleHookFun rule)
 Set integration rule for boundary left-hand side finite element.
 
template<int DIM = -1>
MoFEMErrorCode setBoundaryRhsIntegrationRule (RuleHookFun rule)
 Set integration rule for boundary right-hand side finite element.
 
template<int DIM = -1>
MoFEMErrorCode setSkeletonLhsIntegrationRule (RuleHookFun rule)
 Set integration rule for skeleton left-hand side finite element.
 
template<int DIM = -1>
MoFEMErrorCode setSkeletonRhsIntegrationRule (RuleHookFun rule)
 Set integration rule for skeleton right-hand side finite element.
 
template<int DIM = -1>
MoFEMErrorCode setDomainExplicitRhsIntegrationRule (RuleHookFun rule)
 Set integration rule for domain explicit right-hand side finite element.
 
template<int DIM = -1>
MoFEMErrorCode setBoundaryExplicitRhsIntegrationRule (RuleHookFun rule)
 Set integration rule for boundary explicit right-hand side finite element.
 
template<int DIM = -1>
MoFEMErrorCode setSkeletonExplicitRhsIntegrationRule (RuleHookFun rule)
 Set integration rule for skeleton explicit right-hand side finite element.
 
template<int DIM = -1>
boost::ptr_deque< UserDataOperator > & getOpDomainLhsPipeline ()
 Get the Op Domain Lhs Pipeline object.
 
template<int DIM = -1>
boost::ptr_deque< UserDataOperator > & getOpDomainRhsPipeline ()
 Get the Op Domain Rhs Pipeline object.
 
template<int DIM = -1>
boost::ptr_deque< UserDataOperator > & getOpBoundaryLhsPipeline ()
 Get the Op Boundary Lhs Pipeline object.
 
template<int DIM = -1>
boost::ptr_deque< UserDataOperator > & getOpBoundaryRhsPipeline ()
 Get the Op Boundary Rhs Pipeline object.
 
template<int DIM = -1>
boost::ptr_deque< UserDataOperator > & getOpSkeletonLhsPipeline ()
 Get the Op Skeleton Lhs Pipeline object.
 
template<int DIM = -1>
boost::ptr_deque< UserDataOperator > & getOpSkeletonRhsPipeline ()
 Get the Op Skeleton Rhs Pipeline object.
 
template<int DIM = -1>
boost::ptr_deque< UserDataOperator > & getOpDomainExplicitRhsPipeline ()
 Get the Op Domain Rhs Pipeline object for implicit-explicit G term.
 
template<int DIM = -1>
boost::ptr_deque< UserDataOperator > & getOpBoundaryExplicitRhsPipeline ()
 Get the Op Boundary Rhs Pipeline object for implicit-explicit G term.
 
template<int DIM = -1>
boost::ptr_deque< UserDataOperator > & getOpSkeletonExplicitRhsPipeline ()
 Get the Op Skeleton Rhs Pipeline object for implicit-explicit G term.
 
boost::ptr_deque< UserDataOperator > & getOpMeshsetRhsPipeline ()
 Get the Op Meshset Rhs Pipeline object.
 
boost::ptr_deque< UserDataOperator > & getOpMeshsetLhsPipeline ()
 Get the Op Meshset Lhs Pipeline object.
 
boost::ptr_deque< UserDataOperator > & getOpMeshsetExplicitRhsPipeline ()
 Get the Op Meshset Explicit Rhs Pipeline object.
 
MoFEMErrorCode loopFiniteElements (SmartPetscObj< DM > dm=nullptr)
 Iterate finite elements.
 
SmartPetscObj< KSP > createKSP (SmartPetscObj< DM > dm=nullptr)
 Create KSP (linear) solver.
 
SmartPetscObj< SNES > createSNES (SmartPetscObj< DM > dm=nullptr)
 Create SNES (nonlinear) solver.
 
SmartPetscObj< TS > createTS (const TSType type, SmartPetscObj< DM > dm=nullptr)
 Create TS (time) solver with specified type.
 
SmartPetscObj< TS > createTSEX (SmartPetscObj< DM > dm=nullptr)
 Create TS (time) explicit solver.
 
SmartPetscObj< TS > createTSIM (SmartPetscObj< DM > dm=nullptr)
 Create TS (time) implicit solver.
 
DEPRECATED auto createTS (SmartPetscObj< DM > dm=nullptr)
 
SmartPetscObj< TS > createTSIM2 (SmartPetscObj< DM > dm=nullptr)
 Create TS (time) solver for second order equation in time.
 
DEPRECATED auto createTS2 (SmartPetscObj< DM > dm=nullptr)
 
SmartPetscObj< TS > createTSIMEX (SmartPetscObj< DM > dm=nullptr)
 Create TS (time) implicit-explicit solver.
 
template<>
MoFEMErrorCode setDomainLhsIntegrationRule (PipelineManager::RuleHookFun rule)
 
template<>
MoFEMErrorCode setDomainRhsIntegrationRule (PipelineManager::RuleHookFun rule)
 
template<>
MoFEMErrorCode setBoundaryLhsIntegrationRule (PipelineManager::RuleHookFun rule)
 
template<>
MoFEMErrorCode setBoundaryRhsIntegrationRule (PipelineManager::RuleHookFun rule)
 
template<>
MoFEMErrorCode setSkeletonLhsIntegrationRule (PipelineManager::RuleHookFun rule)
 
template<>
MoFEMErrorCode setSkeletonRhsIntegrationRule (PipelineManager::RuleHookFun rule)
 
template<>
MoFEMErrorCode setDomainExplicitRhsIntegrationRule (PipelineManager::RuleHookFun rule)
 
template<>
MoFEMErrorCode setBoundaryExplicitRhsIntegrationRule (PipelineManager::RuleHookFun rule)
 
template<>
MoFEMErrorCode setSkeletonExplicitRhsIntegrationRule (PipelineManager::RuleHookFun rule)
 
template<>
boost::ptr_deque< PipelineManager::UserDataOperator > & getOpDomainLhsPipeline ()
 
template<>
boost::ptr_deque< PipelineManager::UserDataOperator > & getOpDomainRhsPipeline ()
 
template<>
boost::ptr_deque< PipelineManager::UserDataOperator > & getOpBoundaryLhsPipeline ()
 
template<>
boost::ptr_deque< PipelineManager::UserDataOperator > & getOpBoundaryRhsPipeline ()
 
template<>
boost::ptr_deque< PipelineManager::UserDataOperator > & getOpSkeletonLhsPipeline ()
 
template<>
boost::ptr_deque< PipelineManager::UserDataOperator > & getOpSkeletonRhsPipeline ()
 
template<>
boost::ptr_deque< PipelineManager::UserDataOperator > & getOpDomainExplicitRhsPipeline ()
 
template<>
boost::ptr_deque< PipelineManager::UserDataOperator > & getOpBoundaryExplicitRhsPipeline ()
 
template<>
boost::ptr_deque< PipelineManager::UserDataOperator > & getOpSkeletonExplicitRhsPipeline ()
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface reference to pointer of interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface.
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface.
 
virtual ~UnknownInterface ()=default
 

Private Member Functions

template<int DIM>
boost::shared_ptr< FEMethod > & createDomainFEPipeline (boost::shared_ptr< FEMethod > &fe)
 Create domain finite element pipeline based on dimension.
 
template<int DIM>
boost::shared_ptr< FEMethod > & createBoundaryFEPipeline (boost::shared_ptr< FEMethod > &fe)
 Create boundary finite element pipeline based on dimension.
 
boost::shared_ptr< FEMethod > & createMeshsetFEPipeline (boost::shared_ptr< FEMethod > &fe)
 Create meshset finite element pipeline.
 
template<>
boost::shared_ptr< FEMethod > & createDomainFEPipeline (boost::shared_ptr< FEMethod > &fe)
 
template<>
boost::shared_ptr< FEMethod > & createDomainFEPipeline (boost::shared_ptr< FEMethod > &fe)
 
template<>
boost::shared_ptr< FEMethod > & createDomainFEPipeline (boost::shared_ptr< FEMethod > &fe)
 
template<>
boost::shared_ptr< FEMethod > & createDomainFEPipeline (boost::shared_ptr< FEMethod > &fe)
 
template<>
boost::shared_ptr< FEMethod > & createBoundaryFEPipeline (boost::shared_ptr< FEMethod > &fe)
 
template<>
boost::shared_ptr< FEMethod > & createBoundaryFEPipeline (boost::shared_ptr< FEMethod > &fe)
 
template<>
boost::shared_ptr< FEMethod > & createBoundaryFEPipeline (boost::shared_ptr< FEMethod > &fe)
 
template<>
boost::shared_ptr< FEMethod > & createBoundaryFEPipeline (boost::shared_ptr< FEMethod > &fe)
 

Private Attributes

MoFEM::CorecOre
 Reference to MoFEM core instance.
 
boost::shared_ptr< FEMethodfeDomainRhs
 Element to assemble RHS side by integrating domain.
 
boost::shared_ptr< FEMethodfeDomainLhs
 Element to assemble LHS side by integrating domain.
 
boost::shared_ptr< FEMethodfeBoundaryRhs
 Element to assemble RHS side by integrating boundary.
 
boost::shared_ptr< FEMethodfeBoundaryLhs
 Element to assemble LHS side by integrating boundary.
 
boost::shared_ptr< FEMethodfeSkeletonRhs
 Element to assemble RHS side by integrating skeleton.
 
boost::shared_ptr< FEMethodfeSkeletonLhs
 Element to assemble LHS side by integrating skeleton.
 
boost::shared_ptr< FEMethodfeDomainExplicitRhs
 Element to assemble explicit Rhs for IMEX solver.
 
boost::shared_ptr< FEMethodfeBoundaryExplicitRhs
 
boost::shared_ptr< FEMethodfeSkeletonExplicitRhs
 
boost::shared_ptr< FEMethodfeMeshsetRhs
 
boost::shared_ptr< FEMethodfeMeshsetLhs
 
boost::shared_ptr< FEMethodfeMeshsetExplicitRhs
 

Additional Inherited Members

- Static Public Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version.
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version.
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version.
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version.
 

Detailed Description

PipelineManager interface.

Examples
mofem/atom_tests/child_and_parent.cpp, mofem/atom_tests/dg_projection.cpp, mofem/atom_tests/hanging_node_approx.cpp, mofem/atom_tests/hcurl_check_approx_in_2d.cpp, mofem/atom_tests/hcurl_curl_operator.cpp, mofem/atom_tests/hdiv_check_approx_in_3d.cpp, mofem/atom_tests/hdiv_divergence_operator.cpp, mofem/atom_tests/higher_derivatives.cpp, mofem/atom_tests/operators_tests.cpp, mofem/atom_tests/scalar_check_approximation.cpp, mofem/atom_tests/schur_test_diag_mat.cpp, mofem/atom_tests/simple_interface.cpp, mofem/atom_tests/simple_l2_only.cpp, mofem/atom_tests/tensor_divergence_operator.cpp, mofem/atom_tests/test_broken_space.cpp, mofem/tutorials/adv-0/plastic.cpp, mofem/tutorials/adv-2/thermo_elastic.cpp, mofem/tutorials/adv-3/level_set.cpp, mofem/tutorials/adv-4/dynamic_first_order_con_law.cpp, mofem/tutorials/adv-5/seepage.cpp, mofem/tutorials/clx-0/helmholtz.cpp, mofem/tutorials/fun-0/hello_world.cpp, mofem/tutorials/fun-1/integration.cpp, mofem/tutorials/fun-2/plot_base.cpp, mofem/tutorials/mix-0/mixed_poisson.cpp, mofem/tutorials/mix-1/phase.cpp, mofem/tutorials/scl-0/approximaton.cpp, mofem/tutorials/scl-1/poisson_2d_homogeneous.cpp, mofem/tutorials/scl-10/initial_diffusion.cpp, mofem/tutorials/scl-10/photon_diffusion.cpp, mofem/tutorials/scl-11/poisson_2d_dis_galerkin.cpp, mofem/tutorials/scl-6/heat_equation.cpp, mofem/tutorials/scl-7/wave_equation.cpp, mofem/tutorials/scl-8/radiation.cpp, mofem/tutorials/scl-9/heat_method.cpp, mofem/tutorials/vec-1/eigen_elastic.cpp, mofem/tutorials/vec-10/schur_elastic.cpp, mofem/tutorials/vec-3/nonlinear_dynamic_elastic.cpp, mofem/tutorials/vec-4/approx_sphere.cpp, mofem/tutorials/vec-4/shallow_wave.cpp, mofem/tutorials/vec-5/free_surface.cpp, mofem/tutorials/vec-6/plate.cpp, mofem/tutorials/vec-7/adjoint.cpp, and mofem/users_modules/adolc-plasticity/adolc_plasticity.cpp.

Definition at line 24 of file PipelineManager.hpp.

Member Typedef Documentation

◆ EdgeEle

Definition at line 48 of file PipelineManager.hpp.

◆ FaceEle

Definition at line 47 of file PipelineManager.hpp.

◆ RuleHookFun

Definition at line 44 of file PipelineManager.hpp.

◆ UserDataOperator

Definition at line 43 of file PipelineManager.hpp.

◆ VolEle

Definition at line 46 of file PipelineManager.hpp.

Member Enumeration Documentation

◆ TSType

Enumeration of time solver types.

Enumerator
EX 

Explicit time integration.

IM 

Implicit time integration.

IM2 

Second-order implicit time integration.

IMEX 

Implicit-explicit time integration.

Definition at line 488 of file PipelineManager.hpp.

488 {
489 EX, ///< Explicit time integration
490 IM, ///< Implicit time integration
491 IM2, ///< Second-order implicit time integration
492 IMEX ///< Implicit-explicit time integration
493 };
@ IM
Implicit time integration.
@ IM2
Second-order implicit time integration.
@ IMEX
Implicit-explicit time integration.
@ EX
Explicit time integration.

Constructor & Destructor Documentation

◆ PipelineManager()

MoFEM::PipelineManager::PipelineManager ( const MoFEM::Core core)

Construct PipelineManager.

Parameters
coreReference to MoFEM core instance

Definition at line 60 of file PipelineManager.cpp.

61 : cOre(const_cast<Core &>(core)) {}
CoreTmp< 0 > Core
Definition Core.hpp:1148
MoFEM::Core & cOre
Reference to MoFEM core instance.

Member Function Documentation

◆ createBoundaryFEPipeline() [1/5]

template<int DIM>
boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::createBoundaryFEPipeline ( boost::shared_ptr< FEMethod > &  fe)
inlineprivate

Create boundary finite element pipeline based on dimension.

Template Parameters
DIMProblem dimension (1, 2, 3, or -1 for automatic detection)
Parameters
feReference to finite element method shared pointer
Returns
boost::shared_ptr<FEMethod>& Reference to created finite element

Definition at line 674 of file PipelineManager.hpp.

674 {
675 static_assert(DIM == 1 || DIM == 2 || DIM == 3, "not implemented");
676 fe = boost::make_shared<FEMethod>();
677 return fe;
678}

◆ createBoundaryFEPipeline() [2/5]

template<>
boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::createBoundaryFEPipeline ( boost::shared_ptr< FEMethod > &  fe)
inlineprivate

Definition at line 674 of file PipelineManager.hpp.

682 {
683 if (!fe)
684 fe = boost::make_shared<FaceElementForcesAndSourcesCore>(cOre);
685 return fe;
686}

◆ createBoundaryFEPipeline() [3/5]

template<>
boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::createBoundaryFEPipeline ( boost::shared_ptr< FEMethod > &  fe)
inlineprivate

Definition at line 674 of file PipelineManager.hpp.

690 {
691 if (!fe)
692 fe = boost::make_shared<EdgeEle>(cOre);
693 return fe;
694}

◆ createBoundaryFEPipeline() [4/5]

template<>
boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::createBoundaryFEPipeline ( boost::shared_ptr< FEMethod > &  fe)
inlineprivate

Definition at line 674 of file PipelineManager.hpp.

698 {
699 if (!fe)
700 fe = boost::make_shared<VertexElementForcesAndSourcesCore>(cOre);
701 return fe;
702}

◆ createBoundaryFEPipeline() [5/5]

template<>
boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::createBoundaryFEPipeline ( boost::shared_ptr< FEMethod > &  fe)
inlineprivate

Definition at line 748 of file PipelineManager.hpp.

754 {
755 switch (cOre.getInterface<Simple>()->getDim()) {
756 case 1:
757 return createBoundaryFEPipeline<1>(fe);
758 case 2:
759 return createBoundaryFEPipeline<2>(fe);
760 case 3:
761 return createBoundaryFEPipeline<3>(fe);
762 default:
763 THROW_MESSAGE("Not implemented");
764 }
765}
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.

◆ createDomainFEPipeline() [1/5]

template<int DIM>
boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::createDomainFEPipeline ( boost::shared_ptr< FEMethod > &  fe)
inlineprivate

Create domain finite element pipeline based on dimension.

Template Parameters
DIMProblem dimension (1, 2, 3, or -1 for automatic detection)
Parameters
feReference to finite element method shared pointer
Returns
boost::shared_ptr<FEMethod>& Reference to created finite element

Definition at line 627 of file PipelineManager.hpp.

627 {
628 static_assert(DIM == 1 || DIM == 2 || DIM == 3, "not implemented");
629 fe = boost::make_shared<FEMethod>();
630 return fe;
631}

◆ createDomainFEPipeline() [2/5]

template<>
boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::createDomainFEPipeline ( boost::shared_ptr< FEMethod > &  fe)
inlineprivate

Definition at line 627 of file PipelineManager.hpp.

643 {
644 if (!fe)
645 fe = boost::make_shared<FaceEle>(cOre);
646 return fe;
647}

◆ createDomainFEPipeline() [3/5]

template<>
boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::createDomainFEPipeline ( boost::shared_ptr< FEMethod > &  fe)
inlineprivate

Definition at line 627 of file PipelineManager.hpp.

659 {
660 switch (cOre.getInterface<Simple>()->getDim()) {
661 case 1:
662 return createDomainFEPipeline<1>(fe);
663 case 2:
664 return createDomainFEPipeline<2>(fe);
665 case 3:
666 return createDomainFEPipeline<3>(fe);
667 default:
668 THROW_MESSAGE("Not implemented");
669 }
670}

◆ createDomainFEPipeline() [4/5]

template<>
boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::createDomainFEPipeline ( boost::shared_ptr< FEMethod > &  fe)
inlineprivate

Definition at line 627 of file PipelineManager.hpp.

651 {
652 if (!fe)
653 fe = boost::make_shared<EdgeEle>(cOre);
654 return fe;
655}

◆ createDomainFEPipeline() [5/5]

template<>
boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::createDomainFEPipeline ( boost::shared_ptr< FEMethod > &  fe)
inlineprivate

Definition at line 627 of file PipelineManager.hpp.

635 {
636 if (!fe)
637 fe = boost::make_shared<VolEle>(cOre);
638 return fe;
639}

◆ createMeshsetFEPipeline()

boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::createMeshsetFEPipeline ( boost::shared_ptr< FEMethod > &  fe)
inlineprivate

Create meshset finite element pipeline.

Parameters
feReference to finite element method shared pointer
Returns
boost::shared_ptr<FEMethod>& Reference to created finite element

Definition at line 26 of file PipelineManager.cpp.

26 {
27 if (!fe)
28 fe = boost::make_shared<MeshsetFE>(cOre);
29 return fe;
30}

◆ createTS() [1/2]

SmartPetscObj< TS > MoFEM::PipelineManager::createTS ( const TSType  type,
SmartPetscObj< DM >  dm = nullptr 
)

Create TS (time) solver with specified type.

Parameters
typeType of time solver (EX/IM/IM2/IMEX)
dmOptional DM object (default: nullptr)
Returns
SmartPetscObj<TS> Smart pointer to PETSc TS solver

Definition at line 231 of file PipelineManager.cpp.

232 {
233 switch (type) {
234 case EX:
235 return createTSEX(dm);
236 break;
237 case IM:
238 return createTSIM(dm);
239 break;
240 case IM2:
241 return createTSIM2(dm);
242 break;
243 case IMEX:
244 return createTSIMEX(dm);
245 break;
246 default:
248 "TS solver handling not implemented");
249 break;
250 }
251 return SmartPetscObj<TS>();
252}
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
SmartPetscObj< TS > createTSIM(SmartPetscObj< DM > dm=nullptr)
Create TS (time) implicit solver.
SmartPetscObj< TS > createTSIM2(SmartPetscObj< DM > dm=nullptr)
Create TS (time) solver for second order equation in time.
SmartPetscObj< TS > createTSIMEX(SmartPetscObj< DM > dm=nullptr)
Create TS (time) implicit-explicit solver.
SmartPetscObj< TS > createTSEX(SmartPetscObj< DM > dm=nullptr)
Create TS (time) explicit solver.

◆ createTS() [2/2]

DEPRECATED auto MoFEM::PipelineManager::createTS ( SmartPetscObj< DM >  dm = nullptr)
inline
Deprecated:
Use version with explicit TS solver type

Definition at line 525 of file PipelineManager.hpp.

525 {
526 return createTSIM(dm);
527 }

◆ createTS2()

DEPRECATED auto MoFEM::PipelineManager::createTS2 ( SmartPetscObj< DM >  dm = nullptr)
inline
Deprecated:
Change name. Use createTSIM2 instead.

Definition at line 541 of file PipelineManager.hpp.

541 {
542 return createTSIM2(dm);
543 }

◆ getBoundaryExplicitRhsFE()

boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::getBoundaryExplicitRhsFE ( )
inline

Get boundary explicit right-hand side finite element.

Returns
boost::shared_ptr<FEMethod>& Reference to boundary explicit RHS finite element for IMEX schemes

Definition at line 732 of file PipelineManager.hpp.

732 {
734}
boost::shared_ptr< FEMethod > feBoundaryExplicitRhs

◆ getBoundaryLhsFE()

boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::getBoundaryLhsFE ( )
inline

Get boundary left-hand side finite element.

Returns
boost::shared_ptr<FEMethod>& Reference to boundary LHS finite element for matrix assembly
Examples
mofem/tutorials/clx-0/helmholtz.cpp, and mofem/tutorials/scl-8/radiation.cpp.

Definition at line 712 of file PipelineManager.hpp.

712 {
713 return feBoundaryLhs;
714}
boost::shared_ptr< FEMethod > feBoundaryLhs
Element to assemble LHS side by integrating boundary.

◆ getBoundaryRhsFE()

boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::getBoundaryRhsFE ( )
inline

Get boundary right-hand side finite element.

Returns
boost::shared_ptr<FEMethod>& Reference to boundary RHS finite element for vector assembly
Examples
mofem/atom_tests/child_and_parent.cpp, mofem/atom_tests/hanging_node_approx.cpp, mofem/tutorials/clx-0/helmholtz.cpp, mofem/tutorials/mix-1/phase.cpp, mofem/tutorials/scl-6/heat_equation.cpp, mofem/tutorials/scl-7/wave_equation.cpp, and mofem/tutorials/scl-8/radiation.cpp.

Definition at line 716 of file PipelineManager.hpp.

716 {
717 return feBoundaryRhs;
718}
boost::shared_ptr< FEMethod > feBoundaryRhs
Element to assemble RHS side by integrating boundary.

◆ getCastBoundaryExplicitRhsFE()

template<typename T , int DIM>
auto MoFEM::PipelineManager::getCastBoundaryExplicitRhsFE ( )
inline

Get typed boundary explicit right-hand side finite element.

Template Parameters
TTarget finite element type (default: ForcesAndSourcesCore)
DIMDimension (-1 for automatic detection from problem)
Returns
auto Shared pointer to typed boundary explicit RHS finite element

Definition at line 806 of file PipelineManager.hpp.

806 {
807 return boost::dynamic_pointer_cast<T>(
808 createBoundaryFEPipeline<DIM>(feBoundaryExplicitRhs));
809}

◆ getCastBoundaryLhsFE()

template<typename T , int DIM>
auto MoFEM::PipelineManager::getCastBoundaryLhsFE ( )
inline

Get typed boundary left-hand side finite element.

Template Parameters
TTarget finite element type (default: ForcesAndSourcesCore)
DIMDimension (-1 for automatic detection from problem)
Returns
auto Shared pointer to typed boundary LHS finite element

Definition at line 777 of file PipelineManager.hpp.

777 {
778 return boost::dynamic_pointer_cast<T>(
779 createBoundaryFEPipeline<DIM>(feBoundaryLhs));
780}

◆ getCastBoundaryRhsFE()

template<typename T , int DIM>
auto MoFEM::PipelineManager::getCastBoundaryRhsFE ( )
inline

Get typed boundary right-hand side finite element.

Template Parameters
TTarget finite element type (default: ForcesAndSourcesCore)
DIMDimension (-1 for automatic detection from problem)
Returns
auto Shared pointer to typed boundary RHS finite element

Definition at line 783 of file PipelineManager.hpp.

783 {
784 return boost::dynamic_pointer_cast<T>(
785 createBoundaryFEPipeline<DIM>(feBoundaryRhs));
786}

◆ getCastDomainExplicitRhsFE()

template<typename T , int DIM>
auto MoFEM::PipelineManager::getCastDomainExplicitRhsFE ( )
inline

Get typed domain explicit right-hand side finite element.

Template Parameters
TTarget finite element type (default: ForcesAndSourcesCore)
DIMDimension (-1 for automatic detection from problem)
Returns
auto Shared pointer to typed domain explicit RHS finite element

Definition at line 800 of file PipelineManager.hpp.

800 {
801 return boost::dynamic_pointer_cast<T>(
802 createDomainFEPipeline<DIM>(feDomainExplicitRhs));
803}
boost::shared_ptr< FEMethod > feDomainExplicitRhs
Element to assemble explicit Rhs for IMEX solver.

◆ getCastDomainLhsFE()

template<typename T , int DIM>
auto MoFEM::PipelineManager::getCastDomainLhsFE ( )
inline

Get typed domain left-hand side finite element.

Template Parameters
TTarget finite element type (default: ForcesAndSourcesCore)
DIMDimension (-1 for automatic detection from problem)
Returns
auto Shared pointer to typed domain LHS finite element

Definition at line 767 of file PipelineManager.hpp.

767 {
768 return boost::dynamic_pointer_cast<T>(
769 createDomainFEPipeline<DIM>(feDomainLhs));
770}
boost::shared_ptr< FEMethod > feDomainLhs
Element to assemble LHS side by integrating domain.

◆ getCastDomainRhsFE()

template<typename T , int DIM>
auto MoFEM::PipelineManager::getCastDomainRhsFE ( )
inline

Get typed domain right-hand side finite element.

Template Parameters
TTarget finite element type (default: ForcesAndSourcesCore)
DIMDimension (-1 for automatic detection from problem)
Returns
auto Shared pointer to typed domain RHS finite element

Definition at line 772 of file PipelineManager.hpp.

772 {
773 return boost::dynamic_pointer_cast<T>(
774 createDomainFEPipeline<DIM>(feDomainRhs));
775}
boost::shared_ptr< FEMethod > feDomainRhs
Element to assemble RHS side by integrating domain.

◆ getCastMeshsetExplicitRhsFE()

template<typename T >
auto MoFEM::PipelineManager::getCastMeshsetExplicitRhsFE ( )
inline

Get typed meshset explicit right-hand side finite element.

Template Parameters
TTarget finite element type (default: ForcesAndSourcesCore)
Returns
auto Shared pointer to typed meshset explicit RHS finite element

Definition at line 825 of file PipelineManager.hpp.

825 {
826 return boost::dynamic_pointer_cast<T>(
828}
boost::shared_ptr< FEMethod > feMeshsetExplicitRhs
boost::shared_ptr< FEMethod > & createMeshsetFEPipeline(boost::shared_ptr< FEMethod > &fe)
Create meshset finite element pipeline.

◆ getCastMeshsetLhsFE()

template<typename T >
auto MoFEM::PipelineManager::getCastMeshsetLhsFE ( )
inline

Get typed meshset left-hand side finite element.

Template Parameters
TTarget finite element type (default: ForcesAndSourcesCore)
Returns
auto Shared pointer to typed meshset LHS finite element

Definition at line 821 of file PipelineManager.hpp.

821 {
822 return boost::dynamic_pointer_cast<T>(createMeshsetFEPipeline(feMeshsetLhs));
823}
boost::shared_ptr< FEMethod > feMeshsetLhs

◆ getCastMeshsetRhsFE()

template<typename T >
auto MoFEM::PipelineManager::getCastMeshsetRhsFE ( )
inline

Get typed meshset right-hand side finite element.

Template Parameters
TTarget finite element type (default: ForcesAndSourcesCore)
Returns
auto Shared pointer to typed meshset RHS finite element

Definition at line 817 of file PipelineManager.hpp.

817 {
818 return boost::dynamic_pointer_cast<T>(createMeshsetFEPipeline(feMeshsetRhs));
819}
boost::shared_ptr< FEMethod > feMeshsetRhs

◆ getCastSkeletonExplicitRhsFE()

template<typename T , int DIM>
auto MoFEM::PipelineManager::getCastSkeletonExplicitRhsFE ( )
inline

Get typed skeleton explicit right-hand side finite element.

Template Parameters
TTarget finite element type (default: ForcesAndSourcesCore)
DIMDimension (-1 for automatic detection from problem)
Returns
auto Shared pointer to typed skeleton explicit RHS finite element

Definition at line 812 of file PipelineManager.hpp.

812 {
813 return boost::dynamic_pointer_cast<T>(
814 createBoundaryFEPipeline<DIM>(feSkeletonExplicitRhs));
815}
boost::shared_ptr< FEMethod > feSkeletonExplicitRhs

◆ getCastSkeletonLhsFE()

template<typename T , int DIM>
auto MoFEM::PipelineManager::getCastSkeletonLhsFE ( )
inline

Get typed skeleton left-hand side finite element.

Template Parameters
TTarget finite element type (default: ForcesAndSourcesCore)
DIMDimension (-1 for automatic detection from problem)
Returns
auto Shared pointer to typed skeleton LHS finite element

Definition at line 789 of file PipelineManager.hpp.

789 {
790 return boost::dynamic_pointer_cast<T>(
791 createBoundaryFEPipeline<DIM>(feSkeletonLhs));
792}
boost::shared_ptr< FEMethod > feSkeletonLhs
Element to assemble LHS side by integrating skeleton.

◆ getCastSkeletonRhsFE()

template<typename T , int DIM>
auto MoFEM::PipelineManager::getCastSkeletonRhsFE ( )
inline

Get typed skeleton right-hand side finite element.

Template Parameters
TTarget finite element type (default: ForcesAndSourcesCore)
DIMDimension (-1 for automatic detection from problem)
Returns
auto Shared pointer to typed skeleton RHS finite element

Definition at line 794 of file PipelineManager.hpp.

794 {
795 return boost::dynamic_pointer_cast<T>(
796 createBoundaryFEPipeline<DIM>(feSkeletonRhs));
797}
boost::shared_ptr< FEMethod > feSkeletonRhs
Element to assemble RHS side by integrating skeleton.

◆ getDomainExplicitRhsFE()

boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::getDomainExplicitRhsFE ( )
inline

Get domain explicit right-hand side finite element.

Returns
boost::shared_ptr<FEMethod>& Reference to domain explicit RHS finite element for IMEX schemes

Definition at line 728 of file PipelineManager.hpp.

728 {
729 return feDomainExplicitRhs;
730}

◆ getDomainLhsFE()

boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::getDomainLhsFE ( )
inline

◆ getDomainRhsFE()

boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::getDomainRhsFE ( )
inline

◆ getMeshsetExplicitRhsFE()

boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::getMeshsetExplicitRhsFE ( )
inline

Get meshset explicit right-hand side finite element.

Returns
boost::shared_ptr<FEMethod>& Reference to meshset explicit RHS finite element for IMEX schemes

Definition at line 748 of file PipelineManager.hpp.

748 {
750}

◆ getMeshsetLhsFE()

boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::getMeshsetLhsFE ( )
inline

Get meshset left-hand side finite element.

Returns
boost::shared_ptr<FEMethod>& Reference to meshset LHS finite element for matrix assembly

Definition at line 744 of file PipelineManager.hpp.

744 {
745 return feMeshsetLhs;
746}

◆ getMeshsetRhsFE()

boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::getMeshsetRhsFE ( )
inline

Get meshset right-hand side finite element.

Returns
boost::shared_ptr<FEMethod>& Reference to meshset RHS finite element for vector assembly

Definition at line 740 of file PipelineManager.hpp.

740 {
741 return feMeshsetRhs;
742}

◆ getOpBoundaryExplicitRhsPipeline()

template<>
boost::ptr_deque< PipelineManager::UserDataOperator > & MoFEM::PipelineManager::getOpBoundaryExplicitRhsPipeline ( )
inline

Definition at line 1245 of file PipelineManager.hpp.

1253 {
1254 switch (cOre.getInterface<Simple>()->getDim()) {
1255 case 1:
1256 return getOpBoundaryExplicitRhsPipeline<1>();
1257 case 2:
1258 return getOpBoundaryExplicitRhsPipeline<2>();
1259 case 3:
1260 break;
1261 default:
1262 THROW_MESSAGE("Not implemented");
1263 }
1264 return getOpBoundaryExplicitRhsPipeline<3>();
1265}

◆ getOpBoundaryLhsPipeline()

template<>
boost::ptr_deque< PipelineManager::UserDataOperator > & MoFEM::PipelineManager::getOpBoundaryLhsPipeline ( )
inline

Definition at line 1125 of file PipelineManager.hpp.

1133 {
1134 switch (cOre.getInterface<Simple>()->getDim()) {
1135 case 1:
1136 return getOpBoundaryLhsPipeline<1>();
1137 case 2:
1138 return getOpBoundaryLhsPipeline<2>();
1139 case 3:
1140 break;
1141 default:
1142 THROW_MESSAGE("Not implemented");
1143 }
1144 return getOpBoundaryLhsPipeline<3>();
1145}

◆ getOpBoundaryRhsPipeline()

template<>
boost::ptr_deque< PipelineManager::UserDataOperator > & MoFEM::PipelineManager::getOpBoundaryRhsPipeline ( )
inline

Definition at line 1149 of file PipelineManager.hpp.

1157 {
1158 switch (cOre.getInterface<Simple>()->getDim()) {
1159 case 1:
1160 return getOpBoundaryRhsPipeline<1>();
1161 case 2:
1162 return getOpBoundaryRhsPipeline<2>();
1163 case 3:
1164 break;
1165 default:
1166 THROW_MESSAGE("Not implemented");
1167 }
1168 return getOpBoundaryRhsPipeline<3>();
1169}

◆ getOpDomainExplicitRhsPipeline()

template<>
boost::ptr_deque< PipelineManager::UserDataOperator > & MoFEM::PipelineManager::getOpDomainExplicitRhsPipeline ( )
inline

Definition at line 1221 of file PipelineManager.hpp.

1229 {
1230 switch (cOre.getInterface<Simple>()->getDim()) {
1231 case 1:
1232 return getOpDomainExplicitRhsPipeline<1>();
1233 case 2:
1234 return getOpDomainExplicitRhsPipeline<2>();
1235 case 3:
1236 break;
1237 default:
1238 THROW_MESSAGE("Not implemented");
1239 }
1240 return getOpDomainExplicitRhsPipeline<3>();
1241}

◆ getOpDomainLhsPipeline()

template<>
boost::ptr_deque< PipelineManager::UserDataOperator > & MoFEM::PipelineManager::getOpDomainLhsPipeline ( )
inline

Definition at line 1077 of file PipelineManager.hpp.

1085 {
1086 switch (cOre.getInterface<Simple>()->getDim()) {
1087 case 1:
1088 return getOpDomainLhsPipeline<1>();
1089 case 2:
1090 return getOpDomainLhsPipeline<2>();
1091 case 3:
1092 break;
1093 default:
1094 THROW_MESSAGE("Not implemented");
1095 }
1096 return getOpDomainLhsPipeline<3>();
1097}

◆ getOpDomainRhsPipeline()

template<>
boost::ptr_deque< PipelineManager::UserDataOperator > & MoFEM::PipelineManager::getOpDomainRhsPipeline ( )
inline

Definition at line 1101 of file PipelineManager.hpp.

1109 {
1110 switch (cOre.getInterface<Simple>()->getDim()) {
1111 case 1:
1112 return getOpDomainRhsPipeline<1>();
1113 case 2:
1114 return getOpDomainRhsPipeline<2>();
1115 case 3:
1116 break;
1117 default:
1118 THROW_MESSAGE("Not implemented");
1119 }
1120 return getOpDomainRhsPipeline<3>();
1121}

◆ getOpMeshsetExplicitRhsPipeline()

boost::ptr_deque< PipelineManager::UserDataOperator > & MoFEM::PipelineManager::getOpMeshsetExplicitRhsPipeline ( )

Get the Op Meshset Explicit Rhs Pipeline object.

Returns
boost::ptr_deque<UserDataOperator>&

Definition at line 47 of file PipelineManager.cpp.

47 {
48 return boost::dynamic_pointer_cast<ForcesAndSourcesCore>(
50 ->getOpPtrVector();
51}

◆ getOpMeshsetLhsPipeline()

boost::ptr_deque< PipelineManager::UserDataOperator > & MoFEM::PipelineManager::getOpMeshsetLhsPipeline ( )

Get the Op Meshset Lhs Pipeline object.

Returns
boost::ptr_deque<UserDataOperator>&

Definition at line 40 of file PipelineManager.cpp.

40 {
41 return boost::dynamic_pointer_cast<ForcesAndSourcesCore>(
43 ->getOpPtrVector();
44}

◆ getOpMeshsetRhsPipeline()

boost::ptr_deque< PipelineManager::UserDataOperator > & MoFEM::PipelineManager::getOpMeshsetRhsPipeline ( )

Get the Op Meshset Rhs Pipeline object.

Returns
boost::ptr_deque<UserDataOperator>&
Examples
mofem/atom_tests/simple_interface.cpp, and mofem/atom_tests/simple_l2_only.cpp.

Definition at line 33 of file PipelineManager.cpp.

33 {
34 return boost::dynamic_pointer_cast<ForcesAndSourcesCore>(
36 ->getOpPtrVector();
37}

◆ getOpSkeletonExplicitRhsPipeline()

template<>
boost::ptr_deque< PipelineManager::UserDataOperator > & MoFEM::PipelineManager::getOpSkeletonExplicitRhsPipeline ( )
inline

Definition at line 1269 of file PipelineManager.hpp.

1277 {
1278 switch (cOre.getInterface<Simple>()->getDim()) {
1279 case 1:
1280 return getOpSkeletonExplicitRhsPipeline<1>();
1281 case 2:
1282 return getOpSkeletonExplicitRhsPipeline<2>();
1283 case 3:
1284 break;
1285 default:
1286 THROW_MESSAGE("Not implemented");
1287 }
1288 return getOpSkeletonExplicitRhsPipeline<3>();
1289}

◆ getOpSkeletonLhsPipeline()

template<>
boost::ptr_deque< PipelineManager::UserDataOperator > & MoFEM::PipelineManager::getOpSkeletonLhsPipeline ( )
inline

Definition at line 1173 of file PipelineManager.hpp.

1181 {
1182 switch (cOre.getInterface<Simple>()->getDim()) {
1183 case 1:
1184 return getOpSkeletonLhsPipeline<1>();
1185 case 2:
1186 return getOpSkeletonLhsPipeline<2>();
1187 case 3:
1188 break;
1189 default:
1190 THROW_MESSAGE("Not implemented");
1191 }
1192 return getOpSkeletonLhsPipeline<3>();
1193}

◆ getOpSkeletonRhsPipeline()

template<>
boost::ptr_deque< PipelineManager::UserDataOperator > & MoFEM::PipelineManager::getOpSkeletonRhsPipeline ( )
inline

Definition at line 1197 of file PipelineManager.hpp.

1205 {
1206 switch (cOre.getInterface<Simple>()->getDim()) {
1207 case 1:
1208 return getOpSkeletonRhsPipeline<1>();
1209 case 2:
1210 return getOpSkeletonRhsPipeline<2>();
1211 case 3:
1212 break;
1213 default:
1214 THROW_MESSAGE("Not implemented");
1215 }
1216 return getOpSkeletonRhsPipeline<3>();
1217}

◆ getSkeletonExplicitRhsFE()

boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::getSkeletonExplicitRhsFE ( )
inline

Get skeleton explicit right-hand side finite element.

Returns
boost::shared_ptr<FEMethod>& Reference to skeleton explicit RHS finite element for IMEX schemes

Definition at line 736 of file PipelineManager.hpp.

736 {
738}

◆ getSkeletonLhsFE()

boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::getSkeletonLhsFE ( )
inline

Get skeleton left-hand side finite element.

Returns
boost::shared_ptr<FEMethod>& Reference to skeleton LHS finite element for matrix assembly

Definition at line 720 of file PipelineManager.hpp.

720 {
721 return feSkeletonLhs;
722}

◆ getSkeletonRhsFE()

boost::shared_ptr< FEMethod > & MoFEM::PipelineManager::getSkeletonRhsFE ( )
inline

Get skeleton right-hand side finite element.

Returns
boost::shared_ptr<FEMethod>& Reference to skeleton RHS finite element for vector assembly

Definition at line 724 of file PipelineManager.hpp.

724 {
725 return feSkeletonRhs;
726}

◆ query_interface()

MoFEMErrorCode MoFEM::PipelineManager::query_interface ( boost::typeindex::type_index  type_index,
UnknownInterface **  iface 
) const
virtual

Query interface for type-safe casting.

Parameters
type_indexType index of the requested interface
ifacePointer to interface pointer to be set
Returns
MoFEMErrorCode Error code (0 on success)

Implements MoFEM::UnknownInterface.

Definition at line 54 of file PipelineManager.cpp.

55 {
56 *iface = const_cast<PipelineManager *>(this);
57 return 0;
58}
PipelineManager(const MoFEM::Core &core)
Construct PipelineManager.

◆ setBoundaryExplicitRhsIntegrationRule() [1/2]

template<>
MoFEMErrorCode MoFEM::PipelineManager::setBoundaryExplicitRhsIntegrationRule ( PipelineManager::RuleHookFun  rule)
inline

Definition at line 1030 of file PipelineManager.hpp.

1032 {
1034 switch (cOre.getInterface<Simple>()->getDim()) {
1035 case 1:
1036 MoFEMFunctionReturnHot(setBoundaryExplicitRhsIntegrationRule<1>(rule));
1037 case 2:
1038 MoFEMFunctionReturnHot(setBoundaryExplicitRhsIntegrationRule<2>(rule));
1039 case 3:
1040 MoFEMFunctionReturnHot(setBoundaryExplicitRhsIntegrationRule<3>(rule));
1041 default:
1042 THROW_MESSAGE("Not implemented");
1043 }
1045}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()

◆ setBoundaryExplicitRhsIntegrationRule() [2/2]

template<int DIM>
MoFEMErrorCode MoFEM::PipelineManager::setBoundaryExplicitRhsIntegrationRule ( PipelineManager::RuleHookFun  rule)
inline

Set integration rule for boundary explicit right-hand side finite element.

Template Parameters
DIMDimension (-1 for automatic detection from problem)
Parameters
ruleIntegration rule hook function
Returns
MoFEMErrorCode Error code (0 on success)

Definition at line 1020 of file PipelineManager.hpp.

1021 {
1023 boost::dynamic_pointer_cast<ForcesAndSourcesCore>(
1024 createBoundaryFEPipeline<DIM>(feBoundaryExplicitRhs))
1025 ->getRuleHook = rule;
1027}

◆ setBoundaryLhsIntegrationRule() [1/2]

template<>
MoFEMErrorCode MoFEM::PipelineManager::setBoundaryLhsIntegrationRule ( PipelineManager::RuleHookFun  rule)
inline

Definition at line 895 of file PipelineManager.hpp.

896 {
898 switch (cOre.getInterface<Simple>()->getDim()) {
899 case 1:
900 MoFEMFunctionReturnHot(setBoundaryLhsIntegrationRule<1>(rule));
901 case 2:
902 MoFEMFunctionReturnHot(setBoundaryLhsIntegrationRule<2>(rule));
903 case 3:
904 MoFEMFunctionReturnHot(setBoundaryLhsIntegrationRule<3>(rule));
905 default:
906 THROW_MESSAGE("Not implemented");
907 }
909}

◆ setBoundaryLhsIntegrationRule() [2/2]

template<int DIM>
MoFEMErrorCode MoFEM::PipelineManager::setBoundaryLhsIntegrationRule ( PipelineManager::RuleHookFun  rule)
inline

Set integration rule for boundary left-hand side finite element.

Template Parameters
DIMDimension (-1 for automatic detection from problem)
Parameters
ruleIntegration rule hook function
Returns
MoFEMErrorCode Error code (0 on success)
Examples
mofem/tutorials/clx-0/helmholtz.cpp, and mofem/tutorials/scl-8/radiation.cpp.

Definition at line 885 of file PipelineManager.hpp.

886 {
888 boost::dynamic_pointer_cast<ForcesAndSourcesCore>(
889 createBoundaryFEPipeline<DIM>(feBoundaryLhs))
890 ->getRuleHook = rule;
892}

◆ setBoundaryRhsIntegrationRule() [1/2]

template<>
MoFEMErrorCode MoFEM::PipelineManager::setBoundaryRhsIntegrationRule ( PipelineManager::RuleHookFun  rule)
inline

Definition at line 922 of file PipelineManager.hpp.

923 {
925 switch (cOre.getInterface<Simple>()->getDim()) {
926 case 1:
927 MoFEMFunctionReturnHot(setBoundaryRhsIntegrationRule<1>(rule));
928 case 2:
929 MoFEMFunctionReturnHot(setBoundaryRhsIntegrationRule<2>(rule));
930 case 3:
931 MoFEMFunctionReturnHot(setBoundaryRhsIntegrationRule<3>(rule));
932 default:
933 THROW_MESSAGE("Not implemented");
934 }
936}

◆ setBoundaryRhsIntegrationRule() [2/2]

template<int DIM>
MoFEMErrorCode MoFEM::PipelineManager::setBoundaryRhsIntegrationRule ( PipelineManager::RuleHookFun  rule)
inline

Set integration rule for boundary right-hand side finite element.

Template Parameters
DIMDimension (-1 for automatic detection from problem)
Parameters
ruleIntegration rule hook function
Returns
MoFEMErrorCode Error code (0 on success)
Examples
mofem/atom_tests/child_and_parent.cpp, mofem/atom_tests/hanging_node_approx.cpp, mofem/tutorials/clx-0/helmholtz.cpp, and mofem/tutorials/scl-8/radiation.cpp.

Definition at line 912 of file PipelineManager.hpp.

913 {
915 boost::dynamic_pointer_cast<ForcesAndSourcesCore>(
916 createBoundaryFEPipeline<DIM>(feBoundaryRhs))
917 ->getRuleHook = rule;
919}

◆ setDomainExplicitRhsIntegrationRule() [1/2]

template<>
MoFEMErrorCode MoFEM::PipelineManager::setDomainExplicitRhsIntegrationRule ( PipelineManager::RuleHookFun  rule)
inline

Definition at line 1003 of file PipelineManager.hpp.

1004 {
1006 switch (cOre.getInterface<Simple>()->getDim()) {
1007 case 1:
1008 MoFEMFunctionReturnHot(setDomainExplicitRhsIntegrationRule<1>(rule));
1009 case 2:
1010 MoFEMFunctionReturnHot(setDomainExplicitRhsIntegrationRule<2>(rule));
1011 case 3:
1012 MoFEMFunctionReturnHot(setDomainExplicitRhsIntegrationRule<3>(rule));
1013 default:
1014 THROW_MESSAGE("Not implemented");
1015 }
1017}

◆ setDomainExplicitRhsIntegrationRule() [2/2]

template<int DIM>
MoFEMErrorCode MoFEM::PipelineManager::setDomainExplicitRhsIntegrationRule ( PipelineManager::RuleHookFun  rule)
inline

Set integration rule for domain explicit right-hand side finite element.

Template Parameters
DIMDimension (-1 for automatic detection from problem)
Parameters
ruleIntegration rule hook function
Returns
MoFEMErrorCode Error code (0 on success)

Definition at line 993 of file PipelineManager.hpp.

994 {
996 boost::dynamic_pointer_cast<ForcesAndSourcesCore>(
997 createDomainFEPipeline<DIM>(feDomainExplicitRhs))
998 ->getRuleHook = rule;
1000}

◆ setDomainLhsIntegrationRule() [1/2]

template<>
MoFEMErrorCode MoFEM::PipelineManager::setDomainLhsIntegrationRule ( PipelineManager::RuleHookFun  rule)
inline

Definition at line 841 of file PipelineManager.hpp.

842 {
844 switch (cOre.getInterface<Simple>()->getDim()) {
845 case 1:
846 MoFEMFunctionReturnHot(setDomainLhsIntegrationRule<1>(rule));
847 case 2:
848 MoFEMFunctionReturnHot(setDomainLhsIntegrationRule<2>(rule));
849 case 3:
850 MoFEMFunctionReturnHot(setDomainLhsIntegrationRule<3>(rule));
851 default:
852 THROW_MESSAGE("Not implemented");
853 }
855}
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...

◆ setDomainLhsIntegrationRule() [2/2]

template<int DIM>
MoFEMErrorCode MoFEM::PipelineManager::setDomainLhsIntegrationRule ( PipelineManager::RuleHookFun  rule)
inline

Set integration rule for domain left-hand side finite element.

Template Parameters
DIMDimension (-1 for automatic detection from problem)
Parameters
ruleIntegration rule hook function
Returns
MoFEMErrorCode Error code (0 on success)
Examples
mofem/atom_tests/child_and_parent.cpp, mofem/atom_tests/dg_projection.cpp, mofem/atom_tests/hanging_node_approx.cpp, mofem/atom_tests/higher_derivatives.cpp, mofem/atom_tests/schur_test_diag_mat.cpp, mofem/tutorials/clx-0/helmholtz.cpp, mofem/tutorials/mix-0/mixed_poisson.cpp, mofem/tutorials/scl-0/approximaton.cpp, mofem/tutorials/scl-1/poisson_2d_homogeneous.cpp, mofem/tutorials/scl-11/poisson_2d_dis_galerkin.cpp, mofem/tutorials/scl-8/radiation.cpp, mofem/tutorials/scl-9/heat_method.cpp, mofem/tutorials/vec-4/approx_sphere.cpp, and mofem/tutorials/vec-6/plate.cpp.

Definition at line 831 of file PipelineManager.hpp.

832 {
834 boost::dynamic_pointer_cast<ForcesAndSourcesCore>(
835 createDomainFEPipeline<DIM>(feDomainLhs))
836 ->getRuleHook = rule;
838}

◆ setDomainRhsIntegrationRule() [1/2]

template<>
MoFEMErrorCode MoFEM::PipelineManager::setDomainRhsIntegrationRule ( PipelineManager::RuleHookFun  rule)
inline

Definition at line 868 of file PipelineManager.hpp.

869 {
871 switch (cOre.getInterface<Simple>()->getDim()) {
872 case 1:
873 MoFEMFunctionReturnHot(setDomainRhsIntegrationRule<1>(rule));
874 case 2:
875 MoFEMFunctionReturnHot(setDomainRhsIntegrationRule<2>(rule));
876 case 3:
877 MoFEMFunctionReturnHot(setDomainRhsIntegrationRule<3>(rule));
878 default:
879 THROW_MESSAGE("Not implemented");
880 }
882}

◆ setDomainRhsIntegrationRule() [2/2]

template<int DIM>
MoFEMErrorCode MoFEM::PipelineManager::setDomainRhsIntegrationRule ( PipelineManager::RuleHookFun  rule)
inline

◆ setSkeletonExplicitRhsIntegrationRule() [1/2]

template<>
MoFEMErrorCode MoFEM::PipelineManager::setSkeletonExplicitRhsIntegrationRule ( PipelineManager::RuleHookFun  rule)
inline

Definition at line 1058 of file PipelineManager.hpp.

1060 {
1062 switch (cOre.getInterface<Simple>()->getDim()) {
1063 case 1:
1064 MoFEMFunctionReturnHot(setSkeletonExplicitRhsIntegrationRule<1>(rule));
1065 case 2:
1066 MoFEMFunctionReturnHot(setSkeletonExplicitRhsIntegrationRule<2>(rule));
1067 case 3:
1068 MoFEMFunctionReturnHot(setSkeletonExplicitRhsIntegrationRule<3>(rule));
1069 default:
1070 THROW_MESSAGE("Not implemented");
1071 }
1073}

◆ setSkeletonExplicitRhsIntegrationRule() [2/2]

template<int DIM>
MoFEMErrorCode MoFEM::PipelineManager::setSkeletonExplicitRhsIntegrationRule ( PipelineManager::RuleHookFun  rule)
inline

Set integration rule for skeleton explicit right-hand side finite element.

Template Parameters
DIMDimension (-1 for automatic detection from problem)
Parameters
ruleIntegration rule hook function
Returns
MoFEMErrorCode Error code (0 on success)

Definition at line 1048 of file PipelineManager.hpp.

1049 {
1051 boost::dynamic_pointer_cast<ForcesAndSourcesCore>(
1052 createBoundaryFEPipeline<DIM>(feSkeletonExplicitRhs))
1053 ->getRuleHook = rule;
1055}

◆ setSkeletonLhsIntegrationRule() [1/2]

template<>
MoFEMErrorCode MoFEM::PipelineManager::setSkeletonLhsIntegrationRule ( PipelineManager::RuleHookFun  rule)
inline

Definition at line 949 of file PipelineManager.hpp.

950 {
952 switch (cOre.getInterface<Simple>()->getDim()) {
953 case 1:
954 MoFEMFunctionReturnHot(setSkeletonLhsIntegrationRule<1>(rule));
955 case 2:
956 MoFEMFunctionReturnHot(setSkeletonLhsIntegrationRule<2>(rule));
957 case 3:
958 MoFEMFunctionReturnHot(setSkeletonLhsIntegrationRule<3>(rule));
959 default:
960 THROW_MESSAGE("Not implemented");
961 }
963}

◆ setSkeletonLhsIntegrationRule() [2/2]

template<int DIM>
MoFEMErrorCode MoFEM::PipelineManager::setSkeletonLhsIntegrationRule ( PipelineManager::RuleHookFun  rule)
inline

Set integration rule for skeleton left-hand side finite element.

Template Parameters
DIMDimension (-1 for automatic detection from problem)
Parameters
ruleIntegration rule hook function
Returns
MoFEMErrorCode Error code (0 on success)

Definition at line 939 of file PipelineManager.hpp.

940 {
942 boost::dynamic_pointer_cast<ForcesAndSourcesCore>(
943 createBoundaryFEPipeline<DIM>(feSkeletonLhs))
944 ->getRuleHook = rule;
946}

◆ setSkeletonRhsIntegrationRule() [1/2]

template<>
MoFEMErrorCode MoFEM::PipelineManager::setSkeletonRhsIntegrationRule ( PipelineManager::RuleHookFun  rule)
inline

Definition at line 976 of file PipelineManager.hpp.

977 {
979 switch (cOre.getInterface<Simple>()->getDim()) {
980 case 1:
981 MoFEMFunctionReturnHot(setSkeletonRhsIntegrationRule<1>(rule));
982 case 2:
983 MoFEMFunctionReturnHot(setSkeletonRhsIntegrationRule<2>(rule));
984 case 3:
985 MoFEMFunctionReturnHot(setSkeletonRhsIntegrationRule<3>(rule));
986 default:
987 THROW_MESSAGE("Not implemented");
988 }
990}

◆ setSkeletonRhsIntegrationRule() [2/2]

template<int DIM>
MoFEMErrorCode MoFEM::PipelineManager::setSkeletonRhsIntegrationRule ( PipelineManager::RuleHookFun  rule)
inline

Set integration rule for skeleton right-hand side finite element.

Template Parameters
DIMDimension (-1 for automatic detection from problem)
Parameters
ruleIntegration rule hook function
Returns
MoFEMErrorCode Error code (0 on success)

Definition at line 966 of file PipelineManager.hpp.

967 {
969 boost::dynamic_pointer_cast<ForcesAndSourcesCore>(
970 createBoundaryFEPipeline<DIM>(feSkeletonRhs))
971 ->getRuleHook = rule;
973}

Member Data Documentation

◆ cOre

MoFEM::Core& MoFEM::PipelineManager::cOre
private

Reference to MoFEM core instance.

Definition at line 555 of file PipelineManager.hpp.

◆ feBoundaryExplicitRhs

boost::shared_ptr<FEMethod> MoFEM::PipelineManager::feBoundaryExplicitRhs
private

Element to assemble explicit Rhs for IMEX solver

Definition at line 573 of file PipelineManager.hpp.

◆ feBoundaryLhs

boost::shared_ptr<FEMethod> MoFEM::PipelineManager::feBoundaryLhs
private

Element to assemble LHS side by integrating boundary.

Definition at line 564 of file PipelineManager.hpp.

◆ feBoundaryRhs

boost::shared_ptr<FEMethod> MoFEM::PipelineManager::feBoundaryRhs
private

Element to assemble RHS side by integrating boundary.

Definition at line 562 of file PipelineManager.hpp.

◆ feDomainExplicitRhs

boost::shared_ptr<FEMethod> MoFEM::PipelineManager::feDomainExplicitRhs
private

Element to assemble explicit Rhs for IMEX solver.

Definition at line 571 of file PipelineManager.hpp.

◆ feDomainLhs

boost::shared_ptr<FEMethod> MoFEM::PipelineManager::feDomainLhs
private

Element to assemble LHS side by integrating domain.

Definition at line 560 of file PipelineManager.hpp.

◆ feDomainRhs

boost::shared_ptr<FEMethod> MoFEM::PipelineManager::feDomainRhs
private

Element to assemble RHS side by integrating domain.

Definition at line 558 of file PipelineManager.hpp.

◆ feMeshsetExplicitRhs

boost::shared_ptr<FEMethod> MoFEM::PipelineManager::feMeshsetExplicitRhs
private

Element to assemble explicit RHS side by integrating meshset

Definition at line 583 of file PipelineManager.hpp.

◆ feMeshsetLhs

boost::shared_ptr<FEMethod> MoFEM::PipelineManager::feMeshsetLhs
private

Element to assemble LHS side by integrating meshset

Definition at line 581 of file PipelineManager.hpp.

◆ feMeshsetRhs

boost::shared_ptr<FEMethod> MoFEM::PipelineManager::feMeshsetRhs
private

Element to assemble RHS side by integrating meshset

Definition at line 579 of file PipelineManager.hpp.

◆ feSkeletonExplicitRhs

boost::shared_ptr<FEMethod> MoFEM::PipelineManager::feSkeletonExplicitRhs
private

Element to assemble explicit Rhs for IMEX solver

Definition at line 576 of file PipelineManager.hpp.

◆ feSkeletonLhs

boost::shared_ptr<FEMethod> MoFEM::PipelineManager::feSkeletonLhs
private

Element to assemble LHS side by integrating skeleton.

Definition at line 568 of file PipelineManager.hpp.

◆ feSkeletonRhs

boost::shared_ptr<FEMethod> MoFEM::PipelineManager::feSkeletonRhs
private

Element to assemble RHS side by integrating skeleton.

Definition at line 566 of file PipelineManager.hpp.


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