![]() |
v0.15.0 |
#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 | SetPtsData = FieldEvaluatorInterface::SetPtsData |
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. | |
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 |
double | rho = 0.0 |
double | alpha_damping = 0 |
static char | help [] = "...\n\n" |
[TestOperators] | |
#define EXECUTABLE_DIMENSION 3 |
Definition at line 13 of file plastic.cpp.
Definition at line 57 of file plastic.cpp.
using BoundaryLhsBCs = NaturalBC<BoundaryEleOp>::Assembly<AT>::BiLinearForm<IT> |
Definition at line 175 of file plastic.cpp.
using BoundaryRhsBCs = NaturalBC<BoundaryEleOp>::Assembly<AT>::LinearForm<IT> |
Definition at line 172 of file plastic.cpp.
using DomainEle = ElementsAndOps<SPACE_DIM>::DomainEle |
Definition at line 55 of file plastic.cpp.
using DomainRhsBCs = NaturalBC<DomainEleOp>::Assembly<AT>::LinearForm<IT> |
Definition at line 169 of file plastic.cpp.
using OpBoundaryLhsBCs |
Definition at line 176 of file plastic.cpp.
using OpBoundaryRhsBCs |
Definition at line 173 of file plastic.cpp.
using OpDomainRhsBCs |
Definition at line 170 of file plastic.cpp.
using PostProcEle = PostProcBrokenMeshInMoab<DomainEle> |
Definition at line 59 of file plastic.cpp.
Definition at line 62 of file plastic.cpp.
using SideEle = ElementsAndOps<SPACE_DIM>::SideEle |
Definition at line 61 of file plastic.cpp.
Definition at line 60 of file plastic.cpp.
Isotropic hardening
Definition at line 73 of file plastic.cpp.
Definition at line 78 of file plastic.cpp.
Definition at line 64 of file plastic.cpp.
|
inline |
Kinematic hardening
Definition at line 92 of file plastic.cpp.
|
inline |
Definition at line 106 of file plastic.cpp.
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 1580 of file plastic.cpp.
double alpha_damping = 0 |
Definition at line 145 of file plastic.cpp.
|
constexpr |
Definition at line 44 of file plastic.cpp.
int atom_test = 0 |
double b_iso = 16.93 |
Saturation exponent.
Definition at line 132 of file plastic.cpp.
double C1_k = 0 |
Kinematic hardening.
Definition at line 133 of file plastic.cpp.
double cn0 = 1 |
Definition at line 135 of file plastic.cpp.
double cn1 = 1 |
Definition at line 136 of file plastic.cpp.
|
constexpr |
Definition at line 52 of file plastic.cpp.
PetscBool do_eval_field = PETSC_FALSE |
Evaluate field.
Definition at line 119 of file plastic.cpp.
int ep_order = order - 1 |
int geom_order = 2 |
double H = 129 |
Hardening.
Definition at line 128 of file plastic.cpp.
|
static |
[TestOperators]
Definition at line 1578 of file plastic.cpp.
PetscBool is_large_strains = PETSC_TRUE |
PetscBool is_quasi_static = PETSC_TRUE |
Definition at line 143 of file plastic.cpp.
|
constexpr |
Definition at line 47 of file plastic.cpp.
int order = 2 |
Order displacement.
Definition at line 138 of file plastic.cpp.
double poisson_ratio = 0.29 |
Poisson ratio.
Definition at line 126 of file plastic.cpp.
double Qinf = 265 |
Saturation yield stress.
Definition at line 131 of file plastic.cpp.
double rho = 0.0 |
Definition at line 144 of file plastic.cpp.
double scale = 1. |
Definition at line 123 of file plastic.cpp.
PetscBool set_timer = PETSC_FALSE |
double sigmaY = 450 |
Yield stress.
Definition at line 127 of file plastic.cpp.
Definition at line 42 of file plastic.cpp.
|
constexpr |
Definition at line 40 of file plastic.cpp.
int tau_order = order - 2 |
double visH = 0 |
Viscous hardening.
Definition at line 129 of file plastic.cpp.
double young_modulus = 206913 |
Young modulus.
Definition at line 125 of file plastic.cpp.
double zeta = 5e-2 |
Viscous hardening.
Definition at line 130 of file plastic.cpp.