v0.14.0
Public Member Functions | Public Attributes | List of all members
SmallStrainJ2Plasticity Struct Reference

J2 plasticity (Kinematic Isotropic (Linear) Hardening) More...

#include <users_modules/small_strain_plasticity/src/SmallStrainPlasticityMaterialModels.hpp>

Inheritance diagram for SmallStrainJ2Plasticity:
[legend]
Collaboration diagram for SmallStrainJ2Plasticity:
[legend]

Public Member Functions

 SmallStrainJ2Plasticity ()
 
virtual PetscErrorCode freeHemholtzEnergy ()
 Available energy. More...
 
PetscErrorCode evalF ()
 
virtual PetscErrorCode yieldFunction ()
 Evaluate yield function. More...
 
virtual PetscErrorCode flowPotential ()
 Associated flow potential. More...
 
- Public Member Functions inherited from SmallStrainPlasticity::ClosestPointProjection
 ClosestPointProjection ()
 
virtual PetscErrorCode setActiveVariablesW ()
 
virtual PetscErrorCode setActiveVariablesYH ()
 
virtual PetscErrorCode rEcordW ()
 
virtual PetscErrorCode rEcordY ()
 
virtual PetscErrorCode rEcordH ()
 
virtual PetscErrorCode pLayW ()
 
virtual PetscErrorCode pLayW_NoHessian ()
 
virtual PetscErrorCode pLayY ()
 
virtual PetscErrorCode pLayY_NoGradient ()
 
virtual PetscErrorCode pLayH ()
 
virtual PetscErrorCode createMatAVecR ()
 
virtual PetscErrorCode destroyMatAVecR ()
 
virtual PetscErrorCode evaluatePotentials ()
 
virtual PetscErrorCode cAlculateR (Vec R)
 
virtual PetscErrorCode cAlculateA ()
 
PetscErrorCode snesCreate ()
 
PetscErrorCode snesDestroy ()
 
PetscErrorCode solveColasetProjection ()
 
PetscErrorCode consistentTangent ()
 

Public Attributes

double mu
 
double lambda
 
double H
 
double K
 
double phi
 
double sIgma_y
 
adouble t
 Auxiliary function. More...
 
adouble f
 
- Public Attributes inherited from SmallStrainPlasticity::ClosestPointProjection
VectorDouble sTrain
 
VectorDouble internalVariables
 
VectorDouble internalVariables0
 
VectorDouble plasticStrain0
 
VectorDouble plasticStrain
 
double deltaGamma
 
double tOl
 
int gG
 
VectorAdaptor sTress
 
VectorAdaptor internalFluxes
 
ublas::symmetric_matrix< double, ublas::lower > C
 
ublas::symmetric_matrix< double, ublas::lower > D
 
MatrixDouble partialWStrainPlasticStrain
 
VectorAdaptor partialYSigma
 
VectorAdaptor partialYFlux
 
VectorAdaptor partialHSigma
 
VectorAdaptor partialHFlux
 
ublas::symmetric_matrix< double, ublas::lower > partial2HSigma
 
ublas::symmetric_matrix< double, ublas::lower > partial2HFlux
 
MatrixDouble partial2HSigmaFlux
 
vector< int > tAgs
 
VectorDouble activeVariablesW
 
VectorDouble activeVariablesYH
 
double w
 
double y
 
double h
 
ublas::vector< adoublea_sTrain
 
ublas::vector< adoublea_plasticStrain
 
ublas::vector< adoublea_internalVariables
 
ublas::vector< adoublea_sTress
 
ublas::vector< adoublea_internalFluxes
 
adouble a_w
 
adouble a_y
 
adouble a_h
 
VectorDouble gradientW
 
MatrixDouble hessianW
 
VectorDouble gradientY
 
VectorDouble gradientH
 
MatrixDouble hessianH
 
Mat A
 
ublas::matrix< double, ublas::column_major > dataA
 
Vec R
 
Vec Chi
 
Vec dChi
 
SNES sNes
 
KSP kSp
 
PC pC
 
SNESLineSearch lineSearch
 
MatrixDouble Ep
 
MatrixDouble Lp
 
MatrixDouble Cp
 
MatrixDouble Cep
 
VectorDouble Yp
 

Detailed Description

J2 plasticity (Kinematic Isotropic (Linear) Hardening)

Examples
SmallStrainPlasticityMaterialModels.hpp.

Definition at line 32 of file SmallStrainPlasticityMaterialModels.hpp.

Constructor & Destructor Documentation

◆ SmallStrainJ2Plasticity()

SmallStrainJ2Plasticity::SmallStrainJ2Plasticity ( )
inline

Member Function Documentation

◆ evalF()

PetscErrorCode SmallStrainJ2Plasticity::evalF ( )
inline
Examples
SmallStrainPlasticityMaterialModels.hpp.

Definition at line 107 of file SmallStrainPlasticityMaterialModels.hpp.

107  {
108  PetscFunctionBegin;
109  adouble tR = 0;
110  for(int dd = 0;dd<3;dd++) {
111  tR += a_sTress[dd];
112  }
113  // Fix the constant to mach uniaxial test
114  const double c = 3./2.;
115  tR /= 3.;
116  f = 0;
117  for(int dd = 0;dd<3;dd++) {
118  t = (a_sTress[dd]-tR)-a_internalFluxes[dd+1];
119  f += c*t*t;
120  }
121  for(int dd = 3;dd<6;dd++) {
123  f += c*2.*t*t; // s:s off diagonal terms need to go twice
124  }
125  PetscFunctionReturn(0);
126  }

◆ flowPotential()

virtual PetscErrorCode SmallStrainJ2Plasticity::flowPotential ( )
inlinevirtual

Associated flow potential.

See yieldFunction.

Reimplemented from SmallStrainPlasticity::ClosestPointProjection.

Examples
SmallStrainPlasticityMaterialModels.hpp.

Definition at line 148 of file SmallStrainPlasticityMaterialModels.hpp.

148  {
149  PetscFunctionBegin;
150  ierr = evalF(); CHKERRQ(ierr);
151  a_h = sqrt(f) - a_internalFluxes[0];
152  PetscFunctionReturn(0);
153  }

◆ freeHemholtzEnergy()

virtual PetscErrorCode SmallStrainJ2Plasticity::freeHemholtzEnergy ( )
inlinevirtual

Available energy.

\[ \psi = \frac{1}{2} \lambda \textrm{tr}[\varepsilon]^2 + \mu \varepsilon : \varepsilon + \sigma_y\alpha + \frac{1}{2} \phi H \alpha^2 + \frac{1}{2} (1-\phi)K \beta^2 \]

Isotropic hardening variable \(\alpha\) is in first index of internal variables vector. Kinematic hardening variable is in the remaining indices of internal variables vector.

Reimplemented from SmallStrainPlasticity::ClosestPointProjection.

Examples
SmallStrainPlasticityMaterialModels.hpp.

Definition at line 65 of file SmallStrainPlasticityMaterialModels.hpp.

65  {
66  PetscFunctionBegin;
67  //elastic part
68  adouble tR = 0;
69  for(int dd = 0;dd<3;dd++) {
70  tR += a_sTrain[dd]-a_plasticStrain[dd];
71  }
72  a_w = 0.5*lambda*tR*tR;
73  adouble t;
74  int dd = 0;
75  for(;dd<3;dd++) {
77  a_w += mu*t*t;
78  }
79  for(;dd<6;dd++) {
81  a_w += 0.5*mu*t*t;
82  }
83  //plastic part
84  //isotropic
86  //kinematic
87  for(unsigned int dd = 1;dd<a_internalVariables.size();dd++) {
89  }
90  PetscFunctionReturn(0);
91  }

◆ yieldFunction()

virtual PetscErrorCode SmallStrainJ2Plasticity::yieldFunction ( )
inlinevirtual

Evaluate yield function.

\[ y = \sqrt{f} - \overline{\alpha} \]

where \(f\) is defined in evalF.

Reimplemented from SmallStrainPlasticity::ClosestPointProjection.

Examples
SmallStrainPlasticityMaterialModels.hpp.

Definition at line 136 of file SmallStrainPlasticityMaterialModels.hpp.

136  {
137  PetscFunctionBegin;
138  ierr = evalF(); CHKERRQ(ierr);
139  a_y = sqrt(f) - a_internalFluxes[0];
140  PetscFunctionReturn(0);
141  }

Member Data Documentation

◆ f

adouble SmallStrainJ2Plasticity::f

◆ H

double SmallStrainJ2Plasticity::H

◆ K

double SmallStrainJ2Plasticity::K

◆ lambda

double SmallStrainJ2Plasticity::lambda

◆ mu

double SmallStrainJ2Plasticity::mu

◆ phi

double SmallStrainJ2Plasticity::phi

◆ sIgma_y

double SmallStrainJ2Plasticity::sIgma_y

◆ t

adouble SmallStrainJ2Plasticity::t

Auxiliary function.

\[ \eta=\textrm{dev}[\sigma]-\overline{\beta} \]

\[ f = \frac{3}{2} \eta:\eta \]

This is \(3J_2\).

Examples
SmallStrainPlasticityMaterialModels.hpp.

Definition at line 106 of file SmallStrainPlasticityMaterialModels.hpp.


The documentation for this struct was generated from the following file:
SmallStrainPlasticity::ClosestPointProjection::a_plasticStrain
ublas::vector< adouble > a_plasticStrain
Definition: SmallStrainPlasticity.hpp:517
SmallStrainJ2Plasticity::sIgma_y
double sIgma_y
Definition: SmallStrainPlasticityMaterialModels.hpp:45
SmallStrainPlasticity::ClosestPointProjection
Closest point projection algorithm.
Definition: SmallStrainPlasticity.hpp:482
SmallStrainJ2Plasticity::f
adouble f
Definition: SmallStrainPlasticityMaterialModels.hpp:106
SmallStrainJ2Plasticity::evalF
PetscErrorCode evalF()
Definition: SmallStrainPlasticityMaterialModels.hpp:107
SmallStrainPlasticity::ClosestPointProjection::a_internalFluxes
ublas::vector< adouble > a_internalFluxes
Definition: SmallStrainPlasticity.hpp:518
c
const double c
speed of light (cm/ns)
Definition: initial_diffusion.cpp:39
SmallStrainJ2Plasticity::phi
double phi
Definition: SmallStrainPlasticityMaterialModels.hpp:44
SmallStrainJ2Plasticity::t
adouble t
Auxiliary function.
Definition: SmallStrainPlasticityMaterialModels.hpp:106
SmallStrainJ2Plasticity::mu
double mu
Definition: SmallStrainPlasticityMaterialModels.hpp:40
SmallStrainPlasticity::ClosestPointProjection::a_internalVariables
ublas::vector< adouble > a_internalVariables
Definition: SmallStrainPlasticity.hpp:517
SmallStrainJ2Plasticity::lambda
double lambda
Definition: SmallStrainPlasticityMaterialModels.hpp:41
SmallStrainPlasticity::ClosestPointProjection::a_sTress
ublas::vector< adouble > a_sTress
Definition: SmallStrainPlasticity.hpp:518
SmallStrainPlasticity::ClosestPointProjection::internalVariables
VectorDouble internalVariables
Definition: SmallStrainPlasticity.hpp:486
SmallStrainPlasticity::ClosestPointProjection::a_w
adouble a_w
Definition: SmallStrainPlasticity.hpp:519
adouble
FTensor::dd
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33
SmallStrainJ2Plasticity::H
double H
Definition: SmallStrainPlasticityMaterialModels.hpp:42
SmallStrainPlasticity::ClosestPointProjection::a_h
adouble a_h
Definition: SmallStrainPlasticity.hpp:519
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
SmallStrainJ2Plasticity::K
double K
Definition: SmallStrainPlasticityMaterialModels.hpp:43
SmallStrainPlasticity::ClosestPointProjection::a_sTrain
ublas::vector< adouble > a_sTrain
Definition: SmallStrainPlasticity.hpp:517
SmallStrainPlasticity::ClosestPointProjection::a_y
adouble a_y
Definition: SmallStrainPlasticity.hpp:519