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

Small strain plasticity with paraboloidal yield criterion (Isotropic Hardening) More...

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

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

Public Member Functions

 SmallStrainParaboloidalPlasticity ()
 
virtual PetscErrorCode freeHemholtzEnergy ()
 Available energy. More...
 
PetscErrorCode evalF ()
 
virtual PetscErrorCode yieldFunction ()
 Evaluate yield function. More...
 
virtual PetscErrorCode flowPotential ()
 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 nup
 
double Ht
 
double Hc
 
double sIgma_yt
 
double sIgma_yc
 
double nt
 
double nc
 
adouble I1
 Auxiliary function. More...
 
adouble J2
 
- 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

Small strain plasticity with paraboloidal yield criterion (Isotropic Hardening)

Examples
SmallStrainPlasticityMaterialModels.hpp.

Definition at line 162 of file SmallStrainPlasticityMaterialModels.hpp.

Constructor & Destructor Documentation

◆ SmallStrainParaboloidalPlasticity()

SmallStrainParaboloidalPlasticity::SmallStrainParaboloidalPlasticity ( )
inline
Examples
SmallStrainPlasticityMaterialModels.hpp.

Definition at line 164 of file SmallStrainPlasticityMaterialModels.hpp.

164  :
166  internalVariables.resize(2,false);
167  internalVariables.clear();
168  }

Member Function Documentation

◆ evalF()

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

Definition at line 269 of file SmallStrainPlasticityMaterialModels.hpp.

269  {
270  PetscFunctionBegin;
271  adouble tR = 0;
272  for(int dd = 0;dd<3;dd++) {
273  tR += a_sTress[dd];
274  }
275  I1=tR;
276 
277  tR /= 3.;
278  J2 = 0;
279  adouble t;
280  for(int dd = 0;dd<3;dd++) {
281  t = (a_sTress[dd]-tR);
282  J2 += t*t;
283  }
284  for(int dd = 3;dd<6;dd++) {
285  J2 += 2.*a_sTress[dd]*a_sTress[dd]; // s:s off diagonal terms need to go twice
286  }
287  J2=0.5*J2;
288 
289  PetscFunctionReturn(0);
290  }

◆ flowPotential()

virtual PetscErrorCode SmallStrainParaboloidalPlasticity::flowPotential ( )
inlinevirtual

Flow potential.

\[ \Psi = 6J_2 + 2\alpha I_1 \left(\overline{\alpha_1}-\overline{\alpha_0}\right) - 2\overline{\alpha_0} \,\overline{\alpha_1} \]

\[ \alpha= \frac{1-2\nu_p}{1+\nu_p} \]

Reimplemented from SmallStrainPlasticity::ClosestPointProjection.

Examples
SmallStrainPlasticityMaterialModels.hpp.

Definition at line 327 of file SmallStrainPlasticityMaterialModels.hpp.

327  {
328  PetscFunctionBegin;
329  ierr = evalF(); CHKERRQ(ierr);
330  double alpha= (1-2*nup)/(1+nup); //relation between alpha and plastic Poission's ratio
332  PetscFunctionReturn(0);
333  }

◆ freeHemholtzEnergy()

virtual PetscErrorCode SmallStrainParaboloidalPlasticity::freeHemholtzEnergy ( )
inlinevirtual

Available energy.

\[ \psi = \frac{1}{2} \lambda \textrm{tr}[\varepsilon]^2 + \mu \varepsilon : \varepsilon + \sigma_{yt}\alpha_0 + \frac{1}{2} H_t \alpha_0^2 + \sigma_{yc}\alpha_1 + \frac{1}{2} H_c \alpha_1^2 \]

Energy with linear hardening

Reimplemented from SmallStrainPlasticity::ClosestPointProjection.

Examples
SmallStrainPlasticityMaterialModels.hpp.

Definition at line 224 of file SmallStrainPlasticityMaterialModels.hpp.

224  {
225  PetscFunctionBegin;
226  //elastic part
227  adouble tR = 0;
228  for(int dd = 0;dd<3;dd++) {
229  tR += a_sTrain[dd]-a_plasticStrain[dd];
230  }
231  a_w = 0.5*lambda*tR*tR;
232  adouble t;
233  int dd = 0;
234  for(;dd<3;dd++) {
236  a_w += mu*t*t;
237  }
238  for(;dd<6;dd++) {
240  a_w += 0.5*mu*t*t;
241  }
242  //plastic part
243  //isotropic
244 // a_w += a_internalVariables[0]*sIgma_yt + 0.5*Ht*a_internalVariables[0]*a_internalVariables[0];
245 // a_w += a_internalVariables[1]*sIgma_yc + 0.5*Hc*a_internalVariables[1]*a_internalVariables[1];
248 
249  PetscFunctionReturn(0);
250  }

◆ yieldFunction()

virtual PetscErrorCode SmallStrainParaboloidalPlasticity::yieldFunction ( )
inlinevirtual

Evaluate yield function.

\[ y = 6J_2 + 2I_1\left(\overline{\alpha_1}-\overline{\alpha_0}\right) - 2\overline{\alpha_0} \,\overline{\alpha_1} \]

where

\[ \overline{\alpha_0}=\frac{\partial \psi}{\partial \alpha_0}=\sigma_{yt} + H_t \alpha_0 \]

\[ \overline{\alpha_1}=\frac{\partial \psi}{\partial \alpha_1}=\sigma_{yc} + H_c \alpha_1 \]

Reimplemented from SmallStrainPlasticity::ClosestPointProjection.

Examples
SmallStrainPlasticityMaterialModels.hpp.

Definition at line 309 of file SmallStrainPlasticityMaterialModels.hpp.

309  {
310  PetscFunctionBegin;
311  ierr = evalF(); CHKERRQ(ierr);
313  PetscFunctionReturn(0);
314  }

Member Data Documentation

◆ Hc

double SmallStrainParaboloidalPlasticity::Hc

◆ Ht

double SmallStrainParaboloidalPlasticity::Ht

◆ I1

adouble SmallStrainParaboloidalPlasticity::I1

Auxiliary function.

\[ I_1 = \textrm{tr} (\boldsymbol{\sigma}) \]

\[ \eta=\textrm{dev}[\boldsymbol{\sigma}] \]

\[ J_2 = \frac{1}{2} \eta:\eta \]

Examples
SmallStrainPlasticityMaterialModels.hpp.

Definition at line 268 of file SmallStrainPlasticityMaterialModels.hpp.

◆ J2

adouble SmallStrainParaboloidalPlasticity::J2

◆ lambda

double SmallStrainParaboloidalPlasticity::lambda

◆ mu

double SmallStrainParaboloidalPlasticity::mu

◆ nc

double SmallStrainParaboloidalPlasticity::nc

◆ nt

double SmallStrainParaboloidalPlasticity::nt

◆ nup

double SmallStrainParaboloidalPlasticity::nup

◆ sIgma_yc

double SmallStrainParaboloidalPlasticity::sIgma_yc

◆ sIgma_yt

double SmallStrainParaboloidalPlasticity::sIgma_yt

The documentation for this struct was generated from the following file:
SmallStrainParaboloidalPlasticity::J2
adouble J2
Definition: SmallStrainPlasticityMaterialModels.hpp:268
SmallStrainParaboloidalPlasticity::nt
double nt
Definition: SmallStrainPlasticityMaterialModels.hpp:175
SmallStrainPlasticity::ClosestPointProjection::a_plasticStrain
ublas::vector< adouble > a_plasticStrain
Definition: SmallStrainPlasticity.hpp:517
SmallStrainParaboloidalPlasticity::lambda
double lambda
Definition: SmallStrainPlasticityMaterialModels.hpp:171
SmallStrainPlasticity::ClosestPointProjection
Closest point projection algorithm.
Definition: SmallStrainPlasticity.hpp:482
SmallStrainPlasticity::ClosestPointProjection::a_internalFluxes
ublas::vector< adouble > a_internalFluxes
Definition: SmallStrainPlasticity.hpp:518
SmallStrainParaboloidalPlasticity::sIgma_yt
double sIgma_yt
Definition: SmallStrainPlasticityMaterialModels.hpp:174
SmallStrainParaboloidalPlasticity::nc
double nc
Definition: SmallStrainPlasticityMaterialModels.hpp:175
SmallStrainParaboloidalPlasticity::evalF
PetscErrorCode evalF()
Definition: SmallStrainPlasticityMaterialModels.hpp:269
SmallStrainPlasticity::ClosestPointProjection::a_internalVariables
ublas::vector< adouble > a_internalVariables
Definition: SmallStrainPlasticity.hpp:517
t
constexpr double t
plate stiffness
Definition: plate.cpp:59
SmallStrainParaboloidalPlasticity::I1
adouble I1
Auxiliary function.
Definition: SmallStrainPlasticityMaterialModels.hpp:268
SmallStrainPlasticity::ClosestPointProjection::a_sTress
ublas::vector< adouble > a_sTress
Definition: SmallStrainPlasticity.hpp:518
SmallStrainPlasticity::ClosestPointProjection::internalVariables
VectorDouble internalVariables
Definition: SmallStrainPlasticity.hpp:486
SmallStrainParaboloidalPlasticity::nup
double nup
Definition: SmallStrainPlasticityMaterialModels.hpp:172
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
SmallStrainParaboloidalPlasticity::mu
double mu
Definition: SmallStrainPlasticityMaterialModels.hpp:170
SmallStrainParaboloidalPlasticity::Hc
double Hc
Definition: SmallStrainPlasticityMaterialModels.hpp:173
SmallStrainPlasticity::ClosestPointProjection::a_h
adouble a_h
Definition: SmallStrainPlasticity.hpp:519
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
SmallStrainParaboloidalPlasticity::Ht
double Ht
Definition: SmallStrainPlasticityMaterialModels.hpp:173
SmallStrainParaboloidalPlasticity::sIgma_yc
double sIgma_yc
Definition: SmallStrainPlasticityMaterialModels.hpp:174
SmallStrainPlasticity::ClosestPointProjection::a_sTrain
ublas::vector< adouble > a_sTrain
Definition: SmallStrainPlasticity.hpp:517
SmallStrainPlasticity::ClosestPointProjection::a_y
adouble a_y
Definition: SmallStrainPlasticity.hpp:519