v0.13.2
Loading...
Searching...
No Matches
Classes | Typedefs | Functions | Variables
elastic.cpp File Reference

elastic example More...

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

Go to the source code of this file.

Classes

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

Typedefs

using EntData = EntitiesFieldData::EntData
 
using DomainEle = PipelineManager::ElementsAndOpsByDim< SPACE_DIM >::DomainEle
 
using BoundaryEle = PipelineManager::ElementsAndOpsByDim< SPACE_DIM >::BoundaryEle
 
using DomainEleOp = DomainEle::UserDataOperator
 
using BoundaryEleOp = BoundaryEle::UserDataOperator
 
using OpK = FormsIntegrators< DomainEleOp >::Assembly< A >::BiLinearForm< I >::OpGradSymTensorGrad< 1, SPACE_DIM, SPACE_DIM, 0 >
 
using OpInternalForce = FormsIntegrators< DomainEleOp >::Assembly< A >::LinearForm< I >::OpGradTimesSymTensor< 1, SPACE_DIM, SPACE_DIM >
 
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 OpEssentialLhs = EssentialBC< BoundaryEleOp >::Assembly< A >::BiLinearForm< GAUSS >::OpEssentialLhs< DisplacementCubitBcData, 1, SPACE_DIM >
 
using OpEssentialRhs = EssentialBC< BoundaryEleOp >::Assembly< A >::LinearForm< GAUSS >::OpEssentialRhs< DisplacementCubitBcData, 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 SPACE_DIM
 [Define dimension] More...
 
constexpr AssemblyType A
 
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))
 
static char help [] = "...\n\n"
 [Check] More...
 

Detailed Description

elastic example

Version
0.13.2
Date
2022-09-19

Definition in file elastic.cpp.

Typedef Documentation

◆ BoundaryEle

Definition at line 27 of file elastic.cpp.

◆ BoundaryEleOp

Definition at line 30 of file elastic.cpp.

◆ BoundaryLhsBCs

Definition at line 44 of file elastic.cpp.

◆ BoundaryRhsBCs

Definition at line 42 of file elastic.cpp.

◆ DomainEle

Definition at line 26 of file elastic.cpp.

◆ DomainEleOp

Definition at line 29 of file elastic.cpp.

◆ DomainRhsBCs

Definition at line 40 of file elastic.cpp.

◆ EntData

Definition at line 25 of file elastic.cpp.

◆ OpBoundaryLhsBCs

Definition at line 45 of file elastic.cpp.

◆ OpBoundaryRhsBCs

Definition at line 43 of file elastic.cpp.

◆ OpDomainRhsBCs

Definition at line 41 of file elastic.cpp.

◆ OpEssentialLhs

Definition at line 47 of file elastic.cpp.

◆ OpEssentialRhs

Definition at line 49 of file elastic.cpp.

◆ OpInternalForce

using OpInternalForce = FormsIntegrators<DomainEleOp>::Assembly<A>::LinearForm< I>::OpGradTimesSymTensor<1, SPACE_DIM, SPACE_DIM>
Examples
nonlinear_elastic.cpp.

Definition at line 34 of file elastic.cpp.

◆ OpK

using OpK = FormsIntegrators<DomainEleOp>::Assembly<A>::BiLinearForm< I>::OpGradSymTensorGrad<1, SPACE_DIM, SPACE_DIM, 0>

Definition at line 32 of file elastic.cpp.

◆ PostProcEleBdy

Definition at line 68 of file elastic.cpp.

◆ PostProcEleDomain

Definition at line 66 of file elastic.cpp.

◆ SideEle

Definition at line 67 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 694 of file elastic.cpp.

694 {
695
696 // Initialisation of MoFEM/PETSc and MOAB data structures
697 const char param_file[] = "param_file.petsc";
699
700 auto core_log = logging::core::get();
701 core_log->add_sink(
703
704 try {
705
706 //! [Register MoFEM discrete manager in PETSc]
707 DMType dm_name = "DMMOFEM";
708 CHKERR DMRegister_MoFEM(dm_name);
709 //! [Register MoFEM discrete manager in PETSc
710
711 //! [Create MoAB]
712 moab::Core mb_instance; ///< mesh database
713 moab::Interface &moab = mb_instance; ///< mesh database interface
714 //! [Create MoAB]
715
716 //! [Create MoFEM]
717 MoFEM::Core core(moab); ///< finite element database
718 MoFEM::Interface &m_field = core; ///< finite element database interface
719 //! [Create MoFEM]
720
721 //! [Example]
722 Example ex(m_field);
723 CHKERR ex.runProblem();
724 //! [Example]
725 }
727
729}
std::string param_file
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
#define CHKERR
Inline error check.
Definition: definitions.h:535
static char help[]
[Check]
Definition: elastic.cpp:692
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:47
[Example]
Definition: plastic.cpp:141
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.
Definition: LogManager.cpp:298
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
Definition: LogManager.cpp:344

Variable Documentation

◆ A

constexpr AssemblyType A
constexpr
Initial value:
= (SCHUR_ASSEMBLE)
? AssemblyType::SCHUR
: AssemblyType::PETSC

Definition at line 18 of file elastic.cpp.

◆ bulk_modulus_K

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

◆ help

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

[Check]

Definition at line 692 of file elastic.cpp.

◆ I

constexpr IntegrationType I
constexpr
Initial value:
=
IntegrationType::GAUSS

Definition at line 22 of file elastic.cpp.

◆ poisson_ratio

constexpr double poisson_ratio = 0.3
constexpr

Definition at line 77 of file elastic.cpp.

◆ shear_modulus_G

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

◆ SPACE_DIM

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

[Define dimension]

Definition at line 16 of file elastic.cpp.

◆ young_modulus

constexpr double young_modulus = 1
constexpr

Definition at line 76 of file elastic.cpp.