v0.13.2
Loading...
Searching...
No Matches
Classes | Macros | Typedefs | Functions | Variables
plastic.cpp File Reference
#include <MoFEM.hpp>
#include <MatrixFunction.hpp>
#include <IntegrationRules.hpp>
#include <HenckyOps.hpp>
#include <PlasticOps.hpp>

Go to the source code of this file.

Classes

struct  ElementsAndOps< DIM >
 
struct  ElementsAndOps< 2 >
 
struct  ElementsAndOps< 3 >
 
struct  Example
 [Example] More...
 
struct  Example::PlasticityTimeScale
 
struct  SetUpSchur
 [Push operators to pipeline] More...
 
struct  SetUpSchurImpl
 

Macros

#define EXECUTABLE_DIMENSION   3
 

Typedefs

using EntData = EntitiesFieldData::EntData
 
using DomainEle = ElementsAndOps< SPACE_DIM >::DomainEle
 
using DomainEleOp = DomainEle::UserDataOperator
 
using BoundaryEle = ElementsAndOps< SPACE_DIM >::BoundaryEle
 
using BoundaryEleOp = BoundaryEle::UserDataOperator
 
using PostProcEle = PostProcBrokenMeshInMoab< DomainEle >
 
using AssemblyDomainEleOp = FormsIntegrators< DomainEleOp >::Assembly< A >::OpBase
 
using OpKCauchy = FormsIntegrators< DomainEleOp >::Assembly< A >::BiLinearForm< GAUSS >::OpGradSymTensorGrad< 1, SPACE_DIM, SPACE_DIM, 0 >
 [Only used with Hooke equation (linear material model)] More...
 
using OpInternalForceCauchy = FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< G >::OpGradTimesSymTensor< 1, SPACE_DIM, SPACE_DIM >
 
using OpMass = FormsIntegrators< DomainEleOp >::Assembly< A >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM >
 [Only used with Hooke equation (linear material model)] More...
 
using OpInertiaForce = FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< G >::OpBaseTimesVector< 1, SPACE_DIM, 1 >
 
using OpKPiola = FormsIntegrators< DomainEleOp >::Assembly< A >::BiLinearForm< GAUSS >::OpGradTensorGrad< 1, SPACE_DIM, SPACE_DIM, 1 >
 [Only used for dynamics] More...
 
using OpInternalForcePiola = FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< G >::OpGradTimesTensor< 1, SPACE_DIM, SPACE_DIM >
 
using OpBoundaryMass = FormsIntegrators< BoundaryEleOp >::Assembly< PETSC >::BiLinearForm< G >::OpMass< 1, SPACE_DIM >
 [Only used with Hencky/nonlinear material] More...
 
using OpBoundaryVec = FormsIntegrators< BoundaryEleOp >::Assembly< PETSC >::LinearForm< G >::OpBaseTimesVector< 1, SPACE_DIM, 0 >
 
using OpBoundaryInternal = FormsIntegrators< BoundaryEleOp >::Assembly< PETSC >::LinearForm< G >::OpBaseTimesVector< 1, SPACE_DIM, 1 >
 
using OpScaleL2 = MoFEM::OpScaleBaseBySpaceInverseOfMeasure< DomainEleOp >
 [Essential boundary conditions] More...
 
using DomainNaturalBC = NaturalBC< DomainEleOp >::Assembly< A >::LinearForm< G >
 
using OpBodyForce = DomainNaturalBC::OpFlux< NaturalMeshsetType< BLOCKSET >, 1, SPACE_DIM >
 
using BoundaryNaturalBC = NaturalBC< BoundaryEleOp >::Assembly< A >::LinearForm< G >
 
using OpForce = BoundaryNaturalBC::OpFlux< NaturalMeshsetType< BLOCKSET >, 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 >
 

Functions

double hardening_exp (double tau)
 
double hardening (double tau)
 
double hardening_dtau (double tau)
 
double hardening_dtau2 (double tau)
 
int main (int argc, char *argv[])
 

Variables

constexpr int SPACE_DIM
 
constexpr auto size_symm = (SPACE_DIM * (SPACE_DIM + 1)) / 2
 
constexpr AssemblyType A
 
constexpr IntegrationType G
 
PetscBool is_large_strains = PETSC_TRUE
 
PetscBool set_timer = PETSC_FALSE
 
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 cn0 = 1
 
double cn1 = 1
 
double zeta = 5e-2
 
double Qinf = 265
 
double b_iso = 16.93
 
int order = 2
 Order if fixed. More...
 
int geom_order = 2
 Order if fixed. More...
 
constexpr size_t activ_history_sise = 1
 
static char help [] = "...\n\n"
 [Solve] More...
 

Macro Definition Documentation

◆ EXECUTABLE_DIMENSION

#define EXECUTABLE_DIMENSION   3

Typedef Documentation

◆ AssemblyDomainEleOp

◆ BoundaryEle

Definition at line 44 of file plastic.cpp.

◆ BoundaryEleOp

Definition at line 45 of file plastic.cpp.

◆ BoundaryNaturalBC

Definition at line 85 of file plastic.cpp.

◆ DomainEle

Definition at line 42 of file plastic.cpp.

◆ DomainEleOp

Definition at line 43 of file plastic.cpp.

◆ DomainNaturalBC

Definition at line 81 of file plastic.cpp.

◆ EntData

Definition at line 41 of file plastic.cpp.

◆ OpBodyForce

Definition at line 82 of file plastic.cpp.

◆ OpBoundaryInternal

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

Definition at line 76 of file plastic.cpp.

◆ OpBoundaryMass

[Only used with Hencky/nonlinear material]

[Essential boundary conditions]

Examples
contact.cpp, and helmholtz.cpp.

Definition at line 72 of file plastic.cpp.

◆ OpBoundaryVec

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

Definition at line 74 of file plastic.cpp.

◆ OpEssentialLhs

Definition at line 89 of file plastic.cpp.

◆ OpEssentialRhs

Definition at line 91 of file plastic.cpp.

◆ OpForce

Definition at line 86 of file plastic.cpp.

◆ OpInertiaForce

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

Definition at line 60 of file plastic.cpp.

◆ OpInternalForceCauchy

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

Definition at line 53 of file plastic.cpp.

◆ OpInternalForcePiola

Examples
contact.cpp, and plastic.cpp.

Definition at line 67 of file plastic.cpp.

◆ OpKCauchy

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

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

Examples
contact.cpp, and plastic.cpp.

Definition at line 51 of file plastic.cpp.

◆ OpKPiola

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

[Only used for dynamics]

[Only used with Hencky/nonlinear material]

Examples
contact.cpp, and plastic.cpp.

Definition at line 65 of file plastic.cpp.

◆ OpMass

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

[Only used for dynamics]

Definition at line 58 of file plastic.cpp.

◆ OpScaleL2

[Essential boundary conditions]

Definition at line 79 of file plastic.cpp.

◆ PostProcEle

Definition at line 46 of file plastic.cpp.

Function Documentation

◆ hardening()

double hardening ( double  tau)
inline
Examples
PlasticOpsGeneric.hpp, and plastic.cpp.

Definition at line 122 of file plastic.cpp.

122 {
123 return H * tau + Qinf * (1. - hardening_exp(tau)) + sigmaY;
124}
double Qinf
Definition: plastic.cpp:108
double H
Definition: plastic.cpp:103
double sigmaY
Definition: plastic.cpp:102
double hardening_exp(double tau)
Definition: plastic.cpp:116

◆ hardening_dtau()

double hardening_dtau ( double  tau)
inline
Examples
PlasticOpsGeneric.hpp, and plastic.cpp.

Definition at line 126 of file plastic.cpp.

126 {
127 auto r = [](auto tau) { return H + Qinf * b_iso * hardening_exp(tau); };
128 constexpr double eps = 1e-12;
129 return std::max(r(tau), eps * r(0));
130}
static const double eps
const double r
rate factor
double b_iso
Definition: plastic.cpp:109

◆ hardening_dtau2()

double hardening_dtau2 ( double  tau)
inline
Examples
plastic.cpp.

Definition at line 132 of file plastic.cpp.

132 {
133 return -(Qinf * (b_iso * b_iso)) * hardening_exp(tau);
134}

◆ hardening_exp()

double hardening_exp ( double  tau)
inline
Examples
plastic.cpp.

Definition at line 116 of file plastic.cpp.

116 {
117 return std::exp(
118 std::max(static_cast<double>(std::numeric_limits<float>::min_exponent10),
119 -b_iso * tau));
120}

◆ 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]

[Load mesh]

[Load mesh]

[Example]

[Example]

Definition at line 1186 of file plastic.cpp.

1186 {
1187
1188 // Initialisation of MoFEM/PETSc and MOAB data structures
1189 const char param_file[] = "param_file.petsc";
1190 MoFEM::Core::Initialize(&argc, &argv, param_file, help);
1191
1192 // Add logging channel for example
1193 auto core_log = logging::core::get();
1194 core_log->add_sink(
1196 core_log->add_sink(
1198 LogManager::setLog("EXAMPLE");
1199
1200 MOFEM_LOG_TAG("EXAMPLE", "example");
1201
1202 try {
1203
1204 //! [Register MoFEM discrete manager in PETSc]
1205 DMType dm_name = "DMMOFEM";
1206 CHKERR DMRegister_MoFEM(dm_name);
1207 //! [Register MoFEM discrete manager in PETSc
1208
1209 //! [Create MoAB]
1210 moab::Core mb_instance; ///< mesh database
1211 moab::Interface &moab = mb_instance; ///< mesh database interface
1212 //! [Create MoAB]
1213
1214 //! [Create MoFEM]
1215 MoFEM::Core core(moab); ///< finite element database
1216 MoFEM::Interface &m_field = core; ///< finite element database interface
1217 //! [Create MoFEM]
1218
1219 //! [Load mesh]
1220 Simple *simple = m_field.getInterface<Simple>();
1221 CHKERR simple->getOptions();
1222 CHKERR simple->loadFile();
1223 //! [Load mesh]
1224
1225 //! [Example]
1226 Example ex(m_field);
1227 CHKERR ex.runProblem();
1228 //! [Example]
1229 }
1231
1233}
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:372
#define CHKERR
Inline error check.
Definition: definitions.h:535
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:47
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
Definition: LogManager.cpp:391
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
static char help[]
[Solve]
Definition: plastic.cpp:1184
[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:300
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
Definition: LogManager.cpp:346
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.

Variable Documentation

◆ A

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

Definition at line 35 of file plastic.cpp.

◆ activ_history_sise

constexpr size_t activ_history_sise = 1
constexpr
Examples
plastic.cpp.

Definition at line 114 of file plastic.cpp.

◆ b_iso

double b_iso = 16.93
Examples
plastic.cpp.

Definition at line 109 of file plastic.cpp.

◆ cn0

double cn0 = 1
Examples
PlasticOps.hpp, and plastic.cpp.

Definition at line 105 of file plastic.cpp.

◆ cn1

double cn1 = 1
Examples
PlasticOps.hpp, and plastic.cpp.

Definition at line 106 of file plastic.cpp.

◆ G

constexpr IntegrationType G
constexpr
Initial value:
=
IntegrationType::GAUSS
Examples
plastic.cpp.

Definition at line 38 of file plastic.cpp.

◆ geom_order

int geom_order = 2

Order if fixed.

Examples
plastic.cpp.

Definition at line 112 of file plastic.cpp.

◆ H

double H = 129
Examples
plastic.cpp.

Definition at line 103 of file plastic.cpp.

◆ help

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

[Solve]

Definition at line 1184 of file plastic.cpp.

◆ is_large_strains

PetscBool is_large_strains = PETSC_TRUE
Examples
contact.cpp, and plastic.cpp.

Definition at line 94 of file plastic.cpp.

◆ order

int order = 2

Order if fixed.

Definition at line 111 of file plastic.cpp.

◆ poisson_ratio

double poisson_ratio = 0.29

Definition at line 100 of file plastic.cpp.

◆ Qinf

double Qinf = 265
Examples
plastic.cpp.

Definition at line 108 of file plastic.cpp.

◆ rho

double rho = 0

◆ scale

double scale = 1.

◆ set_timer

PetscBool set_timer = PETSC_FALSE
Examples
plastic.cpp.

Definition at line 95 of file plastic.cpp.

◆ sigmaY

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

Definition at line 102 of file plastic.cpp.

◆ size_symm

constexpr auto size_symm = (SPACE_DIM * (SPACE_DIM + 1)) / 2
constexpr

◆ SPACE_DIM

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

Definition at line 31 of file plastic.cpp.

◆ visH

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

Definition at line 104 of file plastic.cpp.

◆ young_modulus

double young_modulus = 206913

◆ zeta

double zeta = 5e-2
Examples
PlasticOps.hpp, plastic.cpp, and plot_base.cpp.

Definition at line 107 of file plastic.cpp.