v0.15.0
Loading...
Searching...
No Matches
Classes | Namespaces | 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  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
 
#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 SetPtsData = FieldEvaluatorInterface::SetPtsData
 
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[])
 
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
 [Specialisation for assembly]
 
PetscBool is_large_strains = PETSC_TRUE
 Large strains.
 
PetscBool set_timer = PETSC_FALSE
 Set timer.
 
PetscBool do_eval_field = PETSC_FALSE
 Evaluate field.
 
int atom_test = 0
 Atom test.
 
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
 [Operators used for contact]
 
double rho = 0.0
 
double alpha_damping = 0
 
static char help [] = "...\n\n"
 [TestOperators]
 
static char help [] = "...\n\n"
 [TestOperators]
 

Macro Definition Documentation

◆ EXECUTABLE_DIMENSION [1/2]

#define EXECUTABLE_DIMENSION   3

◆ EXECUTABLE_DIMENSION [2/2]

#define EXECUTABLE_DIMENSION   3

Definition at line 13 of file plastic.cpp.

Typedef Documentation

◆ BoundaryEle

Definition at line 57 of file plastic.cpp.

◆ BoundaryEleOp

Definition at line 58 of file plastic.cpp.

◆ BoundaryLhsBCs

Definition at line 175 of file plastic.cpp.

◆ BoundaryRhsBCs

Definition at line 172 of file plastic.cpp.

◆ DomainEle

Definition at line 55 of file plastic.cpp.

◆ DomainEleOp

Definition at line 56 of file plastic.cpp.

◆ DomainRhsBCs

typedef NaturalBC< DomainEleOp >::Assembly< A >::LinearForm< I > DomainRhsBCs

Definition at line 169 of file plastic.cpp.

◆ EntData

Definition at line 54 of file plastic.cpp.

◆ OpBoundaryLhsBCs

Definition at line 176 of file plastic.cpp.

◆ OpBoundaryRhsBCs

Definition at line 173 of file plastic.cpp.

◆ OpDomainRhsBCs

Definition at line 170 of file plastic.cpp.

◆ PostProcEle

Definition at line 59 of file plastic.cpp.

◆ SetPtsData

Definition at line 62 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, mofem/tutorials/adv-0/plastic.cpp, mofem/tutorials/adv-0/src/PlasticOpsGeneric.hpp, plastic.cpp, and thermoplastic.cpp.

Definition at line 73 of file plastic.cpp.

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

◆ iso_hardening_dtau()

double iso_hardening_dtau ( double  tau,
double  H,
double  Qinf,
double  b_iso 
)
inline
Examples
PlasticOpsGeneric.hpp, mofem/tutorials/adv-0/plastic.cpp, mofem/tutorials/adv-0/src/PlasticOpsGeneric.hpp, plastic.cpp, and thermoplastic.cpp.

Definition at line 78 of file plastic.cpp.

79 {
80 auto r = [&](auto tau) {
81 return H + Qinf * b_iso * iso_hardening_exp(tau, b_iso);
82 };
83 constexpr double eps = 1e-12;
84 return std::max(r(tau), eps * r(0));
85}
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
mofem/tutorials/adv-0/plastic.cpp, plastic.cpp, and thermoplastic.cpp.

Definition at line 64 of file plastic.cpp.

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

◆ 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, mofem/tutorials/adv-0/plastic.cpp, mofem/tutorials/adv-0/src/PlasticOpsGeneric.hpp, plastic.cpp, and thermoplastic.cpp.

Definition at line 92 of file plastic.cpp.

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

◆ kinematic_hardening_dplastic_strain()

template<int DIM>
auto kinematic_hardening_dplastic_strain ( double  C1_k)
inline
Examples
mofem/tutorials/adv-0/plastic.cpp, plastic.cpp, and thermoplastic.cpp.

Definition at line 106 of file plastic.cpp.

106 {
107 FTensor::Index<'i', DIM> i;
108 FTensor::Index<'j', DIM> j;
109 FTensor::Index<'k', DIM> k;
110 FTensor::Index<'l', DIM> l;
113 t_diff(i, j, k, l) = C1_k * (t_kd(i, k) ^ t_kd(j, l)) / 4.;
114 return t_diff;
115}
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]

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

1610 {
1611
1612 //! [Register MoFEM discrete manager in PETSc]
1613 DMType dm_name = "DMMOFEM";
1614 CHKERR DMRegister_MoFEM(dm_name);
1615 //! [Register MoFEM discrete manager in PETSc
1616
1617 //! [Create MoAB]
1618 moab::Core mb_instance; ///< mesh database
1619 moab::Interface &moab = mb_instance; ///< mesh database interface
1620 //! [Create MoAB]
1621
1622 //! [Create MoFEM]
1623 MoFEM::Core core(moab); ///< finite element database
1624 MoFEM::Interface &m_field = core; ///< finite element database interface
1625 //! [Create MoFEM]
1626
1627 //! [Load mesh]
1628 Simple *simple = m_field.getInterface<Simple>();
1630 CHKERR simple->loadFile();
1631 //! [Load mesh]
1632
1633 //! [Example]
1634 Example ex(m_field);
1635 CHKERR ex.runProblem();
1636 //! [Example]
1637 }
1639
1641
1642#ifdef ADD_CONTACT
1643 #ifdef ENABLE_PYTHON_BINDING
1644 if (Py_FinalizeEx() < 0) {
1645 exit(120);
1646 }
1647 #endif
1648#endif // ADD_CONTACT
1649
1650 return 0;
1651}
1652
1653struct SetUpSchurImpl : public SetUpSchur {
1654
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
[Example]
Definition plastic.cpp:216
Core (interface) class.
Definition Core.hpp:82
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition Core.cpp:118
Deprecated interface functions.
Simple interface for fast problem set-up.
Definition Simple.hpp:27
MoFEMErrorCode getOptions()
get options
Definition Simple.cpp:180
intrusive_ptr for managing petsc objects
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
[Push operators to pipeline]

Variable Documentation

◆ alpha_damping

double alpha_damping = 0

◆ AT

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

Definition at line 44 of file plastic.cpp.

◆ atom_test

int atom_test = 0

◆ b_iso

double b_iso = 16.93

◆ C1_k

double C1_k = 0

◆ cn0

double cn0 = 1

◆ cn1

double cn1 = 1

◆ CONTACT_SPACE

constexpr FieldSpace CONTACT_SPACE = ElementsAndOps<SPACE_DIM>::CONTACT_SPACE
constexpr

◆ do_eval_field

PetscBool do_eval_field = PETSC_FALSE

◆ ep_order

int ep_order = order - 1

Order of ep files.

Examples
mofem/tutorials/adv-0/plastic.cpp, plastic.cpp, and thermoplastic.cpp.

Definition at line 140 of file plastic.cpp.

◆ geom_order

int geom_order = 2

Order if fixed.

Examples
SolutionMapping.hpp, mofem/tutorials/adv-0/plastic.cpp, plastic.cpp, and thermoplastic.cpp.

Definition at line 141 of file plastic.cpp.

◆ H

double H = 129

◆ help [1/2]

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

[TestOperators]

Definition at line 1578 of file plastic.cpp.

◆ help [2/2]

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

[TestOperators]

Definition at line 1574 of file plastic.cpp.

◆ is_large_strains

PetscBool is_large_strains = PETSC_TRUE

◆ is_quasi_static

PetscBool is_quasi_static = PETSC_TRUE

◆ IT

constexpr IntegrationType IT
constexpr
Initial value:
=
IntegrationType::GAUSS

Definition at line 47 of file plastic.cpp.

◆ order

int order = 2

Order displacement.

Definition at line 138 of file plastic.cpp.

◆ poisson_ratio

double poisson_ratio = 0.29

◆ Qinf

double Qinf = 265

◆ rho

double rho = 0.0

◆ scale

double scale = 1.

◆ set_timer

PetscBool set_timer = PETSC_FALSE

Set timer.

Examples
mofem/tutorials/adv-0/plastic.cpp, plastic.cpp, and thermoplastic.cpp.

Definition at line 118 of file plastic.cpp.

◆ sigmaY

double sigmaY = 450

◆ size_symm

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

◆ SPACE_DIM

constexpr int SPACE_DIM
constexpr
Initial value:
=
3

Definition at line 40 of file plastic.cpp.

◆ tau_order

int tau_order = order - 2

Order of tau files.

Examples
mofem/tutorials/adv-0/plastic.cpp, plastic.cpp, and thermoplastic.cpp.

Definition at line 139 of file plastic.cpp.

◆ visH

double visH = 0

◆ young_modulus

double young_modulus = 206913

◆ zeta

double zeta = 5e-2