v0.15.0
Loading...
Searching...
No Matches
dynamic_first_order_con_law.cpp File Reference
#include <MoFEM.hpp>
#include <MatrixFunction.hpp>

Go to the source code of this file.

Classes

struct  PostProcEleByDim< 2 >
 
struct  PostProcEleByDim< 3 >
 
struct  OpCalculateFStab< DIM_0, DIM_1 >
 
struct  OpCalculatePiola< DIM_0, DIM_1 >
 
struct  OpCalculateDisplacement< DIM >
 
struct  OpCalculatePiolaIncompressibleNH< DIM_0, DIM_1 >
 
struct  OpCalculateDeformationGradient< DIM_0, DIM_1 >
 
struct  TSPrePostProc
 Set of functions called by PETSc solver used to refine and update mesh. More...
 
struct  LinMomTimeScale
 
struct  CommonData
 
struct  Example
 [Example] More...
 
struct  Example::DynamicFirstOrderConsSinusTimeScale
 
struct  Example::DynamicFirstOrderConsConstantTimeScale
 
struct  Monitor
 [Push operators to pipeline] More...
 

Typedefs

using DomainEle = PipelineManager::ElementsAndOpsByDim<SPACE_DIM>::DomainEle
 
using PostProcEle = PostProcBrokenMeshInMoab<DomainEle>
 
using PostProcFaceEle
 
using BoundaryEle
 
using PostProcEleDomain = PostProcEleByDim<SPACE_DIM>::PostProcEleDomain
 
using SideEle = PostProcEleByDim<SPACE_DIM>::SideEle
 
using PostProcEleBdy = PostProcEleByDim<SPACE_DIM>::PostProcEleBdy
 
using OpMassV
 
using OpMassF
 
using OpInertiaForce
 
using OpBodyForce
 
using DomainNaturalBC
 
using OpBodyForceVector
 
using OpGradTimesTensor2
 
using OpRhsTestPiola
 
using OpGradTimesPiola
 
using BoundaryNaturalBC
 
using OpForce = BoundaryNaturalBC::OpFlux<NaturalForceMeshsets, 1, SPACE_DIM>
 

Functions

template<typename T >
double trace (FTensor::Tensor2< T, 2, 2 > &t_stress)
 
template<typename T >
double trace (FTensor::Tensor2< T, 3, 3 > &t_stress)
 
int main (int argc, char *argv[])
 

Variables

constexpr int SPACE_DIM
 
constexpr double omega = 1.
 Save field DOFS on vertices/tags.
 
constexpr double young_modulus = 1.
 
constexpr double poisson_ratio = 0.
 
double bulk_modulus_K = young_modulus / (3. * (1. - 2. * poisson_ratio))
 
double shear_modulus_G = young_modulus / (2. * (1. + poisson_ratio))
 
double mu = young_modulus / (2. * (1. + poisson_ratio))
 
double lamme_lambda
 
static boost::weak_ptr< TSPrePostProctsPrePostProc
 
static char help [] = "...\n\n"
 [Check]
 

Typedef Documentation

◆ BoundaryEle

◆ BoundaryNaturalBC

◆ DomainEle

◆ DomainNaturalBC

◆ OpBodyForce

using OpBodyForce
Initial value:

Definition at line 64 of file dynamic_first_order_con_law.cpp.

◆ OpBodyForceVector

◆ OpForce

◆ OpGradTimesPiola

Initial value:
FormsIntegrators< DomainEleOp >::Assembly< A >::LinearForm< I >::OpGradTimesTensor< 1, FIELD_DIM, SPACE_DIM > OpGradTimesTensor

Definition at line 81 of file dynamic_first_order_con_law.cpp.

◆ OpGradTimesTensor2

◆ OpInertiaForce

◆ OpMassF

◆ OpMassV

◆ OpRhsTestPiola

Initial value:
IntegrationType::GAUSS>::OpBaseTimesVector<1, SPACE_DIM * SPACE_DIM, 1>
Examples
dynamic_first_order_con_law.cpp.

Definition at line 77 of file dynamic_first_order_con_law.cpp.

◆ PostProcEle

◆ PostProcEleBdy

◆ PostProcEleDomain

◆ PostProcFaceEle

◆ SideEle

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 1190 of file dynamic_first_order_con_law.cpp.

1190 {
1191
1192 // Initialisation of MoFEM/PETSc and MOAB data structures
1193 const char param_file[] = "param_file.petsc";
1194 MoFEM::Core::Initialize(&argc, &argv, param_file, help);
1195
1196 // Add logging channel for example
1197 auto core_log = logging::core::get();
1198 core_log->add_sink(
1200 LogManager::setLog("EXAMPLE");
1201 MOFEM_LOG_TAG("EXAMPLE", "example");
1202
1203 try {
1204
1205 //! [Register MoFEM discrete manager in PETSc]
1206 DMType dm_name = "DMMOFEM";
1207 CHKERR DMRegister_MoFEM(dm_name);
1208 //! [Register MoFEM discrete manager in PETSc
1209
1210 //! [Create MoAB]
1211 moab::Core mb_instance; ///< mesh database
1212 moab::Interface &moab = mb_instance; ///< mesh database interface
1213 //! [Create MoAB]
1214
1215 //! [Create MoFEM]
1216 MoFEM::Core core(moab); ///< finite element database
1217 MoFEM::Interface &m_field = core; ///< finite element database interface
1218 //! [Create MoFEM]
1219
1220 //! [Example]
1221 Example ex(m_field);
1222 CHKERR ex.runProblem();
1223 //! [Example]
1224 }
1226
1228}
#define CATCH_ERRORS
Catch errors.
#define CHKERR
Inline error check.
static char help[]
[Check]
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition DMMoFEM.cpp:43
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
[Example]
Definition plastic.cpp:216
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.
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.

◆ trace() [1/2]

template<typename T >
double trace ( FTensor::Tensor2< T, 2, 2 > & t_stress)
inline
Examples
dynamic_first_order_con_law.cpp.

Definition at line 14 of file dynamic_first_order_con_law.cpp.

14 {
15 constexpr double third = boost::math::constants::third<double>();
16 return (t_stress(0, 0) + t_stress(1, 1));
17};
constexpr double third

◆ trace() [2/2]

template<typename T >
double trace ( FTensor::Tensor2< T, 3, 3 > & t_stress)
inline

Definition at line 19 of file dynamic_first_order_con_law.cpp.

19 {
20 return (t_stress(0, 0) + t_stress(1, 1) + t_stress(2, 2));
21};

Variable Documentation

◆ bulk_modulus_K

◆ help

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

[Check]

Definition at line 1188 of file dynamic_first_order_con_law.cpp.

◆ lamme_lambda

double lamme_lambda
Initial value:
((1. + poisson_ratio) * (1. - 2. * poisson_ratio))
constexpr double poisson_ratio
constexpr double young_modulus
Examples
dynamic_first_order_con_law.cpp.

Definition at line 98 of file dynamic_first_order_con_law.cpp.

◆ mu

◆ omega

double omega = 1.
constexpr

Save field DOFS on vertices/tags.

Examples
dynamic_first_order_con_law.cpp, and shallow_wave.cpp.

Definition at line 92 of file dynamic_first_order_con_law.cpp.

◆ poisson_ratio

double poisson_ratio = 0.
constexpr

Definition at line 94 of file dynamic_first_order_con_law.cpp.

◆ shear_modulus_G

◆ SPACE_DIM

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

Definition at line 23 of file dynamic_first_order_con_law.cpp.

◆ tsPrePostProc

boost::weak_ptr<TSPrePostProc> tsPrePostProc
static

◆ young_modulus

double young_modulus = 1.
constexpr

Definition at line 93 of file dynamic_first_order_con_law.cpp.