v0.14.0
Loading...
Searching...
No Matches
elastic.cpp File Reference

elastic example More...

#include <MoFEM.hpp>
#include <ElasticSpring.hpp>
#include <FluidLevel.hpp>
#include <CalculateTraction.hpp>
#include <NaturalDomainBC.hpp>
#include <NaturalBoundaryBC.hpp>

Go to the source code of this file.

Classes

struct  DomainBCs
 [OpInternalForce] More...
 
struct  BoundaryBCs
 
struct  PostProcEleByDim< 2 >
 
struct  PostProcEleByDim< 3 >
 
struct  Example
 [Example] More...
 
struct  SetUpSchur
 [Push operators to pipeline] More...
 
struct  SetUpSchurImpl
 

Typedefs

using DomainEle = PipelineManager::ElementsAndOpsByDim<SPACE_DIM>::DomainEle
 
using BoundaryEle
 
using OpK
 [Define entities]
 
using OpInternalForce
 
using DomainRhsBCs = NaturalBC<DomainEleOp>::Assembly<A>::LinearForm<I>
 
using OpDomainRhsBCs = DomainRhsBCs::OpFlux<DomainBCs, 1, SPACE_DIM>
 
using BoundaryRhsBCs = NaturalBC<BoundaryEleOp>::Assembly<A>::LinearForm<I>
 
using OpBoundaryRhsBCs = BoundaryRhsBCs::OpFlux<BoundaryBCs, 1, SPACE_DIM>
 
using BoundaryLhsBCs = NaturalBC<BoundaryEleOp>::Assembly<A>::BiLinearForm<I>
 
using OpBoundaryLhsBCs = BoundaryLhsBCs::OpFlux<BoundaryBCs, 1, SPACE_DIM>
 
using PostProcEleDomain = PostProcEleByDim<SPACE_DIM>::PostProcEleDomain
 
using SideEle = PostProcEleByDim<SPACE_DIM>::SideEle
 
using PostProcEleBdy = PostProcEleByDim<SPACE_DIM>::PostProcEleBdy
 

Functions

int main (int argc, char *argv[])
 

Variables

constexpr int BASE_DIM = 1
 
constexpr int SPACE_DIM
 [Define dimension]
 
constexpr AssemblyType A
 [Define dimension]
 
constexpr IntegrationType I
 
constexpr double young_modulus = 1
 
constexpr double poisson_ratio = 0.3
 
constexpr double bulk_modulus_K = young_modulus / (3 * (1 - 2 * poisson_ratio))
 
constexpr double shear_modulus_G = young_modulus / (2 * (1 + poisson_ratio))
 
PetscBool is_plane_strain = PETSC_FALSE
 
static char help [] = "...\n\n"
 [Check]
 

Detailed Description

elastic example

Version
0.13.2
Date
2022-09-19

Definition in file elastic.cpp.

Typedef Documentation

◆ BoundaryEle

◆ BoundaryLhsBCs

Definition at line 52 of file elastic.cpp.

◆ BoundaryRhsBCs

Definition at line 50 of file elastic.cpp.

◆ DomainEle

◆ DomainRhsBCs

Definition at line 48 of file elastic.cpp.

◆ OpBoundaryLhsBCs

Definition at line 53 of file elastic.cpp.

◆ OpBoundaryRhsBCs

Definition at line 51 of file elastic.cpp.

◆ OpDomainRhsBCs

Definition at line 49 of file elastic.cpp.

◆ OpInternalForce

Initial value:
I>::OpGradTimesSymTensor<BASE_DIM, SPACE_DIM, SPACE_DIM>
constexpr IntegrationType I
Definition elastic.cpp:25

[OpK] [OpInternalForce]

Definition at line 42 of file elastic.cpp.

◆ OpK

using OpK
Initial value:
I>::OpGradSymTensorGrad<BASE_DIM, SPACE_DIM, SPACE_DIM, 0>

[Define entities]

[OpK]

Examples
PoissonOperators.hpp, eigen_elastic.cpp, and simple_elasticity.cpp.

Definition at line 38 of file elastic.cpp.

◆ PostProcEleBdy

◆ PostProcEleDomain

◆ SideEle

Definition at line 70 of file elastic.cpp.

Function Documentation

◆ main()

int main ( int argc,
char * argv[] )

[Register MoFEM discrete manager in PETSc]

[Register MoFEM discrete manager in PETSc

[Create MoAB]

< mesh database

< mesh database interface

[Create MoAB]

[Create MoFEM]

< finite element database

< finite element database interface

[Create MoFEM]

[Example]

[Example]

Definition at line 881 of file elastic.cpp.

881 {
882
883 // Initialisation of MoFEM/PETSc and MOAB data structures
884 const char param_file[] = "param_file.petsc";
885 MoFEM::Core::Initialize(&argc, &argv, param_file, help);
886
887 auto core_log = logging::core::get();
888 core_log->add_sink(
890
891 core_log->add_sink(
893 LogManager::setLog("FieldEvaluator");
894 MOFEM_LOG_TAG("FieldEvaluator", "field_eval");
895
896 try {
897
898 //! [Register MoFEM discrete manager in PETSc]
899 DMType dm_name = "DMMOFEM";
900 CHKERR DMRegister_MoFEM(dm_name);
901 DMType dm_name_mg = "DMMOFEM_MG";
903 //! [Register MoFEM discrete manager in PETSc
904
905 //! [Create MoAB]
906 moab::Core mb_instance; ///< mesh database
907 moab::Interface &moab = mb_instance; ///< mesh database interface
908 //! [Create MoAB]
909
910 //! [Create MoFEM]
911 MoFEM::Core core(moab); ///< finite element database
912 MoFEM::Interface &m_field = core; ///< finite element database interface
913 //! [Create MoFEM]
914
915 //! [Example]
916 Example ex(m_field);
917 CHKERR ex.runProblem();
918 //! [Example]
919 }
921
923}
#define CATCH_ERRORS
Catch errors.
#define CHKERR
Inline error check.
static char help[]
[Check]
Definition elastic.cpp:879
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.
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
[Example]
Definition plastic.cpp:213
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:112
Deprecated interface functions.
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
static boost::shared_ptr< std::ostream > getStrmSync()
Get the strm sync object.

Variable Documentation

◆ A

AssemblyType A
constexpr
Initial value:
? AssemblyType::BLOCK_SCHUR
#define SCHUR_ASSEMBLE
Definition contact.cpp:18
AssemblyType
[Storage and set boundary conditions]

[Define dimension]

Definition at line 21 of file elastic.cpp.

◆ BASE_DIM

int BASE_DIM = 1
constexpr

Definition at line 15 of file elastic.cpp.

◆ bulk_modulus_K

double bulk_modulus_K = young_modulus / (3 * (1 - 2 * poisson_ratio))
constexpr

Definition at line 81 of file elastic.cpp.

◆ help

char help[] = "...\n\n"
static

[Check]

Definition at line 879 of file elastic.cpp.

◆ I

IntegrationType I
constexpr
Initial value:
=
IntegrationType::GAUSS

Definition at line 25 of file elastic.cpp.

◆ is_plane_strain

PetscBool is_plane_strain = PETSC_FALSE

Definition at line 84 of file elastic.cpp.

◆ poisson_ratio

double poisson_ratio = 0.3
constexpr

Definition at line 80 of file elastic.cpp.

◆ shear_modulus_G

double shear_modulus_G = young_modulus / (2 * (1 + poisson_ratio))
constexpr

Definition at line 82 of file elastic.cpp.

◆ SPACE_DIM

int SPACE_DIM
constexpr
Initial value:
=
#define EXECUTABLE_DIMENSION
Definition plastic.cpp:13

[Define dimension]

Definition at line 18 of file elastic.cpp.

◆ young_modulus

double young_modulus = 1
constexpr

Definition at line 79 of file elastic.cpp.