v0.14.0 |
[J2 2D] More...
#include <users_modules/adolc-plasticity/src/ADOLCPlasticityMaterialModels.hpp>
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< ParamsArray > | defaultParamsArrayPtr = 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_TAPE > | tapesTags |
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< adouble > | a_sTrain |
ublas::vector< adouble > | a_plasticStrain |
ublas::vector< adouble > | a_internalVariables |
ublas::vector< adouble > | a_sTress |
ublas::vector< adouble > | a_internalFluxes |
adouble | a_w |
adouble | a_y |
adouble | a_h |
[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.
using ADOLCPlasticity::ParaboloidalPlasticity::ParamsArray = std::array<double, LAST_PARAM> |
Definition at line 570 of file ADOLCPlasticityMaterialModels.hpp.
Enumerator | |
---|---|
LAMBDA | |
MU | |
NUP | |
SIGMA_YT | |
SIGMA_YC | |
HT | |
HC | |
NT | |
NC | |
LAST_PARAM |
Definition at line 557 of file ADOLCPlasticityMaterialModels.hpp.
|
inline |
Definition at line 533 of file ADOLCPlasticityMaterialModels.hpp.
|
inlinevirtual |
Get model parameters from blocks.
Reimplemented from ADOLCPlasticity::ClosestPointProjection.
Definition at line 641 of file ADOLCPlasticityMaterialModels.hpp.
|
inline |
[Paraboloidal yield function]
Definition at line 681 of file ADOLCPlasticityMaterialModels.hpp.
|
inlinevirtual |
[Paraboloidal yield function]
[Paraboloidal flow potential]
Implements ADOLCPlasticity::ClosestPointProjection.
Definition at line 705 of file ADOLCPlasticityMaterialModels.hpp.
|
inlinevirtual |
[Paraboloidal free energy]
Implements ADOLCPlasticity::ClosestPointProjection.
Definition at line 653 of file ADOLCPlasticityMaterialModels.hpp.
|
inline |
Definition at line 573 of file ADOLCPlasticityMaterialModels.hpp.
|
inlinevirtual |
Set parameters for ADOL-C of constitutive relations.
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.
Definition at line 647 of file ADOLCPlasticityMaterialModels.hpp.
|
inlinevirtual |
Set yield function.
Implements ADOLCPlasticity::ClosestPointProjection.
Definition at line 695 of file ADOLCPlasticityMaterialModels.hpp.
boost::shared_ptr<ParamsArray> ADOLCPlasticity::ParaboloidalPlasticity::defaultParamsArrayPtr = nullptr |
Definition at line 571 of file ADOLCPlasticityMaterialModels.hpp.
double ADOLCPlasticity::ParaboloidalPlasticity::Hc |
Definition at line 553 of file ADOLCPlasticityMaterialModels.hpp.
double ADOLCPlasticity::ParaboloidalPlasticity::Ht |
Definition at line 553 of file ADOLCPlasticityMaterialModels.hpp.
FTensor::Index<'i', 3> ADOLCPlasticity::ParaboloidalPlasticity::i |
Definition at line 538 of file ADOLCPlasticityMaterialModels.hpp.
adouble ADOLCPlasticity::ParaboloidalPlasticity::I1 |
[Paraboloidal free energy]
Definition at line 678 of file ADOLCPlasticityMaterialModels.hpp.
FTensor::Index<'j', 3> ADOLCPlasticity::ParaboloidalPlasticity::j |
Definition at line 539 of file ADOLCPlasticityMaterialModels.hpp.
adouble ADOLCPlasticity::ParaboloidalPlasticity::J2 |
Definition at line 678 of file ADOLCPlasticityMaterialModels.hpp.
double ADOLCPlasticity::ParaboloidalPlasticity::lambda |
Definition at line 550 of file ADOLCPlasticityMaterialModels.hpp.
double ADOLCPlasticity::ParaboloidalPlasticity::mu |
Definition at line 551 of file ADOLCPlasticityMaterialModels.hpp.
double ADOLCPlasticity::ParaboloidalPlasticity::nc |
Definition at line 555 of file ADOLCPlasticityMaterialModels.hpp.
double ADOLCPlasticity::ParaboloidalPlasticity::nt |
Definition at line 555 of file ADOLCPlasticityMaterialModels.hpp.
double ADOLCPlasticity::ParaboloidalPlasticity::nup |
Definition at line 552 of file ADOLCPlasticityMaterialModels.hpp.
double ADOLCPlasticity::ParaboloidalPlasticity::sigmaYc |
Definition at line 554 of file ADOLCPlasticityMaterialModels.hpp.
double ADOLCPlasticity::ParaboloidalPlasticity::sigmaYt |
Definition at line 554 of file ADOLCPlasticityMaterialModels.hpp.
FTensor::Tensor2_symmetric<adouble, 3> ADOLCPlasticity::ParaboloidalPlasticity::tDeviator |
Definition at line 546 of file ADOLCPlasticityMaterialModels.hpp.
FTensor::Tensor2_symmetric<adouble, 3> ADOLCPlasticity::ParaboloidalPlasticity::tElasticStrain |
Definition at line 544 of file ADOLCPlasticityMaterialModels.hpp.
FTensor::Tensor2_symmetric<adouble, 3> ADOLCPlasticity::ParaboloidalPlasticity::tPlasticStrain |
Definition at line 543 of file ADOLCPlasticityMaterialModels.hpp.
adouble ADOLCPlasticity::ParaboloidalPlasticity::tR |
Definition at line 548 of file ADOLCPlasticityMaterialModels.hpp.
FTensor::Tensor2_symmetric<adouble, 3> ADOLCPlasticity::ParaboloidalPlasticity::tStress |
Definition at line 545 of file ADOLCPlasticityMaterialModels.hpp.
FTensor::Tensor2_symmetric<adouble, 3> ADOLCPlasticity::ParaboloidalPlasticity::tTotalStrain |
Definition at line 542 of file ADOLCPlasticityMaterialModels.hpp.
FTensor::Index<'Z', 6> ADOLCPlasticity::ParaboloidalPlasticity::Z |
Definition at line 540 of file ADOLCPlasticityMaterialModels.hpp.