v0.14.0
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>
#include <PlasticNaturalBCs.hpp>

Go to the source code of this file.

Classes

struct  ElementsAndOps< 2 >
 
struct  ElementsAndOps< 3 >
 
struct  Example
 [Example] More...
 
struct  Example::ScaledTimeScale
 
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 SkinPostProcEle = PostProcBrokenMeshInMoab< BoundaryEle >
 
using SideEle = ElementsAndOps< SPACE_DIM >::SideEle
 
using DomainRhsBCs = NaturalBC< DomainEleOp >::Assembly< AT >::LinearForm< IT >
 
using OpDomainRhsBCs = DomainRhsBCs::OpFlux< PlasticOps::DomainBCs, 1, SPACE_DIM >
 
using BoundaryRhsBCs = NaturalBC< BoundaryEleOp >::Assembly< AT >::LinearForm< IT >
 
using OpBoundaryRhsBCs = BoundaryRhsBCs::OpFlux< PlasticOps::BoundaryBCs, 1, SPACE_DIM >
 
using BoundaryLhsBCs = NaturalBC< BoundaryEleOp >::Assembly< AT >::BiLinearForm< IT >
 
using OpBoundaryLhsBCs = BoundaryLhsBCs::OpFlux< PlasticOps::BoundaryBCs, 1, SPACE_DIM >
 

Functions

double iso_hardening_exp (double tau, double b_iso)
 
double iso_hardening (double tau, double H, double Qinf, double b_iso, double sigmaY)
 
double iso_hardening_dtau (double tau, double H, double Qinf, double b_iso)
 
template<typename T , int DIM>
auto kinematic_hardening (FTensor::Tensor2_symmetric< T, DIM > &t_plastic_strain, double C1_k)
 
template<int DIM>
auto kinematic_hardening_dplastic_strain (double C1_k)
 
int main (int argc, char *argv[])
 

Variables

constexpr int SPACE_DIM
 
constexpr auto size_symm = (SPACE_DIM * (SPACE_DIM + 1)) / 2
 
constexpr AssemblyType AT
 
constexpr IntegrationType IT
 
constexpr FieldSpace CONTACT_SPACE = ElementsAndOps<SPACE_DIM>::CONTACT_SPACE
 
PetscBool is_large_strains = PETSC_TRUE
 Large strains. More...
 
PetscBool set_timer = PETSC_FALSE
 Set timer. More...
 
double scale = 1.
 
double young_modulus = 206913
 Young modulus. More...
 
double poisson_ratio = 0.29
 Poisson ratio. More...
 
double sigmaY = 450
 Yield stress. More...
 
double H = 129
 Hardening. More...
 
double visH = 0
 Viscous hardening. More...
 
double zeta = 5e-2
 Viscous hardening. More...
 
double Qinf = 265
 Saturation yield stress. More...
 
double b_iso = 16.93
 Saturation exponent. More...
 
double C1_k = 0
 Kinematic hardening. More...
 
double cn0 = 1
 
double cn1 = 1
 
int order = 2
 Order displacement. More...
 
int tau_order = order - 2
 Order of tau files. More...
 
int ep_order = order - 1
 Order of ep files. More...
 
int geom_order = 2
 Order if fixed. More...
 
PetscBool is_quasi_static = PETSC_TRUE
 
double rho = 0.0
 
double alpha_damping = 0
 
static char help [] = "...\n\n"
 [Solve] More...
 

Macro Definition Documentation

◆ EXECUTABLE_DIMENSION

#define EXECUTABLE_DIMENSION   3

Typedef Documentation

◆ BoundaryEle

Definition at line 57 of file plastic.cpp.

◆ BoundaryEleOp

Definition at line 58 of file plastic.cpp.

◆ BoundaryLhsBCs

Definition at line 220 of file plastic.cpp.

◆ BoundaryRhsBCs

Definition at line 217 of file plastic.cpp.

◆ DomainEle

Definition at line 55 of file plastic.cpp.

◆ DomainEleOp

Definition at line 56 of file plastic.cpp.

◆ DomainRhsBCs

Definition at line 214 of file plastic.cpp.

◆ EntData

Definition at line 54 of file plastic.cpp.

◆ OpBoundaryLhsBCs

using OpBoundaryLhsBCs = BoundaryLhsBCs::OpFlux<PlasticOps::BoundaryBCs, 1, SPACE_DIM>

Definition at line 221 of file plastic.cpp.

◆ OpBoundaryRhsBCs

using OpBoundaryRhsBCs = BoundaryRhsBCs::OpFlux<PlasticOps::BoundaryBCs, 1, SPACE_DIM>

Definition at line 218 of file plastic.cpp.

◆ OpDomainRhsBCs

using OpDomainRhsBCs = DomainRhsBCs::OpFlux<PlasticOps::DomainBCs, 1, SPACE_DIM>

Definition at line 215 of file plastic.cpp.

◆ PostProcEle

Definition at line 59 of file plastic.cpp.

◆ SideEle

Definition at line 61 of file plastic.cpp.

◆ SkinPostProcEle

Definition at line 60 of file plastic.cpp.

Function Documentation

◆ iso_hardening()

double iso_hardening ( double  tau,
double  H,
double  Qinf,
double  b_iso,
double  sigmaY 
)
inline

Isotropic hardening

Examples
PlasticOpsGeneric.hpp, and plastic.cpp.

Definition at line 123 of file plastic.cpp.

124 {
125 return H * tau + Qinf * (1. - iso_hardening_exp(tau, b_iso)) + sigmaY;
126}
double Qinf
Saturation yield stress.
Definition: plastic.cpp:178
double H
Hardening.
Definition: plastic.cpp:175
double iso_hardening_exp(double tau, double b_iso)
Definition: plastic.cpp:114
double b_iso
Saturation exponent.
Definition: plastic.cpp:179
double sigmaY
Yield stress.
Definition: plastic.cpp:174

◆ iso_hardening_dtau()

double iso_hardening_dtau ( double  tau,
double  H,
double  Qinf,
double  b_iso 
)
inline
Examples
PlasticOpsGeneric.hpp, and plastic.cpp.

Definition at line 128 of file plastic.cpp.

129 {
130 auto r = [&](auto tau) {
131 return H + Qinf * b_iso * iso_hardening_exp(tau, b_iso);
132 };
133 constexpr double eps = 1e-12;
134 return std::max(r(tau), eps * r(0));
135}
static const double eps
int r
Definition: sdf.py:8

◆ iso_hardening_exp()

double iso_hardening_exp ( double  tau,
double  b_iso 
)
inline
Examples
plastic.cpp.

Definition at line 114 of file plastic.cpp.

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

◆ kinematic_hardening()

template<typename T , int DIM>
auto kinematic_hardening ( FTensor::Tensor2_symmetric< T, DIM > &  t_plastic_strain,
double  C1_k 
)
inline

Kinematic hardening

Examples
PlasticOpsGeneric.hpp, and plastic.cpp.

Definition at line 142 of file plastic.cpp.

143 {
147 if (C1_k < std::numeric_limits<double>::epsilon()) {
148 t_alpha(i, j) = 0;
149 return t_alpha;
150 }
151 t_alpha(i, j) = C1_k * t_plastic_strain(i, j);
152 return t_alpha;
153}
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'j', 3 > j
double C1_k
Kinematic hardening.
Definition: plastic.cpp:180

◆ kinematic_hardening_dplastic_strain()

template<int DIM>
auto kinematic_hardening_dplastic_strain ( double  C1_k)
inline
Examples
plastic.cpp.

Definition at line 156 of file plastic.cpp.

156 {
163 t_diff(i, j, k, l) = C1_k * (t_kd(i, k) ^ t_kd(j, l)) / 4.;
164 return t_diff;
165}
Kronecker Delta class symmetric.
constexpr auto t_kd
FTensor::Index< 'l', 3 > l
FTensor::Index< 'k', 3 > k

◆ 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 1192 of file plastic.cpp.

1192 {
1193
1194#ifdef ADD_CONTACT
1195#ifdef PYTHON_SFD
1196 Py_Initialize();
1197#endif
1198#endif // ADD_CONTACT
1199
1200 // Initialisation of MoFEM/PETSc and MOAB data structures
1201 const char param_file[] = "param_file.petsc";
1202 MoFEM::Core::Initialize(&argc, &argv, param_file, help);
1203
1204 // Add logging channel for example
1205 auto core_log = logging::core::get();
1206 core_log->add_sink(
1208 core_log->add_sink(
1210 LogManager::setLog("PLASTICITY");
1211 MOFEM_LOG_TAG("PLASTICITY", "Plasticity");
1212
1213#ifdef ADD_CONTACT
1214 core_log->add_sink(
1216 LogManager::setLog("CONTACT");
1217 MOFEM_LOG_TAG("CONTACT", "Contact");
1218#endif // ADD_CONTACT
1219
1220 try {
1221
1222 //! [Register MoFEM discrete manager in PETSc]
1223 DMType dm_name = "DMMOFEM";
1224 CHKERR DMRegister_MoFEM(dm_name);
1225 //! [Register MoFEM discrete manager in PETSc
1226
1227 //! [Create MoAB]
1228 moab::Core mb_instance; ///< mesh database
1229 moab::Interface &moab = mb_instance; ///< mesh database interface
1230 //! [Create MoAB]
1231
1232 //! [Create MoFEM]
1233 MoFEM::Core core(moab); ///< finite element database
1234 MoFEM::Interface &m_field = core; ///< finite element database interface
1235 //! [Create MoFEM]
1236
1237 //! [Load mesh]
1238 Simple *simple = m_field.getInterface<Simple>();
1239 CHKERR simple->getOptions();
1240 CHKERR simple->loadFile();
1241 //! [Load mesh]
1242
1243 //! [Example]
1244 Example ex(m_field);
1245 CHKERR ex.runProblem();
1246 //! [Example]
1247 }
1249
1251
1252#ifdef ADD_CONTACT
1253#ifdef PYTHON_SFD
1254 if (Py_FinalizeEx() < 0) {
1255 exit(120);
1256 }
1257#endif
1258#endif // ADD_CONTACT
1259
1260 return 0;
1261}
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:389
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
static char help[]
[Solve]
Definition: plastic.cpp:1190
[Example]
Definition: plastic.cpp:226
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
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

◆ alpha_damping

double alpha_damping = 0
Examples
plastic.cpp.

Definition at line 192 of file plastic.cpp.

◆ AT

constexpr AssemblyType AT
constexpr
Initial value:
=
(SCHUR_ASSEMBLE) ? AssemblyType::SCHUR
: AssemblyType::PETSC
Examples
plastic.cpp.

Definition at line 44 of file plastic.cpp.

◆ b_iso

double b_iso = 16.93

Saturation exponent.

Examples
PlasticOps.hpp, and plastic.cpp.

Definition at line 179 of file plastic.cpp.

◆ C1_k

double C1_k = 0

Kinematic hardening.

Examples
PlasticOps.hpp, and plastic.cpp.

Definition at line 180 of file plastic.cpp.

◆ cn0

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

Definition at line 182 of file plastic.cpp.

◆ cn1

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

Definition at line 183 of file plastic.cpp.

◆ CONTACT_SPACE

constexpr FieldSpace CONTACT_SPACE = ElementsAndOps<SPACE_DIM>::CONTACT_SPACE
constexpr
Examples
plastic.cpp.

Definition at line 52 of file plastic.cpp.

◆ ep_order

int ep_order = order - 1

Order of ep files.

Examples
plastic.cpp.

Definition at line 187 of file plastic.cpp.

◆ geom_order

int geom_order = 2

Order if fixed.

Examples
plastic.cpp.

Definition at line 188 of file plastic.cpp.

◆ H

double H = 129

Hardening.

Examples
PlasticOps.hpp, and plastic.cpp.

Definition at line 175 of file plastic.cpp.

◆ help

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

[Solve]

Definition at line 1190 of file plastic.cpp.

◆ is_large_strains

PetscBool is_large_strains = PETSC_TRUE

Large strains.

Examples
PlasticOps.hpp, and plastic.cpp.

Definition at line 167 of file plastic.cpp.

◆ is_quasi_static

PetscBool is_quasi_static = PETSC_TRUE
Examples
NonlinearElasticElementInterface.hpp, and plastic.cpp.

Definition at line 190 of file plastic.cpp.

◆ IT

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

Definition at line 47 of file plastic.cpp.

◆ order

int order = 2

Order displacement.

Definition at line 185 of file plastic.cpp.

◆ poisson_ratio

double poisson_ratio = 0.29

Poisson ratio.

Definition at line 173 of file plastic.cpp.

◆ Qinf

double Qinf = 265

Saturation yield stress.

Examples
PlasticOps.hpp, and plastic.cpp.

Definition at line 178 of file plastic.cpp.

◆ rho

double rho = 0.0

◆ scale

double scale = 1.

◆ set_timer

PetscBool set_timer = PETSC_FALSE

Set timer.

Examples
plastic.cpp.

Definition at line 168 of file plastic.cpp.

◆ sigmaY

double sigmaY = 450

Yield stress.

Examples
PlasticOps.hpp, and plastic.cpp.

Definition at line 174 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:13

Definition at line 40 of file plastic.cpp.

◆ tau_order

int tau_order = order - 2

Order of tau files.

Examples
plastic.cpp.

Definition at line 186 of file plastic.cpp.

◆ visH

double visH = 0

Viscous hardening.

Examples
PlasticOps.hpp, and plastic.cpp.

Definition at line 176 of file plastic.cpp.

◆ young_modulus

double young_modulus = 206913

◆ zeta

double zeta = 5e-2

Viscous hardening.

Examples
PlasticOpsGeneric.hpp, plastic.cpp, and plot_base.cpp.

Definition at line 177 of file plastic.cpp.