v0.14.0
Classes | Public Member Functions | Static Public Member Functions | Private Types | Private Attributes | List of all members
MoFEM::UnknownInterface Struct Referenceabstract

base class for all interface classes More...

#include <src/interfaces/UnknownInterface.hpp>

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

Classes

struct  NotKnownClass
 
struct  UIdTypeMap
 

Public Member Functions

virtual MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const =0
 
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface reference to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface. More...
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface. More...
 
virtual ~UnknownInterface ()=default
 

Static Public Member Functions

static MoFEMErrorCode getLibVersion (Version &version)
 Get library version. More...
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version. More...
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version. More...
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version. More...
 

Private Types

typedef multi_index_container< UIdTypeMap, indexed_by< hashed_unique< member< UIdTypeMap, boost::typeindex::type_index, &UIdTypeMap::classIdx > > > > iFaceTypeMap_multiIndex
 Data structure for interfaces Id and class names. More...
 

Private Attributes

iFaceTypeMap_multiIndex iFaceTypeMap
 Maps implementation to interface type name. More...
 

Detailed Description

base class for all interface classes

Examples
EshelbianPlasticity.cpp, and forces_and_sources_testing_users_base.cpp.

Definition at line 34 of file UnknownInterface.hpp.

Member Typedef Documentation

◆ iFaceTypeMap_multiIndex

typedef multi_index_container< UIdTypeMap, indexed_by< hashed_unique<member<UIdTypeMap, boost::typeindex::type_index, &UIdTypeMap::classIdx> > > > MoFEM::UnknownInterface::iFaceTypeMap_multiIndex
private

Data structure for interfaces Id and class names.

Definition at line 297 of file UnknownInterface.hpp.

Constructor & Destructor Documentation

◆ ~UnknownInterface()

virtual MoFEM::UnknownInterface::~UnknownInterface ( )
virtualdefault

Member Function Documentation

◆ getFileVersion()

MoFEMErrorCode MoFEM::UnknownInterface::getFileVersion ( moab::Interface &  moab,
Version version 
)
static

Get database major version.

This is database version. MoFEM can read DataBase from file created by older version. Then library version and database version could be different.

Returns
error code
Examples
mortar_contact_thermal.cpp, and simple_contact_thermal.cpp.

Definition at line 16 of file UnknownInterface.cpp.

17  {
19  const EntityHandle root_meshset = 0;
20  const int def_version[] = {-1, -1, -1};
21  Tag th;
22  rval = moab.tag_get_handle("MOFEM_VERSION", 3, MB_TYPE_INTEGER, th,
23  MB_TAG_CREAT | MB_TAG_MESH, &def_version);
24  int *version_ptr;
25  if (rval == MB_ALREADY_ALLOCATED) {
26  const void *tag_data[1];
27  CHKERR moab.tag_get_by_ptr(th, &root_meshset, 1, tag_data);
28  version_ptr = (int *)tag_data[0];
29  } else {
30  const void *tag_data[1];
31  CHKERR moab.tag_get_by_ptr(th, &root_meshset, 1, tag_data);
32  version_ptr = (int *)tag_data[0];
33  version_ptr[0] = MoFEM_VERSION_MAJOR;
34  version_ptr[1] = MoFEM_VERSION_MINOR;
35  version_ptr[2] = MoFEM_VERSION_BUILD;
36  }
37  version = Version(version_ptr);
39 }

◆ getInterface() [1/5]

template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE MoFEM::UnknownInterface::getInterface ( ) const
inline

Get interface pointer to pointer of interface.

// Create moab database
moab::Core mb_instance;
// Access moab database by interface
moab::Interface &moab = mb_instance;
// Create MoFEM database
MoFEM::Core core(moab);
// Acces MoFEM database by Interface
MoFEM::Interface &m_field = core;
// Get interface
// Get simple interface is simplified version enabling quick and
// easy construction of problem.
Simple *simple_interface = m_field.getInterface<Simple*>();
Returns
IFACE*

Definition at line 162 of file UnknownInterface.hpp.

162  {
163  typedef typename boost::remove_pointer<IFACE>::type IFaceType;
164  IFaceType *iface = NULL;
165  ierr = getInterface<IFACE>(iface);
166  CHKERRABORT(PETSC_COMM_SELF, ierr);
167  return iface;
168  }

◆ getInterface() [2/5]

template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE MoFEM::UnknownInterface::getInterface ( ) const
inline

Get reference to interface.

// Create moab database
moab::Core mb_instance;
// Access moab database by interface
moab::Interface &moab = mb_instance;
// Create MoFEM database
MoFEM::Core core(moab);
// Acces MoFEM database by Interface
MoFEM::Interface &m_field = core;
// Get interface
// Get simple interface is simplified version enabling quick and
// easy construction of problem.
Simple &simple_interface = m_field.getInterface<Simple&>();
Returns
IFACE&

Definition at line 195 of file UnknownInterface.hpp.

195  {
196  typedef typename boost::remove_reference<IFACE>::type IFaceType;
197  IFaceType *iface = NULL;
198  ierr = getInterface<IFaceType>(iface);
199  CHKERRABORT(PETSC_COMM_SELF, ierr);
200  return *iface;
201  }

◆ getInterface() [3/5]

template<class IFACE >
IFACE* MoFEM::UnknownInterface::getInterface ( ) const
inline

Function returning pointer to interface.

// Create moab database
moab::Core mb_instance;
// Access moab database by interface
moab::Interface &moab = mb_instance;
// Create MoFEM database
MoFEM::Core core(moab);
// Acces MoFEM database by Interface
MoFEM::Interface &m_field = core;
// Get interface
// Get simple interface is simplified version enabling quick and
// easy construction of problem.
Simple *simple_interface = m_field.getInterface<Simple>();
Returns
IFACE*

Definition at line 226 of file UnknownInterface.hpp.

226  {
227  IFACE *iface = NULL;
228  ierr = getInterface<IFACE>(iface);
229  CHKERRABORT(PETSC_COMM_SELF, ierr);
230  return iface;
231  }

◆ getInterface() [4/5]

template<class IFACE >
MoFEMErrorCode MoFEM::UnknownInterface::getInterface ( IFACE *&  iface) const
inline

Get interface reference to pointer of interface.

// Create moab database
moab::Core mb_instance;
// Access moab database by interface
moab::Interface &moab = mb_instance;
// Create MoFEM database
MoFEM::Core core(moab);
// Acces MoFEM database by Interface
MoFEM::Interface &m_field = core;
// Get interface
// Get simple interface is simplified version enabling quick and
// easy construction of problem.
Simple *simple_interface;
// Query interface and get pointer to Simple interface
CHKERR m_field.getInterface(simple_interface);
Parameters
ifacereference to a interface pointer
Returns
MoFEMErrorCode
Examples
add_blockset.cpp, add_cubit_meshsets.cpp, adolc_plasticity.cpp, ADOLCPlasticityMaterialModels.hpp, analytical_nonlinear_poisson.cpp, analytical_poisson.cpp, analytical_poisson_field_split.cpp, approx_sphere.cpp, bernstein_bezier_generate_base.cpp, bone_adaptation.cpp, boundary_marker.cpp, build_large_problem.cpp, build_problems.cpp, child_and_parent.cpp, continuity_check_on_contact_prism_side_ele.cpp, continuity_check_on_skeleton_3d.cpp, continuity_check_on_skeleton_with_simple_2d_for_h1.cpp, continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp, continuity_check_on_skeleton_with_simple_2d_for_hdiv.cpp, cubit_bc_test.cpp, dg_projection.cpp, dm_build_partitioned_mesh.cpp, dm_create_subdm.cpp, dm_partitioned_no_field.cpp, dynamic_first_order_con_law.cpp, eigen_elastic.cpp, elasticity.cpp, elasticity_mixed_formulation.cpp, ep.cpp, EshelbianPlasticity.cpp, field_blas.cpp, field_evaluator.cpp, field_to_vertices.cpp, forces_and_sources_testing_edge_element.cpp, forces_and_sources_testing_flat_prism_element.cpp, forces_and_sources_testing_users_base.cpp, free_surface.cpp, hanging_node_approx.cpp, hcurl_check_approx_in_2d.cpp, hcurl_curl_operator.cpp, hcurl_divergence_operator_2d.cpp, hdiv_check_approx_in_3d.cpp, hdiv_divergence_operator.cpp, heat_equation.cpp, heat_method.cpp, hello_world.cpp, helmholtz.cpp, HenckyOps.hpp, higher_derivatives.cpp, level_set.cpp, loop_entities.cpp, lorentz_force.cpp, MagneticElement.hpp, mesh_cut.cpp, mesh_insert_interface_atom.cpp, mesh_smoothing.cpp, meshset_to_vtk.cpp, minimal_surface_area.cpp, mixed_poisson.cpp, mortar_contact.cpp, mortar_contact_thermal.cpp, node_merge.cpp, nonlinear_dynamics.cpp, nonlinear_elastic.cpp, NonlinearElasticElementInterface.hpp, operators_tests.cpp, partition_mesh.cpp, phase.cpp, photon_diffusion.cpp, plastic.cpp, PlasticOps.hpp, plot_base.cpp, poisson_2d_dis_galerkin.cpp, poisson_2d_homogeneous.cpp, prism_elements_from_surface.cpp, prism_polynomial_approximation.cpp, quad_polynomial_approximation.cpp, reaction_diffusion.cpp, remove_entities_from_problem.cpp, remove_entities_from_problem_not_partitioned.cpp, scalar_check_approximation.cpp, schur_test_diag_mat.cpp, seepage.cpp, shallow_wave.cpp, simple_contact.cpp, simple_contact_thermal.cpp, simple_elasticity.cpp, simple_interface.cpp, simple_l2_only.cpp, split_sideset.cpp, tensor_divergence_operator.cpp, test_broken_space.cpp, test_cache_on_entities.cpp, test_jacobian_of_simple_contact_element.cpp, testing_jacobian_of_hook_element.cpp, testing_jacobian_of_hook_scaled_with_density_element.cpp, thermo_elastic.cpp, ThermoElasticOps.hpp, unsaturated_transport.cpp, UnsaturatedFlow.hpp, and wave_equation.cpp.

Definition at line 93 of file UnknownInterface.hpp.

93  {
95  UnknownInterface *ptr;
96  CHKERR query_interface(boost::typeindex::type_id<IFACE>(), &ptr);
97  if (!(iface = dynamic_cast<IFACE *>(ptr)))
98  SETERRQ2(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE,
99  "Cast impossible %s -> %s",
100  (boost::typeindex::type_id_runtime(*ptr).pretty_name()).c_str(),
101  (boost::typeindex::type_id<IFACE>().pretty_name()).c_str());
103  }

◆ getInterface() [5/5]

template<class IFACE >
MoFEMErrorCode MoFEM::UnknownInterface::getInterface ( IFACE **const  iface) const
inline

Get interface pointer to pointer of interface.

// Create moab database
moab::Core mb_instance;
// Access moab database by interface
moab::Interface &moab = mb_instance;
// Create MoFEM database
MoFEM::Core core(moab);
// Acces MoFEM database by Interface
MoFEM::Interface &m_field = core;
// Get interface
// Get simple interface is simplified version enabling quick and
// easy construction of problem.
Simple *simple_interface;
// Query interface and get pointer to Simple interface
CHKERR m_field.getInterface(&simple_interface);
Parameters
ifaceconst pointer to a interface pointer
Returns
MoFEMErrorCode

Definition at line 133 of file UnknownInterface.hpp.

133  {
134  return getInterface<IFACE>(*iface);
135  }

◆ getInterfaceVersion()

MoFEMErrorCode MoFEM::UnknownInterface::getInterfaceVersion ( Version version)
static

Get database major version.

Implementation of particular interface could be different than main lib. For example user could use older interface, to keep back compatibility.

Returns
error code

Definition at line 59 of file UnknownInterface.cpp.

59  {
61  version =
62  Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD);
64 }

◆ getLibVersion()

MoFEMErrorCode MoFEM::UnknownInterface::getLibVersion ( Version version)
static

Get library version.

This is library version.

Returns
error code

Definition at line 9 of file UnknownInterface.cpp.

9  {
11  version =
12  Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD);
14 }

◆ query_interface()

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

Implemented in MoFEM::BaseFunctionUnknownInterface, MoFEM::DMCtx, MoFEM::DofMethod, MoFEM::EntityMethod, MoFEM::FEMethod, MoFEM::CoreTmp< 0 >, MoFEM::BasicMethod, MoFEM::TSMethod, MoFEM::MeshsetsManager, MoFEM::SnesMethod, MoFEM::IntegratedJacobiPolynomial, MoFEM::KernelLobattoPolynomial, MoFEM::KspMethod, MoFEM::IntegratedJacobiPolynomialCtx, MoFEM::FatPrismPolynomialBase, MoFEM::JacobiPolynomial, MoFEM::KernelLobattoPolynomialCtx, MoFEM::LogManager, SomeUserPolynomialBase, MoFEM::LegendrePolynomial, MoFEM::FlatPrismPolynomialBase, MoFEM::LobattoPolynomial, MoFEM::BaseFunctionCtx, MoFEM::Simple, MoFEM::MeshRefinement, MoFEM::BcManager, MoFEM::SeriesRecorder, MoFEM::PipelineManager, MoFEM::FatPrismPolynomialBaseCtx, MoFEM::ISManager, MoFEM::MedInterface, MoFEM::PrismInterface, MoFEM::TetGenInterface, MoFEM::VecManager, MoFEM::EntPolynomialBaseCtx, MoFEM::FlatPrismPolynomialBaseCtx, MoFEM::QuadPolynomialBase, MoFEM::BitRefManager, MoFEM::CommInterface, MoFEM::FieldBlas, MoFEM::FieldEvaluatorInterface, MoFEM::MatrixManager, MoFEM::OperatorsTester, MoFEM::ProblemsManager, MoFEM::EdgePolynomialBase, MoFEM::HexPolynomialBase, MoFEM::CutMeshInterface, MoFEM::Tools, MoFEM::NodeMergerInterface, MoFEM::JacobiPolynomialCtx, MoFEM::LegendrePolynomialCtx, MoFEM::LobattoPolynomialCtx, MoFEM::TetPolynomialBase, MoFEM::TriPolynomialBase, MoFEM::PrismsFromSurfaceInterface, MoFEM::PetscData, MoFEM::BaseFunction, MoFEM::DMMGViaApproxOrdersCtx, FractureMechanics::CPSolvers, and FractureMechanics::CPMeshCut.

Examples
cell_forces.cpp.

◆ registerInterface()

template<class IFACE >
MoFEMErrorCode MoFEM::UnknownInterface::registerInterface ( bool  error_if_registration_failed = true)
inline

Register interface.

Example:

ierr = regSubInterface<Simple>(IDD_MOFEMSimple);
CHKERRABORT(PETSC_COMM_SELF, ierr);
Parameters
uuid
true
Returns
MoFEMErrorCode

Definition at line 54 of file UnknownInterface.hpp.

54  {
56  auto p =
57  iFaceTypeMap.insert(UIdTypeMap(boost::typeindex::type_id<IFACE>()));
58  if (error_if_registration_failed && (!p.second)) {
59  SETERRQ1(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
60  "Registration of interface typeid(IFACE).name() = %s failed",
61  typeid(IFACE).name());
62  }
64  }

◆ setFileVersion()

MoFEMErrorCode MoFEM::UnknownInterface::setFileVersion ( moab::Interface &  moab,
Version  version = Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR,                                MoFEM_VERSION_BUILD) 
)
static

Get database major version.

This is database version. MoFEM can read DataBase from file created by older version. Then library version and database version could be different.

Returns
error code

Definition at line 41 of file UnknownInterface.cpp.

42  {
44  const EntityHandle root_meshset = 0;
45  const int def_version[] = {-1, -1, -1};
46  Tag th;
47  rval = moab.tag_get_handle("MOFEM_VERSION", 3, MB_TYPE_INTEGER, th,
48  MB_TAG_CREAT | MB_TAG_MESH, &def_version);
49  int *version_ptr;
50  const void *tag_data[1];
51  CHKERR moab.tag_get_by_ptr(th, &root_meshset, 1, tag_data);
52  version_ptr = (int *)tag_data[0];
53  version_ptr[0] = version.majorVersion;
54  version_ptr[1] = version.minorVersion;
55  version_ptr[2] = version.buildVersion;
57 }

Member Data Documentation

◆ iFaceTypeMap

iFaceTypeMap_multiIndex MoFEM::UnknownInterface::iFaceTypeMap
mutableprivate

Maps implementation to interface type name.

Definition at line 300 of file UnknownInterface.hpp.


The documentation for this struct was generated from the following files:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
EntityHandle
MoFEM::UnknownInterface::iFaceTypeMap
iFaceTypeMap_multiIndex iFaceTypeMap
Maps implementation to interface type name.
Definition: UnknownInterface.hpp:300
MoFEM::UnknownInterface::query_interface
virtual MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const =0
MoFEM::th
Tag th
Definition: Projection10NodeCoordsOnField.cpp:122
MOFEM_IMPOSSIBLE_CASE
@ MOFEM_IMPOSSIBLE_CASE
Definition: definitions.h:35
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
MoFEM::Exceptions::rval
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
Definition: Exceptions.hpp:74
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
convert.type
type
Definition: convert.py:64
MOFEM_OPERATION_UNSUCCESSFUL
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition: definitions.h:34
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359