v0.13.1
Classes | Macros | Typedefs | Functions | Variables
plastic.cpp File Reference
#include <MoFEM.hpp>
#include <MatrixFunction.hpp>
#include <IntegrationRules.hpp>
#include <HenckyOps.hpp>
#include <PlasticOps.hpp>
#include <OpPostProcElastic.hpp>

Go to the source code of this file.

Classes

struct  ElementsAndOps< DIM >
 
struct  ElementsAndOps< 2 >
 
struct  ElementsAndOps< 3 >
 
struct  Example
 [Example] More...
 

Macros

#define EXECUTABLE_DIMENSION   3
 

Typedefs

using EntData = EntitiesFieldData::EntData
 
using DomainEle = ElementsAndOps< SPACE_DIM >::DomainEle
 
using DomainEleOp = ElementsAndOps< SPACE_DIM >::DomainEleOp
 
using BoundaryEle = ElementsAndOps< SPACE_DIM >::BoundaryEle
 
using BoundaryEleOp = ElementsAndOps< SPACE_DIM >::BoundaryEleOp
 
using PostProcEle = ElementsAndOps< SPACE_DIM >::PostProcEle
 
using AssemblyDomainEleOp = FormsIntegrators< DomainEleOp >::Assembly< PETSC >::OpBase
 
using OpBodyForce = FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< 1, SPACE_DIM >
 [Body force] More...
 
using OpKCauchy = FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradSymTensorGrad< 1, SPACE_DIM, SPACE_DIM, 0 >
 [Body force] More...
 
using OpInternalForceCauchy = FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesSymTensor< 1, SPACE_DIM, SPACE_DIM >
 
using OpMass = FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM >
 [Only used with Hooke equation (linear material model)] More...
 
using OpInertiaForce = FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesVector< 1, SPACE_DIM, 1 >
 
using OpKPiola = FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradTensorGrad< 1, SPACE_DIM, SPACE_DIM, 1 >
 [Only used for dynamics] More...
 
using OpInternalForcePiola = FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesTensor< 1, SPACE_DIM, SPACE_DIM >
 
using OpBoundaryMass = FormsIntegrators< BoundaryEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM >
 [Only used with Hencky/nonlinear material] More...
 
using OpBoundaryVec = FormsIntegrators< BoundaryEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesVector< 1, SPACE_DIM, 0 >
 
using OpBoundaryInternal = FormsIntegrators< BoundaryEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesVector< 1, SPACE_DIM, 1 >
 
using OpScaleL2 = MoFEM::OpScaleBaseBySpaceInverseOfMeasure< DomainEleOp >
 [Essential boundary conditions] More...
 

Functions

long double hardening (long double tau, double temp)
 
long double hardening_dtau (long double tau, double temp)
 
int main (int argc, char *argv[])
 

Variables

constexpr int SPACE_DIM
 
constexpr EntityType boundary_ent = SPACE_DIM == 3 ? MBTRI : MBEDGE
 
PetscBool is_large_strains = PETSC_TRUE
 
double scale = 1.
 
double young_modulus = 206913
 
double poisson_ratio = 0.29
 
double rho = 0
 
double sigmaY = 450
 
double H = 129
 
double visH = 0
 
double cn = 1
 
double Qinf = 265
 
double b_iso = 16.93
 
int order = 2
 
static char help [] = "...\n\n"
 [Solve] More...
 

Macro Definition Documentation

◆ EXECUTABLE_DIMENSION

#define EXECUTABLE_DIMENSION   3

Typedef Documentation

◆ AssemblyDomainEleOp

◆ BoundaryEle

Definition at line 58 of file plastic.cpp.

◆ BoundaryEleOp

Definition at line 59 of file plastic.cpp.

◆ DomainEle

Definition at line 56 of file plastic.cpp.

◆ DomainEleOp

Definition at line 57 of file plastic.cpp.

◆ EntData

Examples
plastic.cpp.

Definition at line 55 of file plastic.cpp.

◆ OpBodyForce

using OpBodyForce = FormsIntegrators<DomainEleOp>::Assembly<PETSC>::LinearForm< GAUSS>::OpSource<1, SPACE_DIM>

[Body force]

Examples
plastic.cpp.

Definition at line 66 of file plastic.cpp.

◆ OpBoundaryInternal

using OpBoundaryInternal = FormsIntegrators<BoundaryEleOp>::Assembly< PETSC>::LinearForm<GAUSS>::OpBaseTimesVector<1, SPACE_DIM, 1>
Examples
plastic.cpp.

Definition at line 96 of file plastic.cpp.

◆ OpBoundaryMass

[Only used with Hencky/nonlinear material]

[Essential boundary conditions]

Examples
plastic.cpp.

Definition at line 92 of file plastic.cpp.

◆ OpBoundaryVec

using OpBoundaryVec = FormsIntegrators<BoundaryEleOp>::Assembly< PETSC>::LinearForm<GAUSS>::OpBaseTimesVector<1, SPACE_DIM, 0>
Examples
plastic.cpp.

Definition at line 94 of file plastic.cpp.

◆ OpInertiaForce

using OpInertiaForce = FormsIntegrators<DomainEleOp>::Assembly< PETSC>::LinearForm<GAUSS>::OpBaseTimesVector<1, SPACE_DIM, 1>

Definition at line 80 of file plastic.cpp.

◆ OpInternalForceCauchy

using OpInternalForceCauchy = FormsIntegrators<DomainEleOp>::Assembly< PETSC>::LinearForm<GAUSS>::OpGradTimesSymTensor<1, SPACE_DIM, SPACE_DIM>
Examples
plastic.cpp.

Definition at line 73 of file plastic.cpp.

◆ OpInternalForcePiola

using OpInternalForcePiola = FormsIntegrators<DomainEleOp>::Assembly< PETSC>::LinearForm<GAUSS>::OpGradTimesTensor<1, SPACE_DIM, SPACE_DIM>
Examples
plastic.cpp.

Definition at line 87 of file plastic.cpp.

◆ OpKCauchy

using OpKCauchy = FormsIntegrators<DomainEleOp>::Assembly<PETSC>::BiLinearForm< GAUSS>::OpGradSymTensorGrad<1, SPACE_DIM, SPACE_DIM, 0>

[Body force]

[Only used with Hooke equation (linear material model)]

Examples
plastic.cpp.

Definition at line 71 of file plastic.cpp.

◆ OpKPiola

using OpKPiola = FormsIntegrators<DomainEleOp>::Assembly<PETSC>::BiLinearForm< GAUSS>::OpGradTensorGrad<1, SPACE_DIM, SPACE_DIM, 1>

[Only used for dynamics]

[Only used with Hencky/nonlinear material]

Examples
plastic.cpp.

Definition at line 85 of file plastic.cpp.

◆ OpMass

using OpMass = FormsIntegrators<DomainEleOp>::Assembly<PETSC>::BiLinearForm< GAUSS>::OpMass<1, SPACE_DIM>

[Only used with Hooke equation (linear material model)]

[Only used for dynamics]

Definition at line 78 of file plastic.cpp.

◆ OpScaleL2

[Essential boundary conditions]

Definition at line 99 of file plastic.cpp.

◆ PostProcEle

Definition at line 60 of file plastic.cpp.

Function Documentation

◆ hardening()

long double hardening ( long double  tau,
double  temp 
)
Examples
plastic.cpp.

Definition at line 116 of file plastic.cpp.

116 {
117 return H * tau + Qinf * (1. - std::exp(-b_iso * tau)) + sigmaY;
118}
double Qinf
Definition: plastic.cpp:112
double H
Definition: plastic.cpp:109
double b_iso
Definition: plastic.cpp:113
double sigmaY
Definition: plastic.cpp:108

◆ hardening_dtau()

long double hardening_dtau ( long double  tau,
double  temp 
)
Examples
plastic.cpp.

Definition at line 120 of file plastic.cpp.

120 {
121 return H + Qinf * b_iso * std::exp(-b_iso * tau);
122}

◆ main()

int main ( int  argc,
char *  argv[] 
)
Examples
plastic.cpp.

Definition at line 969 of file plastic.cpp.

969 {
970
971 // Initialisation of MoFEM/PETSc and MOAB data structures
972 const char param_file[] = "param_file.petsc";
974
975 // Add logging channel for example
976 auto core_log = logging::core::get();
977 core_log->add_sink(
978 LogManager::createSink(LogManager::getStrmWorld(), "EXAMPLE"));
979 LogManager::setLog("EXAMPLE");
980 MOFEM_LOG_TAG("EXAMPLE", "example");
981
982 try {
983
984 //! [Register MoFEM discrete manager in PETSc]
985 DMType dm_name = "DMMOFEM";
986 CHKERR DMRegister_MoFEM(dm_name);
987 //! [Register MoFEM discrete manager in PETSc
988
989 //! [Create MoAB]
990 moab::Core mb_instance; ///< mesh database
991 moab::Interface &moab = mb_instance; ///< mesh database interface
992 //! [Create MoAB]
993
994 //! [Create MoFEM]
995 MoFEM::Core core(moab); ///< finite element database
996 MoFEM::Interface &m_field = core; ///< finite element database insterface
997 //! [Create MoFEM]
998
999 //! [Load mesh]
1000 Simple *simple = m_field.getInterface<Simple>();
1001 CHKERR simple->getOptions();
1002 CHKERR simple->loadFile();
1003 //! [Load mesh]
1004
1005 //! [Example]
1006 Example ex(m_field);
1007 CHKERR ex.runProblem();
1008 //! [Example]
1009 }
1011
1013}
std::string param_file
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition: acoustic.cpp:69
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
#define CHKERR
Inline error check.
Definition: definitions.h:548
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMMoFEM.cpp:59
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:342
CoreTmp< 0 > Core
Definition: Core.hpp:1096
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1965
static char help[]
[Solve]
Definition: plastic.cpp:967
[Example]
Definition: plastic.cpp:130
Core (interface) class.
Definition: Core.hpp:92
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:85
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:125
Deprecated interface functions.
Simple interface for fast problem set-up.
Definition: Simple.hpp:35
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.

Variable Documentation

◆ b_iso

double b_iso = 16.93
Examples
plastic.cpp.

Definition at line 113 of file plastic.cpp.

◆ boundary_ent

constexpr EntityType boundary_ent = SPACE_DIM == 3 ? MBTRI : MBEDGE
constexpr
Examples
plastic.cpp.

Definition at line 54 of file plastic.cpp.

◆ cn

double cn = 1
Examples
ContactOps.hpp, PlasticOps.hpp, and plastic.cpp.

Definition at line 111 of file plastic.cpp.

◆ H

double H = 129
Examples
plastic.cpp.

Definition at line 109 of file plastic.cpp.

◆ help

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

[Solve]

Examples
plastic.cpp.

Definition at line 967 of file plastic.cpp.

◆ is_large_strains

PetscBool is_large_strains = PETSC_TRUE
Examples
plastic.cpp.

Definition at line 101 of file plastic.cpp.

◆ order

int order = 2
Examples
plastic.cpp.

Definition at line 114 of file plastic.cpp.

◆ poisson_ratio

double poisson_ratio = 0.29
Examples
plastic.cpp.

Definition at line 106 of file plastic.cpp.

◆ Qinf

double Qinf = 265
Examples
plastic.cpp.

Definition at line 112 of file plastic.cpp.

◆ rho

double rho = 0

◆ scale

double scale = 1.

◆ sigmaY

double sigmaY = 450
Examples
PlasticOps.hpp, and plastic.cpp.

Definition at line 108 of file plastic.cpp.

◆ SPACE_DIM

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

Definition at line 51 of file plastic.cpp.

◆ visH

double visH = 0
Examples
PlasticOps.hpp, and plastic.cpp.

Definition at line 110 of file plastic.cpp.

◆ young_modulus

double young_modulus = 206913
Examples
Remodeling.cpp, and plastic.cpp.

Definition at line 105 of file plastic.cpp.