v0.14.0
Public Types | Public Member Functions | Public Attributes | List of all members
ADOLCPlasticity::ParaboloidalPlasticity Struct Reference

[J2 2D] More...

#include <users_modules/adolc-plasticity/src/ADOLCPlasticityMaterialModels.hpp>

Inheritance diagram for ADOLCPlasticity::ParaboloidalPlasticity:
[legend]
Collaboration diagram for ADOLCPlasticity::ParaboloidalPlasticity:
[legend]

Public Types

enum  ModelParams {
  LAMBDA, MU, NUP, SIGMA_YT,
  SIGMA_YC, HT, HC, NT,
  NC, LAST_PARAM
}
 
using ParamsArray = std::array< double, LAST_PARAM >
 
- Public Types inherited from ADOLCPlasticity::ClosestPointProjection
enum  TypesTags { W = 0, Y, H, LAST_TAPE }
 

Public Member Functions

 ParaboloidalPlasticity ()
 
MoFEMErrorCode getDefaultMaterialParameters ()
 
MoFEMErrorCode addMatBlockOps (MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string block_name, Sev sev)
 Get model parameters from blocks. More...
 
MoFEMErrorCode setParams (short tag, bool &recalculate_elastic_tangent)
 Set parameters for ADOL-C of constitutive relations. More...
 
MoFEMErrorCode freeHemholtzEnergy ()
 [Paraboloidal free energy] More...
 
MoFEMErrorCode evalF ()
 [Paraboloidal yield function] More...
 
MoFEMErrorCode yieldFunction ()
 Set yield function. More...
 
MoFEMErrorCode flowPotential ()
 [Paraboloidal yield function] More...
 
- Public Member Functions inherited from ADOLCPlasticity::ClosestPointProjection
VectorAdaptor getPlasticStrain ()
 
VectorAdaptor getTotalStrain ()
 
VectorAdaptor getInternalVariables ()
 
VectorAdaptor getActiveVariablesYH ()
 
VectorAdaptor getStress ()
 
VectorAdaptor getInternalFluxes ()
 
 ClosestPointProjection ()
 
MoFEMErrorCode recordW ()
 Record strain energy. More...
 
MoFEMErrorCode recordY ()
 Record yield function. More...
 
MoFEMErrorCode recordH ()
 Record flow potential. More...
 
MoFEMErrorCode playW ()
 
MoFEMErrorCode playW_NoHessian ()
 
MoFEMErrorCode playY ()
 
MoFEMErrorCode playY_NoGradient ()
 
MoFEMErrorCode playH ()
 
MoFEMErrorCode playH_NoHessian ()
 
MoFEMErrorCode createMatAVecR ()
 
MoFEMErrorCode evaluatePotentials ()
 
MoFEMErrorCode playPotentials ()
 
MoFEMErrorCode playPotentials_NoHessian ()
 
MoFEMErrorCode calculateR (Vec R)
 
MoFEMErrorCode calculateA ()
 
MoFEMErrorCode snesCreate ()
 Create nested snes. More...
 
MoFEMErrorCode solveClosetProjection ()
 Solve nonlinear system of equations to find stress, internal fluxes, and Lagrange plastic multiplier. More...
 
MoFEMErrorCode consistentTangent ()
 Calculate consistent tangent matrix. More...
 
MoFEMErrorCode recordTapes ()
 Record tapes. More...
 

Public Attributes

FTensor::Index< 'i', 3 > i
 
FTensor::Index< 'j', 3 > j
 
FTensor::Index< 'Z', 6 > Z
 
FTensor::Tensor2_symmetric< adouble, 3 > tTotalStrain
 
FTensor::Tensor2_symmetric< adouble, 3 > tPlasticStrain
 
FTensor::Tensor2_symmetric< adouble, 3 > tElasticStrain
 
FTensor::Tensor2_symmetric< adouble, 3 > tStress
 
FTensor::Tensor2_symmetric< adouble, 3 > tDeviator
 
adouble tR
 
double lambda
 
double mu
 
double nup
 
double Ht
 
double Hc
 
double sigmaYt
 
double sigmaYc
 
double nt
 
double nc
 
boost::shared_ptr< ParamsArraydefaultParamsArrayPtr = nullptr
 
adouble I1
 [Paraboloidal free energy] More...
 
adouble J2
 
- Public Attributes inherited from ADOLCPlasticity::ClosestPointProjection
int nbInternalVariables
 
boost::function< int(int, int, int)> integrationRule
 
VectorDouble internalVariables0
 
VectorDouble plasticStrain0
 
std::array< int, LAST_TAPEtapesTags
 
VectorAdaptor activeVariablesW
 
VectorAdaptor gradientW
 
double w
 
double y
 
double h
 
double deltaGamma
 
MatrixDouble Ep
 
MatrixDouble Cp
 
MatrixDouble Cep
 
ublas::symmetric_matrix< double, ublas::lower > C
 
ublas::symmetric_matrix< double, ublas::lower > D
 
MatrixDouble hessianW
 
VectorDouble gradientY
 
VectorDouble gradientH
 
MatrixDouble hessianH
 
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
 
SmartPetscObj< Mat > A
 
SmartPetscObj< Vec > R
 
SmartPetscObj< Vec > Chi
 
SmartPetscObj< Vec > dChi
 
ublas::matrix< double, ublas::column_major > dataA
 
SmartPetscObj< SNES > sNes
 
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
 

Detailed Description

[J2 2D]

Paraboloidal yield criterion

See paper for reference [71] and [68]

Linear hardening

\[ \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 \]

Exponential hardening

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

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

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

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

\[ 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 \]

\[ \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} \]

Definition at line 531 of file ADOLCPlasticityMaterialModels.hpp.

Member Typedef Documentation

◆ ParamsArray

Member Enumeration Documentation

◆ ModelParams

Enumerator
LAMBDA 
MU 
NUP 
SIGMA_YT 
SIGMA_YC 
HT 
HC 
NT 
NC 
LAST_PARAM 
Examples
ADOLCPlasticityMaterialModels.hpp.

Definition at line 557 of file ADOLCPlasticityMaterialModels.hpp.

557  {
558  LAMBDA,
559  MU,
560  NUP,
561  SIGMA_YT,
562  SIGMA_YC,
563  HT,
564  HC,
565  NT,
566  NC,
567  LAST_PARAM
568  };

Constructor & Destructor Documentation

◆ ParaboloidalPlasticity()

ADOLCPlasticity::ParaboloidalPlasticity::ParaboloidalPlasticity ( )
inline
Examples
ADOLCPlasticityMaterialModels.hpp.

Definition at line 533 of file ADOLCPlasticityMaterialModels.hpp.

535  activeVariablesW.resize(12 + nbInternalVariables, false);
536  }

Member Function Documentation

◆ addMatBlockOps()

MoFEMErrorCode ADOLCPlasticity::ParaboloidalPlasticity::addMatBlockOps ( MoFEM::Interface m_field,
boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &  pip,
std::string  block_name,
Sev  sev 
)
inlinevirtual

Get model parameters from blocks.

Reimplemented from ADOLCPlasticity::ClosestPointProjection.

Examples
ADOLCPlasticityMaterialModels.hpp.

Definition at line 641 of file ADOLCPlasticityMaterialModels.hpp.

643  {
644  return 0;
645  }

◆ evalF()

MoFEMErrorCode ADOLCPlasticity::ParaboloidalPlasticity::evalF ( )
inline

[Paraboloidal yield function]

Examples
ADOLCPlasticityMaterialModels.hpp.

Definition at line 681 of file ADOLCPlasticityMaterialModels.hpp.

681  {
683  auto t_voight_stress_op = voight_to_stress_op();
684  auto t_vioght_stress = getFTensor1FromPtr<6>(&a_sTress[0]);
685  tStress(i, j) = t_voight_stress_op(i, j, Z) * t_vioght_stress(Z);
686  tR = tStress(i, i);
687  I1 = tR;
689  tR /= 3.;
690  tDeviator(i, j) = tStress(i, j) - tR * t_kd(i, j);
691  J2 = 0.5 * (tDeviator(i, j) * tDeviator(i, j));
693  }

◆ flowPotential()

MoFEMErrorCode ADOLCPlasticity::ParaboloidalPlasticity::flowPotential ( )
inlinevirtual

[Paraboloidal yield function]

[Paraboloidal flow potential]

Implements ADOLCPlasticity::ClosestPointProjection.

Examples
ADOLCPlasticityMaterialModels.hpp.

Definition at line 705 of file ADOLCPlasticityMaterialModels.hpp.

705  {
707  CHKERR evalF();
708  double alpha =
709  (1 - 2 * nup) /
710  (1 + nup); // relation between alpha and plastic poisson's ratio
711  a_h = 6 * J2 +
712  2 * alpha * I1 * (a_internalFluxes[1] - a_internalFluxes[0]) -
715  }

◆ freeHemholtzEnergy()

MoFEMErrorCode ADOLCPlasticity::ParaboloidalPlasticity::freeHemholtzEnergy ( )
inlinevirtual

[Paraboloidal free energy]

Implements ADOLCPlasticity::ClosestPointProjection.

Examples
ADOLCPlasticityMaterialModels.hpp.

Definition at line 653 of file ADOLCPlasticityMaterialModels.hpp.

653  {
655 
656  auto t_voight_strain_op = voight_to_strain_op();
657  auto t_voight_total_strain = getFTensor1FromPtr<6>(&a_sTrain[0]);
658  auto t_voight_plastic_strain = getFTensor1FromPtr<6>(&a_plasticStrain[0]);
659 
660  tTotalStrain(i, j) = t_voight_strain_op(i, j, Z) * t_voight_total_strain(Z);
661  tPlasticStrain(i, j) =
662  t_voight_strain_op(i, j, Z) * t_voight_plastic_strain(Z);
664  tR = tElasticStrain(i, i);
665 
666  a_w = 0.5 * lambda * tR * tR +
667  mu * (tElasticStrain(i, j) * tElasticStrain(i, j));
668 
673 
675  }

◆ getDefaultMaterialParameters()

MoFEMErrorCode ADOLCPlasticity::ParaboloidalPlasticity::getDefaultMaterialParameters ( )
inline
Examples
ADOLCPlasticityMaterialModels.hpp.

Definition at line 573 of file ADOLCPlasticityMaterialModels.hpp.

573  {
575 
576  double young_modulus = 1;
577  double poisson_ratio = 0.25;
578 
579  sigmaYt = 1;
580  sigmaYc = 1;
581  Ht = 0.1;
582  Hc = 0.1;
583  nc = 1.;
584  nt = 1.;
585  nup = 0.;
586 
587  CHKERR PetscOptionsGetScalar(PETSC_NULL, "-young_modulus", &young_modulus,
588  PETSC_NULL);
589  CHKERR PetscOptionsGetScalar(PETSC_NULL, "-poisson_ratio", &poisson_ratio,
590  0);
591  CHKERR PetscOptionsGetScalar(PETSC_NULL, "-sigma_yt", &sigmaYt, PETSC_NULL);
592  CHKERR PetscOptionsGetScalar(PETSC_NULL, "-sigma_yc", &sigmaYc, PETSC_NULL);
593  CHKERR PetscOptionsGetScalar(PETSC_NULL, "-Ht", &Ht, PETSC_NULL);
594  CHKERR PetscOptionsGetScalar(PETSC_NULL, "-Hc", &Hc, PETSC_NULL);
595  CHKERR PetscOptionsGetScalar(PETSC_NULL, "-nt", &nt, PETSC_NULL);
596  CHKERR PetscOptionsGetScalar(PETSC_NULL, "-nc", &nc, PETSC_NULL);
597  CHKERR PetscOptionsGetScalar(PETSC_NULL, "-nup", &nup, PETSC_NULL);
598 
599  MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
600  << "Young modulus: " << young_modulus;
601  MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
602  << "Poisson ratio: " << poisson_ratio;
603 
606 
607  defaultParamsArrayPtr = boost::make_shared<ParamsArray>();
608  (*defaultParamsArrayPtr)[LAMBDA] = lambda;
609  (*defaultParamsArrayPtr)[MU] = mu;
610  (*defaultParamsArrayPtr)[NUP] = nup;
611  (*defaultParamsArrayPtr)[SIGMA_YT] = sigmaYt;
612  (*defaultParamsArrayPtr)[SIGMA_YC] = sigmaYc;
613  (*defaultParamsArrayPtr)[HT] = Ht;
614  (*defaultParamsArrayPtr)[HC] = Hc;
615  (*defaultParamsArrayPtr)[NT] = nt;
616  (*defaultParamsArrayPtr)[NC] = nc;
617 
618  MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
619  << "LAMBDA: " << (*defaultParamsArrayPtr)[LAMBDA];
620  MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
621  << "MU: " << (*defaultParamsArrayPtr)[MU];
622  MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
623  << "NUP: " << (*defaultParamsArrayPtr)[NUP];
624  MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
625  << "SIGMA_YT: " << (*defaultParamsArrayPtr)[SIGMA_YT];
626  MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
627  << "SIGMA_YC: " << (*defaultParamsArrayPtr)[SIGMA_YC];
628  MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
629  << "HT: " << (*defaultParamsArrayPtr)[HT];
630  MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
631  << "HC: " << (*defaultParamsArrayPtr)[HC];
632  MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
633  << "NT: " << (*defaultParamsArrayPtr)[NT];
634  MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
635  << "NC: " << (*defaultParamsArrayPtr)[NC];
636 
638  };

◆ setParams()

MoFEMErrorCode ADOLCPlasticity::ParaboloidalPlasticity::setParams ( short  tag,
bool recalculate_elastic_tangent 
)
inlinevirtual

Set parameters for ADOL-C of constitutive relations.

Parameters
tag[in] - tag of the tape
recalculate_elastic_tangent[out] - if setParam set it to true, tangent matrix for elastic domain should be recalculated

Reimplemented from ADOLCPlasticity::ClosestPointProjection.

Examples
ADOLCPlasticityMaterialModels.hpp.

Definition at line 647 of file ADOLCPlasticityMaterialModels.hpp.

647  {
648  recalculate_elastic_tangent = true;
649  return 0;
650  }

◆ yieldFunction()

MoFEMErrorCode ADOLCPlasticity::ParaboloidalPlasticity::yieldFunction ( )
inlinevirtual

Set yield function.

Implements ADOLCPlasticity::ClosestPointProjection.

Examples
ADOLCPlasticityMaterialModels.hpp.

Definition at line 695 of file ADOLCPlasticityMaterialModels.hpp.

695  {
697  CHKERR evalF();
698  a_y = 6 * J2 + 2 * I1 * (a_internalFluxes[1] - a_internalFluxes[0]) -
701  }

Member Data Documentation

◆ defaultParamsArrayPtr

boost::shared_ptr<ParamsArray> ADOLCPlasticity::ParaboloidalPlasticity::defaultParamsArrayPtr = nullptr

◆ Hc

double ADOLCPlasticity::ParaboloidalPlasticity::Hc

◆ Ht

double ADOLCPlasticity::ParaboloidalPlasticity::Ht

◆ i

FTensor::Index<'i', 3> ADOLCPlasticity::ParaboloidalPlasticity::i

◆ I1

adouble ADOLCPlasticity::ParaboloidalPlasticity::I1

[Paraboloidal free energy]

Examples
ADOLCPlasticityMaterialModels.hpp.

Definition at line 678 of file ADOLCPlasticityMaterialModels.hpp.

◆ j

FTensor::Index<'j', 3> ADOLCPlasticity::ParaboloidalPlasticity::j

◆ J2

adouble ADOLCPlasticity::ParaboloidalPlasticity::J2

◆ lambda

double ADOLCPlasticity::ParaboloidalPlasticity::lambda

◆ mu

double ADOLCPlasticity::ParaboloidalPlasticity::mu

◆ nc

double ADOLCPlasticity::ParaboloidalPlasticity::nc

◆ nt

double ADOLCPlasticity::ParaboloidalPlasticity::nt

◆ nup

double ADOLCPlasticity::ParaboloidalPlasticity::nup

◆ sigmaYc

double ADOLCPlasticity::ParaboloidalPlasticity::sigmaYc

◆ sigmaYt

double ADOLCPlasticity::ParaboloidalPlasticity::sigmaYt

◆ tDeviator

FTensor::Tensor2_symmetric<adouble, 3> ADOLCPlasticity::ParaboloidalPlasticity::tDeviator

◆ tElasticStrain

FTensor::Tensor2_symmetric<adouble, 3> ADOLCPlasticity::ParaboloidalPlasticity::tElasticStrain

◆ tPlasticStrain

FTensor::Tensor2_symmetric<adouble, 3> ADOLCPlasticity::ParaboloidalPlasticity::tPlasticStrain

◆ tR

adouble ADOLCPlasticity::ParaboloidalPlasticity::tR

◆ tStress

FTensor::Tensor2_symmetric<adouble, 3> ADOLCPlasticity::ParaboloidalPlasticity::tStress

◆ tTotalStrain

FTensor::Tensor2_symmetric<adouble, 3> ADOLCPlasticity::ParaboloidalPlasticity::tTotalStrain

◆ Z

FTensor::Index<'Z', 6> ADOLCPlasticity::ParaboloidalPlasticity::Z

The documentation for this struct was generated from the following file:
ADOLCPlasticity::ClosestPointProjection::a_h
adouble a_h
Definition: ADOLCPlasticity.hpp:330
ADOLCPlasticity::ParaboloidalPlasticity::defaultParamsArrayPtr
boost::shared_ptr< ParamsArray > defaultParamsArrayPtr
Definition: ADOLCPlasticityMaterialModels.hpp:571
young_modulus
double young_modulus
Young modulus.
Definition: plastic.cpp:172
ADOLCPlasticity::ParaboloidalPlasticity::tElasticStrain
FTensor::Tensor2_symmetric< adouble, 3 > tElasticStrain
Definition: ADOLCPlasticityMaterialModels.hpp:544
ADOLCPlasticity::ParaboloidalPlasticity::SIGMA_YC
@ SIGMA_YC
Definition: ADOLCPlasticityMaterialModels.hpp:562
ADOLCPlasticity::ClosestPointProjection::ClosestPointProjection
ClosestPointProjection()
Definition: ClosetPointProjection.cpp:15
ADOLCPlasticity::ClosestPointProjection::a_plasticStrain
ublas::vector< adouble > a_plasticStrain
Definition: ADOLCPlasticity.hpp:325
ADOLCPlasticity::ParaboloidalPlasticity::HT
@ HT
Definition: ADOLCPlasticityMaterialModels.hpp:563
ADOLCPlasticity::ParaboloidalPlasticity::nt
double nt
Definition: ADOLCPlasticityMaterialModels.hpp:555
ADOLCPlasticity::ParaboloidalPlasticity::tR
adouble tR
Definition: ADOLCPlasticityMaterialModels.hpp:548
ADOLCPlasticity::ParaboloidalPlasticity::j
FTensor::Index< 'j', 3 > j
Definition: ADOLCPlasticityMaterialModels.hpp:539
ADOLCPlasticity::ParaboloidalPlasticity::Ht
double Ht
Definition: ADOLCPlasticityMaterialModels.hpp:553
ADOLCPlasticity::ParaboloidalPlasticity::sigmaYt
double sigmaYt
Definition: ADOLCPlasticityMaterialModels.hpp:554
ADOLCPlasticity::ParaboloidalPlasticity::Hc
double Hc
Definition: ADOLCPlasticityMaterialModels.hpp:553
ADOLCPlasticity::ClosestPointProjection::a_sTrain
ublas::vector< adouble > a_sTrain
Definition: ADOLCPlasticity.hpp:324
ADOLCPlasticity::ParaboloidalPlasticity::NUP
@ NUP
Definition: ADOLCPlasticityMaterialModels.hpp:560
ADOLCPlasticity::ParaboloidalPlasticity::mu
double mu
Definition: ADOLCPlasticityMaterialModels.hpp:551
ADOLCPlasticity::ParaboloidalPlasticity::NT
@ NT
Definition: ADOLCPlasticityMaterialModels.hpp:565
ADOLCPlasticity::ParaboloidalPlasticity::nc
double nc
Definition: ADOLCPlasticityMaterialModels.hpp:555
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
ADOLCPlasticity::ParaboloidalPlasticity::evalF
MoFEMErrorCode evalF()
[Paraboloidal yield function]
Definition: ADOLCPlasticityMaterialModels.hpp:681
ADOLCPlasticity::ParaboloidalPlasticity::I1
adouble I1
[Paraboloidal free energy]
Definition: ADOLCPlasticityMaterialModels.hpp:678
ADOLCPlasticity::ParaboloidalPlasticity::i
FTensor::Index< 'i', 3 > i
Definition: ADOLCPlasticityMaterialModels.hpp:538
poisson_ratio
double poisson_ratio
Poisson ratio.
Definition: plastic.cpp:173
ADOLCPlasticity::voight_to_strain_op
FTensor::Dg< double, 3, 6 > voight_to_strain_op()
Op convert Vight strain vector to strain tensor.
Definition: ADOLCPlasticity.hpp:29
ADOLCPlasticity::ParaboloidalPlasticity::Z
FTensor::Index< 'Z', 6 > Z
Definition: ADOLCPlasticityMaterialModels.hpp:540
ADOLCPlasticity::ParaboloidalPlasticity::tPlasticStrain
FTensor::Tensor2_symmetric< adouble, 3 > tPlasticStrain
Definition: ADOLCPlasticityMaterialModels.hpp:543
ADOLCPlasticity::ClosestPointProjection::a_w
adouble a_w
Definition: ADOLCPlasticity.hpp:328
ADOLCPlasticity::ParaboloidalPlasticity::nup
double nup
Definition: ADOLCPlasticityMaterialModels.hpp:552
ADOLCPlasticity::ClosestPointProjection::a_internalFluxes
ublas::vector< adouble > a_internalFluxes
Definition: ADOLCPlasticity.hpp:327
t_kd
constexpr auto t_kd
Definition: free_surface.cpp:137
ADOLCPlasticity::ClosestPointProjection::nbInternalVariables
int nbInternalVariables
Definition: ADOLCPlasticity.hpp:141
ADOLCPlasticity::ClosestPointProjection::activeVariablesW
VectorAdaptor activeVariablesW
Definition: ADOLCPlasticity.hpp:173
ADOLCPlasticity::ParaboloidalPlasticity::tDeviator
FTensor::Tensor2_symmetric< adouble, 3 > tDeviator
Definition: ADOLCPlasticityMaterialModels.hpp:546
ADOLCPlasticity::ParaboloidalPlasticity::LAST_PARAM
@ LAST_PARAM
Definition: ADOLCPlasticityMaterialModels.hpp:567
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
ADOLCPlasticity::ParaboloidalPlasticity::tTotalStrain
FTensor::Tensor2_symmetric< adouble, 3 > tTotalStrain
Definition: ADOLCPlasticityMaterialModels.hpp:542
ADOLCPlasticity::ClosestPointProjection::a_internalVariables
ublas::vector< adouble > a_internalVariables
Definition: ADOLCPlasticity.hpp:326
ADOLCPlasticity::ClosestPointProjection::a_y
adouble a_y
Definition: ADOLCPlasticity.hpp:329
ADOLCPlasticity::ParaboloidalPlasticity::LAMBDA
@ LAMBDA
Definition: ADOLCPlasticityMaterialModels.hpp:558
ADOLCPlasticity::voight_to_stress_op
FTensor::Dg< double, 3, 6 > voight_to_stress_op()
Op convert Vight stress vector to stress tensor.
Definition: ADOLCPlasticity.hpp:49
ADOLCPlasticity::ClosestPointProjection::a_sTress
ublas::vector< adouble > a_sTress
Definition: ADOLCPlasticity.hpp:327
ADOLCPlasticity::ParaboloidalPlasticity::MU
@ MU
Definition: ADOLCPlasticityMaterialModels.hpp:559
ADOLCPlasticity::ParaboloidalPlasticity::sigmaYc
double sigmaYc
Definition: ADOLCPlasticityMaterialModels.hpp:554
FTensor::Kronecker_Delta_symmetric
Kronecker Delta class symmetric.
Definition: Kronecker_Delta.hpp:49
ADOLCPlasticity::ParaboloidalPlasticity::SIGMA_YT
@ SIGMA_YT
Definition: ADOLCPlasticityMaterialModels.hpp:561
MoFEM::PetscOptionsGetScalar
PetscErrorCode PetscOptionsGetScalar(PetscOptions *, const char pre[], const char name[], PetscScalar *dval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:162
ADOLCPlasticity::ParaboloidalPlasticity::NC
@ NC
Definition: ADOLCPlasticityMaterialModels.hpp:566
ADOLCPlasticity::ParaboloidalPlasticity::HC
@ HC
Definition: ADOLCPlasticityMaterialModels.hpp:564
ADOLCPlasticity::ParaboloidalPlasticity::tStress
FTensor::Tensor2_symmetric< adouble, 3 > tStress
Definition: ADOLCPlasticityMaterialModels.hpp:545
ADOLCPlasticity::ParaboloidalPlasticity::lambda
double lambda
Definition: ADOLCPlasticityMaterialModels.hpp:550
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
ADOLCPlasticity::ParaboloidalPlasticity::J2
adouble J2
Definition: ADOLCPlasticityMaterialModels.hpp:678