v0.15.0
Loading...
Searching...
No Matches
mofem/tutorials/vec-0/elastic.cpp

elastic example

elastic example

Version
0.15.0
Date
2025-10-07
/**
* @file elastic.cpp
* @example mofem/tutorials/vec-0/elastic.cpp
* @brief elastic example
* @version 0.15.0
* @date 2025-10-07
*
* @copyright Anonymous authors (c) 2025 under the MIT license (see LICENSE for
* details)
*/
#include <MoFEM.hpp>
using namespace MoFEM;
constexpr int BASE_DIM = 1; //< Dimension of the base functions
//! [Define dimension]
constexpr int SPACE_DIM =
EXECUTABLE_DIMENSION; //< Space dimension of problem, mesh
//! [Define dimension]
? AssemblyType::BLOCK_SCHUR
: AssemblyType::PETSC; //< selected assembly type
constexpr IntegrationType I =
IntegrationType::GAUSS; //< selected integration type
//! [Define entities]
using BoundaryEle =
using DomainEleOp = DomainEle::UserDataOperator;
using BoundaryEleOp = BoundaryEle::UserDataOperator;
//! [Define entities]
struct DomainBCs {};
struct BoundaryBCs {};
using OpDomainRhsBCs = DomainRhsBCs::OpFlux<DomainBCs, 1, SPACE_DIM>;
using OpBoundaryRhsBCs = BoundaryRhsBCs::OpFlux<BoundaryBCs, 1, SPACE_DIM>;
using OpBoundaryLhsBCs = BoundaryLhsBCs::OpFlux<BoundaryBCs, 1, SPACE_DIM>;
template <int DIM> struct PostProcEleByDim;
template <> struct PostProcEleByDim<2> {
};
template <> struct PostProcEleByDim<3> {
};
static char help[] = "...\n\n";
int main(int argc, char *argv[]) {
// Initialisation of MoFEM/PETSc and MOAB data structures
const char param_file[] = "param_file.petsc";
MoFEM::Core::Initialize(&argc, &argv, param_file, help);
auto core_log = logging::core::get();
core_log->add_sink(
LogManager::createSink(LogManager::getStrmWorld(), "TIMER"));
core_log->add_sink(
LogManager::createSink(LogManager::getStrmSync(), "FieldEvaluator"));
LogManager::setLog("FieldEvaluator");
MOFEM_LOG_TAG("FieldEvaluator", "field_eval");
try {
//! [Register MoFEM discrete manager in PETSc]
DMType dm_name = "DMMOFEM";
DMType dm_name_mg = "DMMOFEM_MG";
//! [Register MoFEM discrete manager in PETSc
//! [Create MoAB]
moab::Core mb_instance; ///< mesh database
moab::Interface &moab = mb_instance; ///< mesh database interface
//! [Create MoAB]
//! [Create MoFEM]
MoFEM::Core core(moab); ///< finite element database
MoFEM::Interface &m_field = core; ///< finite element database interface
//! [Create MoFEM]
//! [Example]
ElasticExample ex(m_field);
CHKERR ex.runProblem();
//! [Example]
}
}
Implementation of elastic example class.
static char help[]
int main()
constexpr int SPACE_DIM
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
#define SCHUR_ASSEMBLE
Definition contact.cpp:18
#define CATCH_ERRORS
Catch errors.
#define CHKERR
Inline error check.
constexpr int BASE_DIM
PostProcEleByDim< SPACE_DIM >::PostProcEleDomain PostProcEleDomain
PostProcEleByDim< SPACE_DIM >::PostProcEleBdy PostProcEleBdy
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition DMMoFEM.cpp:43
MoFEMErrorCode DMRegister_MGViaApproxOrders(const char sname[])
Register DM for Multi-Grid via approximation orders.
IntegrationType
Form integrator integration types.
AssemblyType
[Storage and set boundary conditions]
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
constexpr IntegrationType I
constexpr AssemblyType A
NaturalBC< BoundaryEleOp >::Assembly< AT >::BiLinearForm< IT > BoundaryLhsBCs
Definition plastic.cpp:175
NaturalBC< BoundaryEleOp >::Assembly< AT >::LinearForm< IT > BoundaryRhsBCs
Definition plastic.cpp:172
#define EXECUTABLE_DIMENSION
Definition plastic.cpp:13
BoundaryLhsBCs::OpFlux< PlasticOps::BoundaryBCs, 1, SPACE_DIM > OpBoundaryLhsBCs
Definition plastic.cpp:177
NaturalBC< DomainEleOp >::Assembly< AT >::LinearForm< IT > DomainRhsBCs
Definition plastic.cpp:169
DomainRhsBCs::OpFlux< PlasticOps::DomainBCs, 1, SPACE_DIM > OpDomainRhsBCs
Definition plastic.cpp:171
BoundaryRhsBCs::OpFlux< PlasticOps::BoundaryBCs, 1, SPACE_DIM > OpBoundaryRhsBCs
Definition plastic.cpp:174
ElementsAndOps< SPACE_DIM >::SideEle SideEle
Definition plastic.cpp:61
Boundary conditions marker.
Definition elastic.cpp:39
[Define entities]
Definition elastic.cpp:38
Core (interface) class.
Definition Core.hpp:82
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
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition Core.cpp:118
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
Assembly methods.
Definition Natural.hpp:65
Template struct for dimension-specific finite element types.