v0.14.0
Loading...
Searching...
No Matches
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  PlasticOps::AddHOOps< 2, 3, 3 >
 
struct  PlasticOps::AddHOOps< 1, 2, 2 >
 
struct  PlasticOps::AddHOOps< 3, 3, 3 >
 
struct  PlasticOps::AddHOOps< 2, 2, 2 >
 
struct  Example
 [Example] More...
 
struct  Example::ScaledTimeScale
 
struct  SetUpSchur
 [Push operators to pipeline] More...
 
struct  SetUpSchurImpl
 

Namespaces

namespace  PlasticOps
 

Macros

#define EXECUTABLE_DIMENSION   3
 

Typedefs

using DomainEle = ElementsAndOps<SPACE_DIM>::DomainEle
 
using BoundaryEle = ElementsAndOps<SPACE_DIM>::BoundaryEle
 
using PostProcEle = PostProcBrokenMeshInMoab<DomainEle>
 
using SkinPostProcEle = PostProcBrokenMeshInMoab<BoundaryEle>
 
using SideEle = ElementsAndOps<SPACE_DIM>::SideEle
 
using DomainRhsBCs = NaturalBC<DomainEleOp>::Assembly<AT>::LinearForm<IT>
 
using OpDomainRhsBCs
 
using BoundaryRhsBCs = NaturalBC<BoundaryEleOp>::Assembly<AT>::LinearForm<IT>
 
using OpBoundaryRhsBCs
 
using BoundaryLhsBCs = NaturalBC<BoundaryEleOp>::Assembly<AT>::BiLinearForm<IT>
 
using OpBoundaryLhsBCs
 

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[])
 
template<int FE_DIM, int PROBLEM_DIM, int SPACE_DIM>
MoFEMErrorCode PlasticOps::scaleL2 (boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pipeline, std::string geom_field_name)
 

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.
 
PetscBool set_timer = PETSC_FALSE
 Set timer.
 
double scale = 1.
 
double young_modulus = 206913
 Young modulus.
 
double poisson_ratio = 0.29
 Poisson ratio.
 
double sigmaY = 450
 Yield stress.
 
double H = 129
 Hardening.
 
double visH = 0
 Viscous hardening.
 
double zeta = 5e-2
 Viscous hardening.
 
double Qinf = 265
 Saturation yield stress.
 
double b_iso = 16.93
 Saturation exponent.
 
double C1_k = 0
 Kinematic hardening.
 
double cn0 = 1
 
double cn1 = 1
 
int order = 2
 Order displacement.
 
int tau_order = order - 2
 Order of tau files.
 
int ep_order = order - 1
 Order of ep files.
 
int geom_order = 2
 Order if fixed.
 
PetscBool is_quasi_static = PETSC_TRUE
 
double rho = 0.0
 
double alpha_damping = 0
 
static char help [] = "...\n\n"
 [TestOperators]
 

Macro Definition Documentation

◆ EXECUTABLE_DIMENSION

Typedef Documentation

◆ BoundaryEle

Definition at line 57 of file plastic.cpp.

◆ BoundaryLhsBCs

Definition at line 171 of file plastic.cpp.

◆ BoundaryRhsBCs

Definition at line 168 of file plastic.cpp.

◆ DomainEle

Definition at line 55 of file plastic.cpp.

◆ DomainRhsBCs

Definition at line 165 of file plastic.cpp.

◆ OpBoundaryLhsBCs

◆ OpBoundaryRhsBCs

◆ OpDomainRhsBCs

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

73 {
74 return H * tau + Qinf * (1. - iso_hardening_exp(tau, b_iso)) + sigmaY;
75}
double Qinf
Saturation yield stress.
Definition plastic.cpp:127
double H
Hardening.
Definition plastic.cpp:124
double iso_hardening_exp(double tau, double b_iso)
Definition plastic.cpp:63
double b_iso
Saturation exponent.
Definition plastic.cpp:128
double sigmaY
Yield stress.
Definition plastic.cpp:123

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

78 {
79 auto r = [&](auto tau) {
80 return H + Qinf * b_iso * iso_hardening_exp(tau, b_iso);
81 };
82 constexpr double eps = 1e-12;
83 return std::max(r(tau), eps * r(0));
84}
constexpr double eps
Definition HenckyOps.hpp:13
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 63 of file plastic.cpp.

63 {
64 return std::exp(
65 std::max(static_cast<double>(std::numeric_limits<float>::min_exponent10),
66 -b_iso * tau));
67}

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

92 {
93 FTensor::Index<'i', DIM> i;
94 FTensor::Index<'j', DIM> j;
96 if (C1_k < std::numeric_limits<double>::epsilon()) {
97 t_alpha(i, j) = 0;
98 return t_alpha;
99 }
100 t_alpha(i, j) = C1_k * t_plastic_strain(i, j);
101 return t_alpha;
102}
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'j', 3 > j
double C1_k
Kinematic hardening.
Definition plastic.cpp:129

◆ kinematic_hardening_dplastic_strain()

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

Definition at line 105 of file plastic.cpp.

105 {
106 FTensor::Index<'i', DIM> i;
107 FTensor::Index<'j', DIM> j;
108 FTensor::Index<'k', DIM> k;
109 FTensor::Index<'l', DIM> l;
112 t_diff(i, j, k, l) = C1_k * (t_kd(i, k) ^ t_kd(j, l)) / 4.;
113 return t_diff;
114}
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 1390 of file plastic.cpp.

1390 {
1391
1392#ifdef ADD_CONTACT
1393 #ifdef PYTHON_SDF
1394 Py_Initialize();
1395 np::initialize();
1396 #endif
1397#endif // ADD_CONTACT
1398
1399 // Initialisation of MoFEM/PETSc and MOAB data structures
1400 const char param_file[] = "param_file.petsc";
1401 MoFEM::Core::Initialize(&argc, &argv, param_file, help);
1402
1403 // Add logging channel for example
1404 auto core_log = logging::core::get();
1405 core_log->add_sink(
1407 core_log->add_sink(
1409 LogManager::setLog("PLASTICITY");
1410 MOFEM_LOG_TAG("PLASTICITY", "Plasticity");
1411
1412#ifdef ADD_CONTACT
1413 core_log->add_sink(
1415 LogManager::setLog("CONTACT");
1416 MOFEM_LOG_TAG("CONTACT", "Contact");
1417#endif // ADD_CONTACT
1418
1419 try {
1420
1421 //! [Register MoFEM discrete manager in PETSc]
1422 DMType dm_name = "DMMOFEM";
1423 CHKERR DMRegister_MoFEM(dm_name);
1424 //! [Register MoFEM discrete manager in PETSc
1425
1426 //! [Create MoAB]
1427 moab::Core mb_instance; ///< mesh database
1428 moab::Interface &moab = mb_instance; ///< mesh database interface
1429 //! [Create MoAB]
1430
1431 //! [Create MoFEM]
1432 MoFEM::Core core(moab); ///< finite element database
1433 MoFEM::Interface &m_field = core; ///< finite element database interface
1434 //! [Create MoFEM]
1435
1436 //! [Load mesh]
1437 Simple *simple = m_field.getInterface<Simple>();
1438 CHKERR simple->getOptions();
1439 CHKERR simple->loadFile();
1440 //! [Load mesh]
1441
1442 //! [Example]
1443 Example ex(m_field);
1444 CHKERR ex.runProblem();
1445 //! [Example]
1446 }
1448
1450
1451#ifdef ADD_CONTACT
1452 #ifdef PYTHON_SDF
1453 if (Py_FinalizeEx() < 0) {
1454 exit(120);
1455 }
1456 #endif
1457#endif // ADD_CONTACT
1458
1459 return 0;
1460}
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition acoustic.cpp:69
#define CATCH_ERRORS
Catch errors.
#define CHKERR
Inline error 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.
static char help[]
[TestOperators]
Definition plastic.cpp:1388
[Example]
Definition plastic.cpp:213
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.
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
Simple interface for fast problem set-up.
Definition Simple.hpp:27
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.

Variable Documentation

◆ alpha_damping

double alpha_damping = 0
Examples
plastic.cpp.

Definition at line 141 of file plastic.cpp.

◆ AT

AssemblyType AT
constexpr
Initial value:
=
(SCHUR_ASSEMBLE) ? AssemblyType::BLOCK_SCHUR
#define SCHUR_ASSEMBLE
Definition contact.cpp:18
AssemblyType
[Storage and set boundary conditions]

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

◆ C1_k

double C1_k = 0

Kinematic hardening.

Examples
PlasticOps.hpp, and plastic.cpp.

Definition at line 129 of file plastic.cpp.

◆ cn0

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

Definition at line 131 of file plastic.cpp.

◆ cn1

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

Definition at line 132 of file plastic.cpp.

◆ CONTACT_SPACE

FieldSpace CONTACT_SPACE = ElementsAndOps<SPACE_DIM>::CONTACT_SPACE
constexpr
Examples
adolc_plasticity.cpp, and 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 136 of file plastic.cpp.

◆ geom_order

int geom_order = 2

Order if fixed.

Examples
plastic.cpp.

Definition at line 137 of file plastic.cpp.

◆ H

◆ help

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

[TestOperators]

Definition at line 1388 of file plastic.cpp.

◆ is_large_strains

PetscBool is_large_strains = PETSC_TRUE

Large strains.

Examples
PlasticOps.hpp, and plastic.cpp.

Definition at line 116 of file plastic.cpp.

◆ is_quasi_static

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

Definition at line 139 of file plastic.cpp.

◆ IT

IntegrationType IT
constexpr
Initial value:
=
IntegrationType::GAUSS

Definition at line 47 of file plastic.cpp.

◆ order

int order = 2

Order displacement.

Definition at line 134 of file plastic.cpp.

◆ poisson_ratio

◆ Qinf

double Qinf = 265

Saturation yield stress.

Examples
PlasticOps.hpp, and plastic.cpp.

Definition at line 127 of file plastic.cpp.

◆ rho

◆ scale

◆ set_timer

PetscBool set_timer = PETSC_FALSE

Set timer.

Examples
plastic.cpp.

Definition at line 117 of file plastic.cpp.

◆ sigmaY

double sigmaY = 450

Yield stress.

Examples
ADOLCPlasticityMaterialModels.hpp, PlasticOps.hpp, and plastic.cpp.

Definition at line 123 of file plastic.cpp.

◆ size_symm

◆ SPACE_DIM

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

◆ visH

double visH = 0

Viscous hardening.

Examples
PlasticOps.hpp, and plastic.cpp.

Definition at line 125 of file plastic.cpp.

◆ young_modulus

◆ zeta

double zeta = 5e-2

Viscous hardening.

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

Definition at line 126 of file plastic.cpp.