v0.14.0
Macros | Enumerations | Functions | Variables
definitions.h File Reference

useful compiler directives and definitions More...

Go to the source code of this file.

Macros

#define DEPRECATED
 
#define MYPCOMM_INDEX   0
 default communicator number PCOMM More...
 
#define MAX_CORE_TMP   1
 maximal number of cores More...
 
#define BITREFEDGES_SIZE   32
 number refined edges More...
 
#define BITREFLEVEL_SIZE   64
 max number of refinements More...
 
#define BITFIELDID_SIZE   32
 max number of fields More...
 
#define BITFEID_SIZE   32
 max number of finite elements More...
 
#define BITPROBLEMID_SIZE   32
 max number of problems More...
 
#define BITINTERFACEUID_SIZE   32
 
#define MB_TYPE_WIDTH   4
 
#define MB_ID_WIDTH   (8 * sizeof(EntityHandle) - MB_TYPE_WIDTH)
 
#define MB_TYPE_MASK   ((EntityHandle)0xF << MB_ID_WIDTH)
 
#define MB_START_ID   ((EntityID)1)
 All entity id's currently start at 1. More...
 
#define MB_END_ID   ((EntityID)MB_ID_MASK)
 Last id is the complement of the MASK. More...
 
#define MB_ID_MASK   (~MB_TYPE_MASK)
 
#define MAX_DOFS_ON_ENTITY   512
 Maximal number of DOFs on entity. More...
 
#define MAX_PROCESSORS_NUMBER   1024
 Maximal number of processors. More...
 
#define DOF_UID_MASK   (MAX_DOFS_ON_ENTITY - 1)
 Mask for DOF number on entity form UId. More...
 
#define ENTITY_UID_MASK   (~DOF_UID_MASK)
 
#define NOT_USED(x)   ((void)(x))
 
#define BARRIER_PCOMM_RANK_START(PCMB)
 set barrier start Run code in sequence, starting from process 0, and ends on last process. More...
 
#define BARRIER_RANK_START(PCMB)
 
#define BARRIER_PCOMM_RANK_END(PCMB)
 set barrier start Run code in sequence, starting from process 0, and ends on last process. More...
 
#define BARRIER_RANK_END(PCMB)
 
#define BARRIER_MOFEM_RANK_START(MOFEM)
 set barrier start Run code in sequence, starting from process 0, and ends on last process. More...
 
#define BARRIER_MOFEM_RANK_END(MOFEM)
 set barrier start Run code in sequence, starting from process 0, and ends on last process. More...
 
#define MoFEMFunctionBegin
 First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions should be MoFEMFunctionReturn(0);. More...
 
#define CATCH_ERRORS
 Catch errors. More...
 
#define MoFEMFunctionReturn(a)
 Last executable line of each PETSc function used for error handling. Replaces return() More...
 
#define MoFEMFunctionBeginHot   PetscFunctionBeginHot
 First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions should be MoFEMFunctionReturn(0); Use of this function allows for lighter profiling by default. More...
 
#define MoFEMFunctionReturnHot(a)   PetscFunctionReturn(a)
 Last executable line of each PETSc function used for error handling. Replaces return() More...
 
#define CHKERRQ_PETSC(n)   CHKERRQ(n)
 
#define CHKERRQ_MOAB(a)
 check error code of MoAB function More...
 
#define CHKERRG(n)
 Check error code of MoFEM/MOAB/PETSc function. More...
 
#define CHKERR   MoFEM::ErrorChecker<__LINE__>() <<
 Inline error check. More...
 
#define MOAB_THROW(err)
 Check error code of MoAB function and throw MoFEM exception. More...
 
#define THROW_MESSAGE(msg)
 Throw MoFEM exception. More...
 
#define CHK_MOAB_THROW(err, msg)
 Check error code of MoAB function and throw MoFEM exception. More...
 
#define CHK_THROW_MESSAGE(err, msg)
 Check and throw MoFEM exception. More...
 
#define SSTR(x)   toString(x)
 Convert number to string. More...
 
#define TENSOR1_VEC_PTR(VEC)   &VEC[0], &VEC[1], &VEC[2]
 
#define SYMMETRIC_TENSOR4_MAT_PTR(MAT)
 
#define TENSOR4_MAT_PTR(MAT)   &MAT(0, 0), MAT.size2()
 
#define TENSOR2_MAT_PTR(MAT)
 
#define SYMMETRIC_TENSOR2_MAT_PTR(MAT)   &MAT(0, 0), &MAT(0, 1), &MAT(0, 2), &MAT(0, 3), &MAT(0, 4), &MAT(0, 5)
 
#define SYMMETRIC_TENSOR2_VEC_PTR(VEC)   &VEC[0], &VEC[1], &VEC[2], &VEC[3], &VEC[4], &VEC[5]
 

Enumerations

enum  MoFEMErrorCodes {
  MOFEM_SUCCESS = 0, MOFEM_DATA_INCONSISTENCY = 100, MOFEM_NOT_IMPLEMENTED = 101, MOFEM_NOT_FOUND = 102,
  MOFEM_OPERATION_UNSUCCESSFUL = 103, MOFEM_IMPOSSIBLE_CASE = 104, MOFEM_INVALID_DATA = 105, MOFEM_NOT_INSTALLED = 106,
  MOFEM_MOFEMEXCEPTION_THROW = 107, MOFEM_STD_EXCEPTION_THROW = 108, MOFEM_ATOM_TEST_INVALID = 109, MOFEM_MOAB_ERROR = 110
}
 Error handling. More...
 
enum  FieldApproximationBase {
  NOBASE = 0, AINSWORTH_LEGENDRE_BASE, AINSWORTH_LOBATTO_BASE, AINSWORTH_BERNSTEIN_BEZIER_BASE,
  DEMKOWICZ_JACOBI_BASE, USER_BASE, LASTBASE
}
 approximation base More...
 
enum  FieldSpace {
  NOSPACE = 0, NOFIELD = 1, H1, HCURL,
  HDIV, L2, LASTSPACE
}
 approximation spaces More...
 
enum  FieldContinuity { CONTINUOUS = 0, DISCONTINUOUS = 1, LASTCONTINUITY }
 Field continuity. More...
 
enum  MoFEMTypes { MF_ZERO = 0, MF_EXCL = 1 << 0, MF_EXIST = 1 << 1, MF_NOT_THROW = 1 << 2 }
 Those types control how functions respond on arguments, f.e. error handling. More...
 
enum  CoordinateTypes {
  CARTESIAN, POLAR, CYLINDRICAL, SPHERICAL,
  LAST_COORDINATE_SYSTEM
}
 Coodinate system. More...
 
enum  RowColData { ROW = 0, COL, DATA, LASTROWCOLDATA }
 RowColData. More...
 
enum  ByWhat {
  BYROW = 1 << 0, BYCOL = 1 << 1, BYDATA = 1 << 2, BYROWDATA = 1 << 0 | 1 << 2,
  BYCOLDATA = 1 << 1 | 1 << 2, BYROWCOL = 1 << 0 | 1 << 1, BYALL = 1 << 0 | 1 << 1 | 1 << 2
}
 
enum  CubitBC {
  UNKNOWNSET = 0, NODESET = 1 << 0, SIDESET = 1 << 1, BLOCKSET = 1 << 2,
  MATERIALSET = 1 << 3, DISPLACEMENTSET = 1 << 4, FORCESET = 1 << 5, PRESSURESET = 1 << 6,
  VELOCITYSET = 1 << 7, ACCELERATIONSET = 1 << 8, TEMPERATURESET = 1 << 9, HEATFLUXSET = 1 << 10,
  INTERFACESET = 1 << 11, UNKNOWNNAME = 1 << 12, MAT_ELASTICSET = 1 << 13, MAT_INTERFSET = 1 << 14,
  MAT_THERMALSET = 1 << 15, BODYFORCESSET = 1 << 16, MAT_MOISTURESET = 1 << 17, DIRICHLET_BC = 1 << 18,
  NEUMANN_BC = 1 << 19, LASTSET_BC = 1 << 20
}
 Types of sets and boundary conditions. More...
 
enum  HVecFormatting { HVEC0 = 0, HVEC1, HVEC2 }
 Format in rows of vectorial base functions. More...
 
enum  HVecDiffFormatting {
  HVEC0_0 = 0, HVEC1_0, HVEC2_0, HVEC0_1,
  HVEC1_1, HVEC2_1, HVEC0_2, HVEC1_2,
  HVEC2_2
}
 Format in rows of vectorial base gradients of base functions. More...
 
enum  VERBOSITY_LEVELS {
  DEFAULT_VERBOSITY = -1, QUIET = 0, VERBOSE, VERY_VERBOSE,
  NOISY, VERY_NOISY
}
 Verbosity levels. More...
 

Functions

DEPRECATED void macro_is_deprecated_using_deprecated_function ()
 Is used to indicate that macro is deprecated Do nothing just triggers error at the compilation. More...
 

Variables

const static char *const MoFEMErrorCodesNames []
 
const static char *const ApproximationBaseNames []
 
const static char *const FieldSpaceNames []
 
const static char *const FieldContinuityNames []
 
const static char *const CoordinateTypesNames []
 Coordinate system names. More...
 
const static char *const CubitBCNames []
 Names of types of sets and boundary conditions. More...
 

Detailed Description

useful compiler directives and definitions

Definition in file definitions.h.

Macro Definition Documentation

◆ BARRIER_MOFEM_RANK_END

#define BARRIER_MOFEM_RANK_END (   MOFEM)
Value:
{ \
for (int i = (MOFEM)->get_comm_rank(); i < (MOFEM)->get_comm_size(); i++) \
MPI_Barrier((MOFEM)->get_comm()); \
};

set barrier start Run code in sequence, starting from process 0, and ends on last process.

It can be only used for testing. Do not use that function as a part of these code.

Definition at line 323 of file definitions.h.

◆ BARRIER_MOFEM_RANK_START

#define BARRIER_MOFEM_RANK_START (   MOFEM)
Value:
{ \
for (int i = 0; i < (MOFEM)->get_comm_rank(); i++) \
MPI_Barrier((MOFEM)->get_comm()); \
};

set barrier start Run code in sequence, starting from process 0, and ends on last process.

It can be only used for testing. Do not use that function as a part of these code.

Definition at line 310 of file definitions.h.

◆ BARRIER_PCOMM_RANK_END

#define BARRIER_PCOMM_RANK_END (   PCMB)
Value:
{ \
for (unsigned int i = PCMB->proc_config().proc_rank(); \
i < PCMB->proc_config().proc_size(); i++) \
MPI_Barrier(PCMB->proc_config().proc_comm()); \
};

set barrier start Run code in sequence, starting from process 0, and ends on last process.

It can be only used for testing. Do not use that function as a part of these code.

Definition at line 286 of file definitions.h.

◆ BARRIER_PCOMM_RANK_START

#define BARRIER_PCOMM_RANK_START (   PCMB)
Value:
{ \
for (unsigned int i = 0; i < PCMB->proc_config().proc_rank(); i++) \
MPI_Barrier(PCMB->proc_config().proc_comm()); \
};

set barrier start Run code in sequence, starting from process 0, and ends on last process.

It can be only used for testing. Do not use that function as a part of these code.

Definition at line 264 of file definitions.h.

◆ BARRIER_RANK_END

#define BARRIER_RANK_END (   PCMB)
Value:
{ \
macro_is_deprecated_using_deprecated_function(); \
for (unsigned int i = PCMB->proc_config().proc_rank(); \
i < PCMB->proc_config().proc_size(); i++) \
MPI_Barrier(PCMB->proc_config().proc_comm()); \
};
Deprecated:
Do use this macro, instead use BARRIER_PCOMM_RANK_START

Definition at line 295 of file definitions.h.

◆ BARRIER_RANK_START

#define BARRIER_RANK_START (   PCMB)
Value:
{ \
macro_is_deprecated_using_deprecated_function(); \
for (unsigned int i = 0; i < PCMB->proc_config().proc_rank(); i++) \
MPI_Barrier(PCMB->proc_config().proc_comm()); \
};
Deprecated:
Do use this macro, instead use BARRIER_PCOMM_RANK_START

Definition at line 272 of file definitions.h.

◆ BITFEID_SIZE

#define BITFEID_SIZE   32

max number of finite elements

Definition at line 234 of file definitions.h.

◆ BITFIELDID_SIZE

#define BITFIELDID_SIZE   32

max number of fields

Definition at line 233 of file definitions.h.

◆ BITINTERFACEUID_SIZE

#define BITINTERFACEUID_SIZE   32

Definition at line 236 of file definitions.h.

◆ BITPROBLEMID_SIZE

#define BITPROBLEMID_SIZE   32

max number of problems

Definition at line 235 of file definitions.h.

◆ BITREFEDGES_SIZE

#define BITREFEDGES_SIZE   32

number refined edges

Definition at line 231 of file definitions.h.

◆ BITREFLEVEL_SIZE

#define BITREFLEVEL_SIZE   64

max number of refinements

Examples
free_surface.cpp, and hanging_node_approx.cpp.

Definition at line 232 of file definitions.h.

◆ CATCH_ERRORS

#define CATCH_ERRORS
Value:
catch (MoFEMExceptionInitial const &ex) { \
return PetscError(PETSC_COMM_SELF, ex.lINE, PETSC_FUNCTION_NAME, __FILE__, \
ex.errorCode, PETSC_ERROR_INITIAL, ex.what()); \
} \
catch (MoFEMExceptionRepeat const &ex) { \
return PetscError(PETSC_COMM_SELF, ex.lINE, PETSC_FUNCTION_NAME, __FILE__, \
ex.errorCode, PETSC_ERROR_REPEAT, " "); \
} \
catch (MoFEMException const &ex) { \
SETERRQ(PETSC_COMM_SELF, ex.errorCode, ex.errorMessage); \
} \
catch (boost::bad_weak_ptr & ex) { \
std::string message("Boost bad weak ptr: " + std::string(ex.what()) + \
" at " + boost::lexical_cast<std::string>(__LINE__) + \
" : " + std::string(__FILE__) + " in " + \
std::string(PETSC_FUNCTION_NAME)); \
SETERRQ(PETSC_COMM_SELF, MOFEM_STD_EXCEPTION_THROW, message.c_str()); \
} \
catch (std::out_of_range & ex) { \
std::string message("Std out of range error: " + std::string(ex.what()) + \
" at " + boost::lexical_cast<std::string>(__LINE__) + \
" : " + std::string(__FILE__) + " in " + \
std::string(PETSC_FUNCTION_NAME)); \
SETERRQ(PETSC_COMM_SELF, MOFEM_STD_EXCEPTION_THROW, message.c_str()); \
} \
catch (std::exception const &ex) { \
std::string message("Std error: " + std::string(ex.what()) + " at " + \
boost::lexical_cast<std::string>(__LINE__) + " : " + \
std::string(__FILE__) + " in " + \
std::string(PETSC_FUNCTION_NAME)); \
SETERRQ(PETSC_COMM_SELF, MOFEM_STD_EXCEPTION_THROW, message.c_str()); \
}

Catch errors.

Usage in main functions

int main(int argc, char *argv[]) {
MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
try {
// More code here
}
}
Examples
add_blockset.cpp, add_cubit_meshsets.cpp, adolc_plasticity.cpp, 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, cell_forces.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, delete_ho_nodes.cpp, dg_projection.cpp, dm_build_partitioned_mesh.cpp, dm_create_subdm.cpp, dm_partitioned_no_field.cpp, dynamic_first_order_con_law.cpp, edge_and_bubble_shape_functions_on_quad.cpp, eigen_elastic.cpp, elasticity.cpp, elasticity_mixed_formulation.cpp, ep.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, gauss_points_on_outer_product.cpp, hanging_node_approx.cpp, hcurl_check_approx_in_2d.cpp, hcurl_curl_operator.cpp, hcurl_divergence_operator_2d.cpp, hdiv_divergence_operator.cpp, heat_equation.cpp, heat_method.cpp, hello_world.cpp, helmholtz.cpp, hertz_surface.cpp, higher_derivatives.cpp, level_set.cpp, log.cpp, loop_entities.cpp, lorentz_force.cpp, magnetostatic.cpp, matrix_function.cpp, 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, navier_stokes.cpp, node_merge.cpp, nonlinear_dynamics.cpp, nonlinear_elastic.cpp, operators_tests.cpp, partition_mesh.cpp, petsc_smart_ptr_objects.cpp, phase.cpp, photon_diffusion.cpp, plastic.cpp, plate.cpp, 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, unsaturated_transport.cpp, wave_equation.cpp, and wavy_surface.cpp.

Definition at line 385 of file definitions.h.

◆ CHK_MOAB_THROW

#define CHK_MOAB_THROW (   err,
  msg 
)
Value:
{ \
if (PetscUnlikely(static_cast<int>(MB_SUCCESS) != (err))) \
{ \
std::string str; \
throw MoFEMException( \
("MOAB error (" + boost::lexical_cast<std::string>((err)) + ") " + \
boost::lexical_cast<std::string>((msg)) + " at line " + \
boost::lexical_cast<std::string>(__LINE__) + " : " + \
std::string(__FILE__)) \
.c_str()); \
} \
}

Check error code of MoAB function and throw MoFEM exception.

Parameters
errMoABErrorCode
msgerror message
Examples
ContactOps.hpp, EshelbianPlasticity.cpp, free_surface.cpp, level_set.cpp, plastic.cpp, and schur_test_diag_mat.cpp.

Definition at line 589 of file definitions.h.

◆ CHK_THROW_MESSAGE

#define CHK_THROW_MESSAGE (   err,
  msg 
)
Value:
{ \
if (PetscUnlikely((err) != MOFEM_SUCCESS)) \
THROW_MESSAGE(msg); \
}

Check and throw MoFEM exception.

Parameters
errerror code
msgmessage
Examples
adolc_plasticity.cpp, ADOLCPlasticity.hpp, ADOLCPlasticityMaterialModels.hpp, ContactOps.hpp, free_surface.cpp, level_set.cpp, nonlinear_dynamics.cpp, plastic.cpp, PlasticOps.hpp, seepage.cpp, test_broken_space.cpp, and thermo_elastic.cpp.

Definition at line 609 of file definitions.h.

◆ CHKERR

#define CHKERR   MoFEM::ErrorChecker<__LINE__>() <<

Inline error check.

// Call other functions
CHKERR fun_moab();
CHKERR fun_petsc();
CHKERR fun_mofem();
// Throw error
SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY, "Some error message");
}
int main(int argc, char *argv[]) {
// Initailise MoFEM and Petsc
MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
try {
moab::Core mb_instance; // MoAB database
moab::Interface &moab = mb_instance;
MoFEM::Core core(moab); // MOFEM database
MoFEM::CoreInterface &m_field = core;
CHKERR foo(); // Call function
}
}
Examples
add_blockset.cpp, add_cubit_meshsets.cpp, adolc_plasticity.cpp, ADOLCPlasticity.hpp, ADOLCPlasticityMaterialModels.hpp, analytical_nonlinear_poisson.cpp, analytical_poisson.cpp, analytical_poisson_field_split.cpp, approx_sphere.cpp, AuxPoissonFunctions.hpp, bernstein_bezier_generate_base.cpp, bone_adaptation.cpp, boundary_marker.cpp, build_large_problem.cpp, build_problems.cpp, cell_forces.cpp, child_and_parent.cpp, ContactOps.hpp, 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, delete_ho_nodes.cpp, dg_projection.cpp, dm_build_partitioned_mesh.cpp, dm_create_subdm.cpp, dm_partitioned_no_field.cpp, dynamic_first_order_con_law.cpp, edge_and_bubble_shape_functions_on_quad.cpp, eigen_elastic.cpp, elasticity.cpp, elasticity_mixed_formulation.cpp, ElasticityMixedFormulation.hpp, ep.cpp, EshelbianOperators.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, gauss_points_on_outer_product.cpp, hanging_node_approx.cpp, hcurl_check_approx_in_2d.cpp, hcurl_curl_operator.cpp, hcurl_divergence_operator_2d.cpp, hdiv_divergence_operator.cpp, heat_equation.cpp, heat_method.cpp, hello_world.cpp, helmholtz.cpp, hertz_surface.cpp, higher_derivatives.cpp, HookeElement.cpp, HookeElement.hpp, HookeInternalStressElement.hpp, level_set.cpp, log.cpp, loop_entities.cpp, lorentz_force.cpp, MagneticElement.hpp, magnetostatic.cpp, matrix_function.cpp, 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, navier_stokes.cpp, NavierStokesElement.cpp, NavierStokesElement.hpp, NeoHookean.hpp, node_merge.cpp, nonlinear_dynamics.cpp, nonlinear_elastic.cpp, NonlinearElasticElementInterface.hpp, operators_tests.cpp, partition_mesh.cpp, petsc_smart_ptr_objects.cpp, phase.cpp, photon_diffusion.cpp, plastic.cpp, PlasticOps.hpp, PlasticOpsMonitor.hpp, plate.cpp, plot_base.cpp, poisson_2d_dis_galerkin.cpp, poisson_2d_homogeneous.cpp, PoissonDiscontinousGalerkin.hpp, PoissonOperators.hpp, prism_elements_from_surface.cpp, prism_polynomial_approximation.cpp, quad_polynomial_approximation.cpp, reaction_diffusion.cpp, Remodeling.cpp, Remodeling.hpp, 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, wave_equation.cpp, and wavy_surface.cpp.

Definition at line 548 of file definitions.h.

◆ CHKERRG

#define CHKERRG (   n)
Value:
if ((boost::is_same<BOOST_TYPEOF((n)), \
MoFEMErrorCodeGeneric<PetscErrorCode>>::value)) { \
CHKERRQ_PETSC((n)); \
} else if (boost::is_same<BOOST_TYPEOF((n)), \
MoFEMErrorCodeGeneric<moab::ErrorCode>>::value) { \
CHKERRQ_MOAB((n)); \
}

Check error code of MoFEM/MOAB/PETSc function.

Parameters
aMoFEMErrorCode
rval = fun_moab(); CHKERRG(rval);
ierr = fun_petsc(); CHKERRG(ierr);
merr = fun_mofem(); CHKERRG(merr);
Note
Function detect type of errocode using specialized template function getErrorType, i.e. condition is evaluated at compilation time.
Examples
analytical_nonlinear_poisson.cpp, analytical_poisson.cpp, analytical_poisson_field_split.cpp, delete_ho_nodes.cpp, elasticity.cpp, EshelbianPlasticity.cpp, field_to_vertices.cpp, hertz_surface.cpp, lorentz_force.cpp, magnetostatic.cpp, mesh_cut.cpp, mesh_smoothing.cpp, minimal_surface_area.cpp, simple_elasticity.cpp, split_sideset.cpp, unsaturated_transport.cpp, UnsaturatedFlow.hpp, and wavy_surface.cpp.

Definition at line 496 of file definitions.h.

◆ CHKERRQ_MOAB

#define CHKERRQ_MOAB (   a)
Value:
if (PetscUnlikely(MB_SUCCESS != (a))) { \
std::string error_str = (unsigned)(a) <= (unsigned)MB_FAILURE \
? moab::ErrorCodeStr[a] \
: "INVALID ERROR CODE"; \
std::string str("MOAB error (" + boost::lexical_cast<std::string>((a)) + \
") " + error_str + " at line " + \
boost::lexical_cast<std::string>(__LINE__) + " : " + \
std::string(__FILE__)); \
SETERRQ(PETSC_COMM_SELF, MOFEM_MOAB_ERROR, str.c_str()); \
}

check error code of MoAB function

Parameters
aMoABErrorCode

Definition at line 467 of file definitions.h.

◆ CHKERRQ_PETSC

#define CHKERRQ_PETSC (   n)    CHKERRQ(n)

Definition at line 462 of file definitions.h.

◆ DEPRECATED

#define DEPRECATED
Examples
PlasticOpsGeneric.hpp, and ThermoElasticOps.hpp.

Definition at line 17 of file definitions.h.

◆ DOF_UID_MASK

#define DOF_UID_MASK   (MAX_DOFS_ON_ENTITY - 1)

Mask for DOF number on entity form UId.

Definition at line 251 of file definitions.h.

◆ ENTITY_UID_MASK

#define ENTITY_UID_MASK   (~DOF_UID_MASK)

Definition at line 253 of file definitions.h.

◆ MAX_CORE_TMP

#define MAX_CORE_TMP   1

maximal number of cores

Definition at line 230 of file definitions.h.

◆ MAX_DOFS_ON_ENTITY

#define MAX_DOFS_ON_ENTITY   512

Maximal number of DOFs on entity.

Examples
EshelbianPlasticity.cpp, level_set.cpp, remove_entities_from_problem.cpp, and tensor_divergence_operator.cpp.

Definition at line 249 of file definitions.h.

◆ MAX_PROCESSORS_NUMBER

#define MAX_PROCESSORS_NUMBER   1024

Maximal number of processors.

Definition at line 250 of file definitions.h.

◆ MB_END_ID

#define MB_END_ID   ((EntityID)MB_ID_MASK)

Last id is the complement of the MASK.

Definition at line 245 of file definitions.h.

◆ MB_ID_MASK

#define MB_ID_MASK   (~MB_TYPE_MASK)

Definition at line 247 of file definitions.h.

◆ MB_ID_WIDTH

#define MB_ID_WIDTH   (8 * sizeof(EntityHandle) - MB_TYPE_WIDTH)

Definition at line 240 of file definitions.h.

◆ MB_START_ID

#define MB_START_ID   ((EntityID)1)

All entity id's currently start at 1.

Definition at line 244 of file definitions.h.

◆ MB_TYPE_MASK

#define MB_TYPE_MASK   ((EntityHandle)0xF << MB_ID_WIDTH)

Definition at line 241 of file definitions.h.

◆ MB_TYPE_WIDTH

#define MB_TYPE_WIDTH   4

Definition at line 239 of file definitions.h.

◆ MOAB_THROW

#define MOAB_THROW (   err)
Value:
{ \
if (PetscUnlikely(MB_SUCCESS != (err))) { \
std::string error_str = (unsigned)(err) <= (unsigned)MB_FAILURE \
? moab::ErrorCodeStr[err] \
: "INVALID ERROR CODE"; \
throw MoFEMException(MOFEM_MOAB_ERROR, \
("MOAB error (" + \
boost::lexical_cast<std::string>((err)) + ") " + \
error_str + " at line " + \
boost::lexical_cast<std::string>(__LINE__) + \
" : " + std::string(__FILE__)) \
.c_str()); \
} \
}

Check error code of MoAB function and throw MoFEM exception.

Parameters
errMoABErrorCode
Examples
nonlinear_dynamics.cpp, phase.cpp, and plate.cpp.

Definition at line 554 of file definitions.h.

◆ MoFEMFunctionBegin

#define MoFEMFunctionBegin
Value:
PetscFunctionBegin; \
try {

First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions should be MoFEMFunctionReturn(0);.

\node Not collective

Example

PetscErrorCode fun() {
int something;
}
Examples
add_blockset.cpp, adolc_plasticity.cpp, ADOLCPlasticity.hpp, ADOLCPlasticityMaterialModels.hpp, analytical_poisson_field_split.cpp, approx_sphere.cpp, AuxPoissonFunctions.hpp, bernstein_bezier_generate_base.cpp, boundary_marker.cpp, build_large_problem.cpp, child_and_parent.cpp, ContactOps.hpp, 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, dg_projection.cpp, dm_build_partitioned_mesh.cpp, dm_partitioned_no_field.cpp, dynamic_first_order_con_law.cpp, eigen_elastic.cpp, elasticity.cpp, ElasticityMixedFormulation.hpp, EshelbianOperators.cpp, EshelbianPlasticity.cpp, field_blas.cpp, field_evaluator.cpp, forces_and_sources_testing_edge_element.cpp, forces_and_sources_testing_users_base.cpp, free_surface.cpp, gauss_points_on_outer_product.cpp, hanging_node_approx.cpp, hcurl_check_approx_in_2d.cpp, hcurl_curl_operator.cpp, hcurl_divergence_operator_2d.cpp, hdiv_divergence_operator.cpp, heat_equation.cpp, heat_method.cpp, helmholtz.cpp, higher_derivatives.cpp, HookeElement.cpp, HookeElement.hpp, HookeInternalStressElement.hpp, level_set.cpp, log.cpp, loop_entities.cpp, lorentz_force.cpp, MagneticElement.hpp, mesh_cut.cpp, mesh_smoothing.cpp, mixed_poisson.cpp, mortar_contact_thermal.cpp, navier_stokes.cpp, NavierStokesElement.cpp, NavierStokesElement.hpp, NeoHookean.hpp, nonlinear_dynamics.cpp, nonlinear_elastic.cpp, NonlinearElasticElementInterface.hpp, operators_tests.cpp, petsc_smart_ptr_objects.cpp, phase.cpp, photon_diffusion.cpp, plastic.cpp, PlasticOps.hpp, PlasticOpsGeneric.hpp, PlasticOpsLargeStrains.hpp, PlasticOpsMonitor.hpp, PlasticOpsSmallStrains.hpp, plate.cpp, plot_base.cpp, poisson_2d_dis_galerkin.cpp, poisson_2d_homogeneous.cpp, poisson_2d_homogeneous.hpp, PoissonDiscontinousGalerkin.hpp, PoissonOperators.hpp, prism_elements_from_surface.cpp, prism_polynomial_approximation.cpp, quad_polynomial_approximation.cpp, reaction_diffusion.cpp, Remodeling.cpp, Remodeling.hpp, remove_entities_from_problem.cpp, remove_entities_from_problem_not_partitioned.cpp, scalar_check_approximation.cpp, schur_test_diag_mat.cpp, seepage.cpp, SeepageOps.hpp, 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_scaled_with_density_element.cpp, thermo_elastic.cpp, ThermoElasticOps.hpp, UnsaturatedFlow.hpp, and wave_equation.cpp.

Definition at line 359 of file definitions.h.

◆ MoFEMFunctionBeginHot

#define MoFEMFunctionBeginHot   PetscFunctionBeginHot

◆ MoFEMFunctionReturn

#define MoFEMFunctionReturn (   a)
Value:
} \
CATCH_ERRORS \
PetscFunctionReturn(a)

Last executable line of each PETSc function used for error handling. Replaces return()

Parameters
aerror code
Note
MoFEMFunctionReturn has to be used with MoFEMFunctionBegin and can be used only at the end of the function. If is need to return function in earlier use MoFEMFunctionReturnHot
Examples
add_blockset.cpp, adolc_plasticity.cpp, ADOLCPlasticity.hpp, ADOLCPlasticityMaterialModels.hpp, analytical_poisson_field_split.cpp, approx_sphere.cpp, AuxPoissonFunctions.hpp, bernstein_bezier_generate_base.cpp, boundary_marker.cpp, build_large_problem.cpp, child_and_parent.cpp, ContactOps.hpp, 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, dg_projection.cpp, dm_build_partitioned_mesh.cpp, dm_partitioned_no_field.cpp, dynamic_first_order_con_law.cpp, eigen_elastic.cpp, elasticity.cpp, ElasticityMixedFormulation.hpp, EshelbianOperators.cpp, EshelbianPlasticity.cpp, field_blas.cpp, field_evaluator.cpp, forces_and_sources_testing_edge_element.cpp, forces_and_sources_testing_users_base.cpp, free_surface.cpp, gauss_points_on_outer_product.cpp, hanging_node_approx.cpp, hcurl_check_approx_in_2d.cpp, hcurl_curl_operator.cpp, hcurl_divergence_operator_2d.cpp, hdiv_divergence_operator.cpp, heat_equation.cpp, heat_method.cpp, hello_world.cpp, helmholtz.cpp, higher_derivatives.cpp, HookeElement.cpp, HookeElement.hpp, HookeInternalStressElement.hpp, level_set.cpp, log.cpp, loop_entities.cpp, lorentz_force.cpp, MagneticElement.hpp, mesh_cut.cpp, mesh_smoothing.cpp, mixed_poisson.cpp, mortar_contact_thermal.cpp, navier_stokes.cpp, NavierStokesElement.cpp, NavierStokesElement.hpp, NeoHookean.hpp, nonlinear_dynamics.cpp, nonlinear_elastic.cpp, NonlinearElasticElementInterface.hpp, operators_tests.cpp, petsc_smart_ptr_objects.cpp, phase.cpp, photon_diffusion.cpp, plastic.cpp, PlasticOps.hpp, PlasticOpsGeneric.hpp, PlasticOpsLargeStrains.hpp, PlasticOpsMonitor.hpp, PlasticOpsSmallStrains.hpp, plate.cpp, plot_base.cpp, poisson_2d_dis_galerkin.cpp, poisson_2d_homogeneous.cpp, poisson_2d_homogeneous.hpp, PoissonDiscontinousGalerkin.hpp, PoissonOperators.hpp, prism_elements_from_surface.cpp, prism_polynomial_approximation.cpp, quad_polynomial_approximation.cpp, reaction_diffusion.cpp, Remodeling.cpp, Remodeling.hpp, remove_entities_from_problem.cpp, remove_entities_from_problem_not_partitioned.cpp, scalar_check_approximation.cpp, schur_test_diag_mat.cpp, seepage.cpp, SeepageOps.hpp, 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_scaled_with_density_element.cpp, thermo_elastic.cpp, ThermoElasticOps.hpp, UnsaturatedFlow.hpp, and wave_equation.cpp.

Definition at line 429 of file definitions.h.

◆ MoFEMFunctionReturnHot

#define MoFEMFunctionReturnHot (   a)    PetscFunctionReturn(a)

Last executable line of each PETSc function used for error handling. Replaces return()

Parameters
aerror code
Examples
adolc_plasticity.cpp, ADOLCPlasticityMaterialModels.hpp, analytical_poisson_field_split.cpp, approx_sphere.cpp, boundary_marker.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, ElasticityMixedFormulation.hpp, EshelbianOperators.cpp, EshelbianPlasticity.cpp, field_blas.cpp, forces_and_sources_testing_flat_prism_element.cpp, forces_and_sources_testing_users_base.cpp, hanging_node_approx.cpp, hcurl_curl_operator.cpp, hcurl_divergence_operator_2d.cpp, hdiv_divergence_operator.cpp, heat_equation.cpp, HookeElement.cpp, HookeElement.hpp, HookeInternalStressElement.hpp, loop_entities.cpp, MagneticElement.hpp, mesh_smoothing.cpp, mortar_contact_thermal.cpp, NavierStokesElement.cpp, NeoHookean.hpp, nonlinear_dynamics.cpp, nonlinear_elastic.cpp, NonlinearElasticElementInterface.hpp, photon_diffusion.cpp, plastic.cpp, PlasticOps.hpp, plate.cpp, plot_base.cpp, PoissonDiscontinousGalerkin.hpp, PoissonOperators.hpp, prism_elements_from_surface.cpp, Remodeling.cpp, Remodeling.hpp, seepage.cpp, shallow_wave.cpp, simple_contact_thermal.cpp, simple_elasticity.cpp, simple_interface.cpp, testing_jacobian_of_hook_scaled_with_density_element.cpp, thermo_elastic.cpp, UnsaturatedFlow.hpp, and wave_equation.cpp.

Definition at line 460 of file definitions.h.

◆ MYPCOMM_INDEX

#define MYPCOMM_INDEX   0

◆ NOT_USED

#define NOT_USED (   x)    ((void)(x))

Definition at line 255 of file definitions.h.

◆ SSTR

#define SSTR (   x)    toString(x)

Convert number to string.

Parameters
xnumber

Definition at line 619 of file definitions.h.

◆ SYMMETRIC_TENSOR2_MAT_PTR

#define SYMMETRIC_TENSOR2_MAT_PTR (   MAT)    &MAT(0, 0), &MAT(0, 1), &MAT(0, 2), &MAT(0, 3), &MAT(0, 4), &MAT(0, 5)

Definition at line 637 of file definitions.h.

◆ SYMMETRIC_TENSOR2_VEC_PTR

#define SYMMETRIC_TENSOR2_VEC_PTR (   VEC)    &VEC[0], &VEC[1], &VEC[2], &VEC[3], &VEC[4], &VEC[5]

Definition at line 640 of file definitions.h.

◆ SYMMETRIC_TENSOR4_MAT_PTR

#define SYMMETRIC_TENSOR4_MAT_PTR (   MAT)
Value:
&MAT(0, 0), &MAT(0, 1), &MAT(0, 2), &MAT(0, 3), &MAT(0, 4), &MAT(0, 5), \
&MAT(1, 0), &MAT(1, 1), &MAT(1, 2), &MAT(1, 3), &MAT(1, 4), &MAT(1, 5), \
&MAT(2, 0), &MAT(2, 1), &MAT(2, 2), &MAT(2, 3), &MAT(2, 4), &MAT(2, 5), \
&MAT(3, 0), &MAT(3, 1), &MAT(3, 2), &MAT(3, 3), &MAT(3, 4), &MAT(3, 5), \
&MAT(4, 0), &MAT(4, 1), &MAT(4, 2), &MAT(4, 3), &MAT(4, 4), &MAT(4, 5), \
&MAT(5, 0), &MAT(5, 1), &MAT(5, 2), &MAT(5, 3), &MAT(5, 4), &MAT(5, 5)

Definition at line 623 of file definitions.h.

◆ TENSOR1_VEC_PTR

#define TENSOR1_VEC_PTR (   VEC)    &VEC[0], &VEC[1], &VEC[2]

Definition at line 621 of file definitions.h.

◆ TENSOR2_MAT_PTR

#define TENSOR2_MAT_PTR (   MAT)
Value:
&MAT(0, 0), &MAT(1, 0), &MAT(2, 0), &MAT(3, 0), &MAT(4, 0), &MAT(5, 0), \
&MAT(6, 0), &MAT(7, 0), &MAT(8, 0)

Definition at line 633 of file definitions.h.

◆ TENSOR4_MAT_PTR

#define TENSOR4_MAT_PTR (   MAT)    &MAT(0, 0), MAT.size2()

Definition at line 631 of file definitions.h.

◆ THROW_MESSAGE

#define THROW_MESSAGE (   msg)
Value:
{ \
("MoFEM error " + boost::lexical_cast<std::string>((msg)) + \
" at line " + boost::lexical_cast<std::string>(__LINE__) + " : " + \
std::string(__FILE__)) \
.c_str()); \
}

Throw MoFEM exception.

Parameters
msgmessage
Examples
matrix_function.cpp, and phase.cpp.

Definition at line 574 of file definitions.h.

Enumeration Type Documentation

◆ ByWhat

enum ByWhat

Controls adjency multi_index container (e.g. BYROW is adjacenciecy by field on on rows), see MoFEM::FieldEntityEntFiniteElementAdjacencyMap

Enumerator
BYROW 
BYCOL 
BYDATA 
BYROWDATA 
BYCOLDATA 
BYROWCOL 
BYALL 

Definition at line 143 of file definitions.h.

143  {
144  BYROW = 1 << 0,
145  BYCOL = 1 << 1,
146  BYDATA = 1 << 2,
147  BYROWDATA = 1 << 0 | 1 << 2,
148  BYCOLDATA = 1 << 1 | 1 << 2,
149  BYROWCOL = 1 << 0 | 1 << 1,
150  BYALL = 1 << 0 | 1 << 1 | 1 << 2
151 };

◆ CoordinateTypes

Coodinate system.

Enumerator
CARTESIAN 
POLAR 
CYLINDRICAL 
SPHERICAL 
LAST_COORDINATE_SYSTEM 

Definition at line 127 of file definitions.h.

127  {
128  CARTESIAN,
129  POLAR,
130  CYLINDRICAL,
131  SPHERICAL,
133 };

◆ CubitBC

enum CubitBC

Types of sets and boundary conditions.

Enumerator
UNKNOWNSET 
NODESET 
SIDESET 
BLOCKSET 
MATERIALSET 
DISPLACEMENTSET 
FORCESET 
PRESSURESET 
VELOCITYSET 
ACCELERATIONSET 
TEMPERATURESET 
HEATFLUXSET 
INTERFACESET 
UNKNOWNNAME 
MAT_ELASTICSET 

block name is "MAT_ELASTIC"

MAT_INTERFSET 
MAT_THERMALSET 

block name is "MAT_THERMAL"

BODYFORCESSET 

block name is "BODY_FORCES"

MAT_MOISTURESET 

block name is "MAT_MOISTURE"

DIRICHLET_BC 
NEUMANN_BC 
LASTSET_BC 

Definition at line 157 of file definitions.h.

157  {
158  UNKNOWNSET = 0,
159  NODESET = 1 << 0,
160  SIDESET = 1 << 1,
161  BLOCKSET = 1 << 2,
162  MATERIALSET = 1 << 3,
163  DISPLACEMENTSET = 1 << 4,
164  FORCESET = 1 << 5,
165  PRESSURESET = 1 << 6,
166  VELOCITYSET = 1 << 7,
167  ACCELERATIONSET = 1 << 8,
168  TEMPERATURESET = 1 << 9,
169  HEATFLUXSET = 1 << 10,
170  INTERFACESET = 1 << 11,
171  UNKNOWNNAME = 1 << 12,
172  MAT_ELASTICSET = 1 << 13, ///< block name is "MAT_ELASTIC"
173  MAT_INTERFSET = 1 << 14,
174  MAT_THERMALSET = 1 << 15, ///< block name is "MAT_THERMAL"
175  BODYFORCESSET = 1 << 16, ///< block name is "BODY_FORCES"
176  MAT_MOISTURESET = 1 << 17, ///< block name is "MAT_MOISTURE"
177  DIRICHLET_BC = 1 << 18,
178  NEUMANN_BC = 1 << 19,
179  LASTSET_BC = 1 << 20
180 };

◆ FieldApproximationBase

approximation base

Enumerator
NOBASE 
AINSWORTH_LEGENDRE_BASE 

Ainsworth Cole (Legendre) approx. base [1].

AINSWORTH_LOBATTO_BASE 

Like AINSWORTH_LEGENDRE_BASE but with Lobatto base instead Legendre [13]

AINSWORTH_BERNSTEIN_BEZIER_BASE 

See [3] and [2]

DEMKOWICZ_JACOBI_BASE 

Construction of base is by Demkowicz [29]

USER_BASE 

user implemented approximation base

LASTBASE 

Definition at line 58 of file definitions.h.

58  {
59  NOBASE = 0,
61  1, ///< Ainsworth Cole (Legendre) approx. base \cite NME:NME847
62  AINSWORTH_LOBATTO_BASE, ///< Like AINSWORTH_LEGENDRE_BASE but with Lobatto
63  ///< base instead Legendre \cite beriot2015efficient
64  AINSWORTH_BERNSTEIN_BEZIER_BASE, ///< See \cite ainsworth2011bernstein and
65  ///< \cite ainsworth2018bernstein
66  DEMKOWICZ_JACOBI_BASE, ///< Construction of base is by Demkowicz \cite
67  ///< fuentes2015orientation
68  USER_BASE, ///< user implemented approximation base
69  LASTBASE
70 };

◆ FieldContinuity

Field continuity.

Enumerator
CONTINUOUS 

Regular field.

DISCONTINUOUS 

Broken continuity (No effect on L2 space)

LASTCONTINUITY 

Definition at line 99 of file definitions.h.

99  {
100  CONTINUOUS = 0, ///< Regular field
101  DISCONTINUOUS = 1, ///< Broken continuity (No effect on L2 space)
103 };

◆ FieldSpace

enum FieldSpace

approximation spaces

Enumerator
NOSPACE 
NOFIELD 

scalar or vector of scalars describe (no true field)

H1 

continuous field

HCURL 

field with continuous tangents

HDIV 

field with continuous normal traction

L2 

field with C-1 continuity

LASTSPACE 

FieldSpace in [ 0, LASTSPACE )

Definition at line 82 of file definitions.h.

82  {
83  NOSPACE = 0,
84  NOFIELD = 1, ///< scalar or vector of scalars describe (no true field)
85  H1, ///< continuous field
86  HCURL, ///< field with continuous tangents
87  HDIV, ///< field with continuous normal traction
88  L2, ///< field with C-1 continuity
89  LASTSPACE ///< FieldSpace in [ 0, LASTSPACE )
90 };

◆ HVecDiffFormatting

Format in rows of vectorial base gradients of base functions.

Enumerator
HVEC0_0 
HVEC1_0 
HVEC2_0 
HVEC0_1 
HVEC1_1 
HVEC2_1 
HVEC0_2 
HVEC1_2 
HVEC2_2 

Definition at line 204 of file definitions.h.

204  {
205  HVEC0_0 = 0,
206  HVEC1_0,
207  HVEC2_0,
208  HVEC0_1,
209  HVEC1_1,
210  HVEC2_1,
211  HVEC0_2,
212  HVEC1_2,
213  HVEC2_2
214 };

◆ HVecFormatting

Format in rows of vectorial base functions.

Enumerator
HVEC0 
HVEC1 
HVEC2 

Definition at line 199 of file definitions.h.

199 { HVEC0 = 0, HVEC1, HVEC2 };

◆ MoFEMErrorCodes

Error handling.

This is complementary to PETSC error codes. The numerical values for these are defined in include/petscerror.h. The names are defined in err.c

MoAB error messages are defined in moab/Types.hpp

Enumerator
MOFEM_SUCCESS 
MOFEM_DATA_INCONSISTENCY 
MOFEM_NOT_IMPLEMENTED 
MOFEM_NOT_FOUND 
MOFEM_OPERATION_UNSUCCESSFUL 
MOFEM_IMPOSSIBLE_CASE 
MOFEM_INVALID_DATA 
MOFEM_NOT_INSTALLED 
MOFEM_MOFEMEXCEPTION_THROW 
MOFEM_STD_EXCEPTION_THROW 
MOFEM_ATOM_TEST_INVALID 
MOFEM_MOAB_ERROR 

Definition at line 29 of file definitions.h.

29  {
30  MOFEM_SUCCESS = 0,
33  MOFEM_NOT_FOUND = 102,
36  MOFEM_INVALID_DATA = 105,
37  MOFEM_NOT_INSTALLED = 106,
41  MOFEM_MOAB_ERROR = 110
42 };

◆ MoFEMTypes

enum MoFEMTypes

Those types control how functions respond on arguments, f.e. error handling.

Enumerator
MF_ZERO 
MF_EXCL 
MF_EXIST 
MF_NOT_THROW 

Definition at line 110 of file definitions.h.

110  {
111  MF_ZERO = 0,
112  MF_EXCL = 1 << 0,
113  MF_EXIST = 1 << 1,
114  MF_NOT_THROW = 1 << 2
115 };

◆ RowColData

enum RowColData

RowColData.

Enumerator
ROW 
COL 
DATA 
LASTROWCOLDATA 

Definition at line 136 of file definitions.h.

136 { ROW = 0, COL, DATA, LASTROWCOLDATA };

◆ VERBOSITY_LEVELS

Verbosity levels.

Enumerator
DEFAULT_VERBOSITY 
QUIET 
VERBOSE 
VERY_VERBOSE 
NOISY 
VERY_NOISY 

Definition at line 219 of file definitions.h.

219  {
220  DEFAULT_VERBOSITY = -1,
221  QUIET = 0,
222  VERBOSE,
223  VERY_VERBOSE,
224  NOISY,
225  VERY_NOISY
226 };

Function Documentation

◆ macro_is_deprecated_using_deprecated_function()

DEPRECATED void macro_is_deprecated_using_deprecated_function ( )

Is used to indicate that macro is deprecated Do nothing just triggers error at the compilation.

Definition at line 11 of file Core.cpp.

11 {}

Variable Documentation

◆ ApproximationBaseNames

const static char* const ApproximationBaseNames[]
static
Initial value:
= {
"NOBASE",
"AINSWORTH_LEGENDRE_BASE",
"AINSWORTH_LOBATTO_BASE",
"AINSWORTH_BERNSTEIN_BEZIER_BASE",
"DEMKOWICZ_JACOBI_BASE",
"USER_BASE",
"LASTBASE"}
Examples
forces_and_sources_testing_edge_element.cpp, hanging_node_approx.cpp, and plastic.cpp.

Definition at line 72 of file definitions.h.

◆ CoordinateTypesNames

const static char* const CoordinateTypesNames[]
static
Initial value:
= {"Cartesian", "Polar",
"Cylindrical", "Spherical"}

Coordinate system names.

Definition at line 121 of file definitions.h.

◆ CubitBCNames

const static char* const CubitBCNames[]
static
Initial value:
= {
"UNKNOWNSET", "NODESET", "SIDESET", "BLOCKSET",
"MATERIALSET", "DISPLACEMENTSET", "FORCESET", "PRESSURESET",
"VELOCITYSET", "ACCELERATIONSET", "TEMPERATURESET", "HEATFLUXSET",
"INTERFACESET", "UNKNOWNNAME", "MAT_ELASTICSET", "MAT_INTERFSET",
"MAT_THERMALSET", "BODYFORCESSET", "MAT_MOISTURESET", "DIRICHLET_BC",
"NEUMANN_BC", "LASTSET_BC"}

Names of types of sets and boundary conditions.

Definition at line 188 of file definitions.h.

◆ FieldContinuityNames

const static char* const FieldContinuityNames[]
static
Initial value:
= {"CONTINUOUS",
"DISCONTINUOUS"}

Definition at line 105 of file definitions.h.

◆ FieldSpaceNames

const static char* const FieldSpaceNames[]
static
Initial value:
= {
"NOSPACE", "NOFIELD", "H1", "HCURL", "HDIV", "L2", "LASTSPACE"}
Examples
hanging_node_approx.cpp.

Definition at line 92 of file definitions.h.

◆ MoFEMErrorCodesNames

const static char* const MoFEMErrorCodesNames[]
static
Initial value:
= {
"MOFEM_SUCCESS",
"MOFEM_DATA_INCONSISTENCY",
"MOFEM_NOT_IMPLEMENTED",
"MOFEM_NOT_FOUND",
"MOFEM_OPERATION_UNSUCCESSFUL",
"MOFEM_IMPOSSIBLE_CASE",
"MOFEM_INVALID_DATA",
"MOFEM_MOFEMEXCEPTION_THROW",
"MOFEM_STD_EXCEPTION_THROW",
"MOFEM_ATOM_TEST_INVALID",
"MOFEM_MOAB_ERROR"}

Definition at line 44 of file definitions.h.

NOSPACE
@ NOSPACE
Definition: definitions.h:83
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
DEFAULT_VERBOSITY
@ DEFAULT_VERBOSITY
Definition: definitions.h:220
SIDESET
@ SIDESET
Definition: definitions.h:160
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
MOFEM_STD_EXCEPTION_THROW
@ MOFEM_STD_EXCEPTION_THROW
Definition: definitions.h:39
LASTBASE
@ LASTBASE
Definition: definitions.h:69
CARTESIAN
@ CARTESIAN
Definition: definitions.h:128
HVEC0_2
@ HVEC0_2
Definition: definitions.h:211
NOISY
@ NOISY
Definition: definitions.h:224
PRESSURESET
@ PRESSURESET
Definition: definitions.h:165
NOBASE
@ NOBASE
Definition: definitions.h:59
MoFEM::Exceptions::MoFEMException
Exception to catch.
Definition: Exceptions.hpp:20
L2
@ L2
field with C-1 continuity
Definition: definitions.h:88
help
static char help[]
Definition: activate_deactivate_dofs.cpp:13
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
LASTCONTINUITY
@ LASTCONTINUITY
Definition: definitions.h:102
HVEC0_1
@ HVEC0_1
Definition: definitions.h:208
VERY_NOISY
@ VERY_NOISY
Definition: definitions.h:225
MOFEM_SUCCESS
@ MOFEM_SUCCESS
Definition: definitions.h:30
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
DATA
@ DATA
Definition: definitions.h:136
HVEC0_0
@ HVEC0_0
Definition: definitions.h:205
HVEC1_1
@ HVEC1_1
Definition: definitions.h:209
MOFEM_IMPOSSIBLE_CASE
@ MOFEM_IMPOSSIBLE_CASE
Definition: definitions.h:35
HVEC1
@ HVEC1
Definition: definitions.h:199
USER_BASE
@ USER_BASE
user implemented approximation base
Definition: definitions.h:68
UNKNOWNNAME
@ UNKNOWNNAME
Definition: definitions.h:171
ROW
@ ROW
Definition: definitions.h:136
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2002
HVEC2_1
@ HVEC2_1
Definition: definitions.h:210
NODESET
@ NODESET
Definition: definitions.h:159
MATERIALSET
@ MATERIALSET
Definition: definitions.h:162
MoFEM::Exceptions::rval
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
Definition: Exceptions.hpp:74
VERBOSE
@ VERBOSE
Definition: definitions.h:222
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
BYALL
@ BYALL
Definition: definitions.h:150
HVEC2_2
@ HVEC2_2
Definition: definitions.h:213
HVEC1_2
@ HVEC1_2
Definition: definitions.h:212
a
constexpr double a
Definition: approx_sphere.cpp:30
POLAR
@ POLAR
Definition: definitions.h:129
FORCESET
@ FORCESET
Definition: definitions.h:164
MAT_THERMALSET
@ MAT_THERMALSET
block name is "MAT_THERMAL"
Definition: definitions.h:174
BODYFORCESSET
@ BODYFORCESSET
block name is "BODY_FORCES"
Definition: definitions.h:175
ACCELERATIONSET
@ ACCELERATIONSET
Definition: definitions.h:167
MAT_ELASTICSET
@ MAT_ELASTICSET
block name is "MAT_ELASTIC"
Definition: definitions.h:172
DISPLACEMENTSET
@ DISPLACEMENTSET
Definition: definitions.h:163
BYROWDATA
@ BYROWDATA
Definition: definitions.h:147
COL
@ COL
Definition: definitions.h:136
SPHERICAL
@ SPHERICAL
Definition: definitions.h:131
LASTSPACE
@ LASTSPACE
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:89
INTERFACESET
@ INTERFACESET
Definition: definitions.h:170
MOFEM_NOT_INSTALLED
@ MOFEM_NOT_INSTALLED
Definition: definitions.h:37
AINSWORTH_LOBATTO_BASE
@ AINSWORTH_LOBATTO_BASE
Definition: definitions.h:62
MOFEM_OPERATION_UNSUCCESSFUL
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition: definitions.h:34
BYDATA
@ BYDATA
Definition: definitions.h:146
UNKNOWNSET
@ UNKNOWNSET
Definition: definitions.h:158
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
MOFEM_MOFEMEXCEPTION_THROW
@ MOFEM_MOFEMEXCEPTION_THROW
Definition: definitions.h:38
MOFEM_NOT_FOUND
@ MOFEM_NOT_FOUND
Definition: definitions.h:33
main
int main(int argc, char *argv[])
Definition: activate_deactivate_dofs.cpp:15
VELOCITYSET
@ VELOCITYSET
Definition: definitions.h:166
AINSWORTH_BERNSTEIN_BEZIER_BASE
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
Definition: definitions.h:64
convert.n
n
Definition: convert.py:82
BYCOLDATA
@ BYCOLDATA
Definition: definitions.h:148
LASTSET_BC
@ LASTSET_BC
Definition: definitions.h:179
MoFEM::CoreTmp< 0 >::Initialize
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:72
MF_ZERO
@ MF_ZERO
Definition: definitions.h:111
BYCOL
@ BYCOL
Definition: definitions.h:145
CYLINDRICAL
@ CYLINDRICAL
Definition: definitions.h:130
DEMKOWICZ_JACOBI_BASE
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:66
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1140
BLOCKSET
@ BLOCKSET
Definition: definitions.h:161
DISCONTINUOUS
@ DISCONTINUOUS
Broken continuity (No effect on L2 space)
Definition: definitions.h:101
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
HCURL
@ HCURL
field with continuous tangents
Definition: definitions.h:86
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
LAST_COORDINATE_SYSTEM
@ LAST_COORDINATE_SYSTEM
Definition: definitions.h:132
TEMPERATURESET
@ TEMPERATURESET
Definition: definitions.h:168
HVEC2_0
@ HVEC2_0
Definition: definitions.h:207
DIRICHLET_BC
@ DIRICHLET_BC
Definition: definitions.h:177
MAT_INTERFSET
@ MAT_INTERFSET
Definition: definitions.h:173
BYROW
@ BYROW
Definition: definitions.h:144
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
NEUMANN_BC
@ NEUMANN_BC
Definition: definitions.h:178
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
HEATFLUXSET
@ HEATFLUXSET
Definition: definitions.h:169
QUIET
@ QUIET
Definition: definitions.h:221
LASTROWCOLDATA
@ LASTROWCOLDATA
Definition: definitions.h:136
MoFEM::CoreInterface
Interface.
Definition: Interface.hpp:27
HVEC2
@ HVEC2
Definition: definitions.h:199
MF_NOT_THROW
@ MF_NOT_THROW
Definition: definitions.h:114
fun
auto fun
Function to approximate.
Definition: dg_projection.cpp:36
MF_EXCL
@ MF_EXCL
Definition: definitions.h:112
CONTINUOUS
@ CONTINUOUS
Regular field.
Definition: definitions.h:100
MF_EXIST
@ MF_EXIST
Definition: definitions.h:113
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
HDIV
@ HDIV
field with continuous normal traction
Definition: definitions.h:87
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
MAT_MOISTURESET
@ MAT_MOISTURESET
block name is "MAT_MOISTURE"
Definition: definitions.h:176
CHKERRG
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:496
MOFEM_MOAB_ERROR
@ MOFEM_MOAB_ERROR
Definition: definitions.h:41
VERY_VERBOSE
@ VERY_VERBOSE
Definition: definitions.h:223
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MOFEM_INVALID_DATA
@ MOFEM_INVALID_DATA
Definition: definitions.h:36
BYROWCOL
@ BYROWCOL
Definition: definitions.h:149
HVEC1_0
@ HVEC1_0
Definition: definitions.h:206
HVEC0
@ HVEC0
Definition: definitions.h:199
NOFIELD
@ NOFIELD
scalar or vector of scalars describe (no true field)
Definition: definitions.h:84