Finite-strain J2 plasticity
Theory
This model is a finite-strain hyperelastic-plastic J2 plasticity formulation with isotropic hardening. The strain energy potential is based on the compressible Neo-Hookean model (Pence–Gou, variant B).
We assume the multiplicative split of the deformation gradient into elastic and plastic parts,
The Pence–Gou potential (variant B) applied to \(\boldsymbol{C}^{\rm e}\) reads:
From the potential we obtain the second Piola–Kirchhoff stress (elastic configuration)
and the Mandel stress in the intermediate configuration:
Plastic admissibility is enforced by a J2 yield function in Mandel stress space with a stress-like hardening variable \(\beta_{\rm p}(\alpha_{\rm p})\):
where \(\boldsymbol{M}^{\rm dev} = \boldsymbol{M} - \frac{1}{3}\operatorname{tr}(\boldsymbol{M})\boldsymbol{I}\) is the deviatoric part of the Mandel stress.
The hardening function is:
with a strain-like hardening internal variable \(\alpha_{\rm p}\).
An associated flow is used:
The Karush-Kuhn-Tucker (KKT) conditions for plastic loading are:
The Kirchhoff stress follows by push-forward from \(\boldsymbol{S}\) with \(\boldsymbol{F}^{\rm e}\):
Stress update at a material point proceeds by computing the trial elastic split \(\boldsymbol{F}^{\rm e,tr}=\boldsymbol{F}\,\big(\boldsymbol{F}^{\rm p,old}\big)^{-1}\) and evaluating the yield function \(f\). If plastic yielding occurs (\(f > 0\)) a fully implicit return-mapping is solved for the unknown vector
A residual system \(\boldsymbol{R}(\boldsymbol{X})=\boldsymbol{0}\) enforces (i) elastic–plastic kinematics, (ii) hardening update, and (iii) consistency \(f=0\). The plastic update of \(\boldsymbol{F}^{\rm p}\) uses the matrix exponential of the flow direction in Mandel-space:
and the strain-like hardening variable is updated as:
Upon convergence of the return-mapping algorithm, the plastic variables are updated according to:
The elastic deformation gradient follows from the multiplicative split as \(\boldsymbol{F}^{\rm e,new}=\boldsymbol{F}\,(\boldsymbol{F}^{\rm p,new})^{-1}\), and the Kirchhoff stress and consistent algorithmic tangent are then computed from the hyperelastic constitutive law using \(\boldsymbol{F}^{\rm e,new}\).
The consistent tangent \(\frac{\partial\boldsymbol{\tau}}{\partial\boldsymbol{F}}\) is computed using the chain rule as:
where the individual terms are:
For plastic step, \(\frac{\partial\boldsymbol{X}}{\partial\boldsymbol{F}}\) is computed by solving the linear system \(\frac{\partial\boldsymbol{R}}{\partial\boldsymbol{X}} \frac{\partial\boldsymbol{X}}{\partial\boldsymbol{F}} = -\frac{\partial\boldsymbol{R}}{\partial\boldsymbol{F}}\) using the converged Jacobian from the return mapping. Then \(\frac{\partial\boldsymbol{F}^{\rm e}}{\partial\boldsymbol{F}}\) is extracted from the first 9 rows of \(\frac{\partial\boldsymbol{X}}{\partial\boldsymbol{F}}\).
Stress update algorithm at a quadrature point for current step \(n+1\)
Notation: \((\cdot)^{\rm old} := (\cdot)_n\), \((\cdot)^{\rm new} := (\cdot)_{n+1}\).
Input (known quantities):
Current deformation gradient: \(\boldsymbol{F}\)
Plastic deformation gradient from previous step: \(\boldsymbol{F}^{\rm p,old}\)
Strain-like hardening variable from previous step: \(\alpha_{\rm p}^{\rm old}\)
Trial elastic state:
Compute trial elastic deformation gradient: \(\boldsymbol{F}^{\rm e,tr}=\boldsymbol{F}\,\big(\boldsymbol{F}^{\rm p,old}\big)^{-1}\)
Compute trial Mandel stress: \(\boldsymbol{M}^{\rm tr} = \boldsymbol{C}^{\rm e,tr} \boldsymbol{S}^{\rm tr}\)
Compute stress-like hardening variable: \(\beta_{\rm p} = f_{\rm y}^\infty + (f_{\rm y}^0 - f_{\rm y}^\infty) e^{-\eta \alpha_{\rm p}^{\rm old}} + H \alpha_{\rm p}^{\rm old}\)
Evaluate yield function: \(f^{\rm tr} = \frac{1}{f_{\rm y}^0}\left(||\boldsymbol{M}^{\rm tr,dev}|| - \sqrt{\frac{2}{3}} \beta_{\rm p}\right)\)
If \(f^{\rm tr} \leq 0\) (elastic step):
Accept the trial state:
\(\boldsymbol{F}^{\rm e,new} = \boldsymbol{F}^{\rm e,tr}\)
\(\boldsymbol{F}^{\rm p,new} = \boldsymbol{F}^{\rm p,old}\)
\(\alpha_{\rm p}^{\rm new} = \alpha_{\rm p}^{\rm old}\)
Compute Kirchhoff stress \(\boldsymbol{\tau}\) and consistent tangent \(\frac{\partial\boldsymbol{\tau}}{\partial\boldsymbol{F}}\) using \(\boldsymbol{F}^{\rm e,new}\)
Return
Return mapping:
If \(f^{\rm tr} > 0\) (plastic step):
Initialize unknowns: \(\boldsymbol{X} = \{\boldsymbol{F}^{\rm e}_{11}, \boldsymbol{F}^{\rm e}_{12}, \ldots, \boldsymbol{F}^{\rm e}_{33}, \alpha_{\rm p}, \Delta\lambda\}^{\mathsf{T}}\)
Solve nonlinear system \(\boldsymbol{R}(\boldsymbol{X}) = \boldsymbol{0}\) by Newton-Raphson iteration:
Residual equations:
\(\boldsymbol{R}_1\): Elastic-plastic kinematics: \(\boldsymbol{F}^{\rm e}\,\Delta\boldsymbol{F}^{\rm p} - \boldsymbol{F}\,\big(\boldsymbol{F}^{\rm p,old}\big)^{-1} = \boldsymbol{0}\)
\(\boldsymbol{R}_2\): Hardening evolution: \(\alpha_{\rm p} - \alpha_{\rm p}^{\rm old} - \frac{\Delta\lambda}{f_{\rm y}^0} \sqrt{\frac{2}{3}} = 0\)
\(\boldsymbol{R}_3\): Consistency condition: \(f(\boldsymbol{M}, \beta_{\rm p}) = 0\)
Newton-Raphson iteration: \(\boldsymbol{X}^{(k+1)} = \boldsymbol{X}^{(k)} - \left[\frac{\partial \boldsymbol{R}}{\partial \boldsymbol{X}}\right]^{-1} \boldsymbol{R}(\boldsymbol{X}^{(k)})\)
Convergence criteria: \(||\boldsymbol{R}|| < \text{TOL}\) and \(||\Delta\boldsymbol{X}|| < \text{TOL}\)
Upon convergence, update plastic variables:
\(\Delta\boldsymbol{F}^{\rm p} = \exp\left(\Delta\lambda \frac{\partial f}{\partial \boldsymbol{M}}\right)\)
\(\boldsymbol{F}^{\rm p,new} = \Delta\boldsymbol{F}^{\rm p} \boldsymbol{F}^{\rm p,old}\)
\(\alpha_{\rm p}^{\rm new} = \alpha_{\rm p}^{\rm old} + \frac{\Delta\lambda}{f_{\rm y}^0} \sqrt{\frac{2}{3}}\)
Compute updated elastic deformation gradient: \(\boldsymbol{F}^{\rm e,new} = \boldsymbol{F}\,(\boldsymbol{F}^{\rm p,new})^{-1}\)
Compute Kirchhoff stress \(\boldsymbol{\tau}\) and consistent tangent \(\frac{\partial\boldsymbol{\tau}}{\partial\boldsymbol{F}}\) using \(\boldsymbol{F}^{\rm e,new}\)
Return
Primary reference: A. Dummer, M. Neuner, P. Gamnitzer, G. Hofstetter (2024). Robust and efficient implementation of finite strain generalized continuum models for material failure: Analytical, numerical, and automatic differentiation with hyper-dual numbers. Computer Methods in Applied Mechanics and Engineering 426:116987.
Implementation
-
class FiniteStrainJ2Plasticity : public MarmotMaterialFiniteStrain
Inheritance diagram for Marmot::Materials::FiniteStrainJ2Plasticity:
![digraph {
graph [bgcolor="#00000000"]
node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2]
edge [color="#1414CE"]
"1" [label="Marmot::Materials::FiniteStrainJ2Plasticity" tooltip="Marmot::Materials::FiniteStrainJ2Plasticity" fillcolor="#BFBFBF"]
"3" [label="MarmotMaterial" tooltip="MarmotMaterial"]
"2" [label="MarmotMaterialFiniteStrain" tooltip="MarmotMaterialFiniteStrain"]
"1" -> "2" [dir=forward tooltip="public-inheritance"]
"2" -> "3" [dir=forward tooltip="public-inheritance"]
}](../../_images/graphviz-a74f16411f7b436a7131e27ba2373521ae7c7301.png)
Collaboration diagram for Marmot::Materials::FiniteStrainJ2Plasticity:
![digraph {
graph [bgcolor="#00000000"]
node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2]
edge [color="#1414CE"]
"1" [label="Marmot::Materials::FiniteStrainJ2Plasticity" tooltip="Marmot::Materials::FiniteStrainJ2Plasticity" fillcolor="#BFBFBF"]
"3" [label="MarmotMaterial" tooltip="MarmotMaterial"]
"2" [label="MarmotMaterialFiniteStrain" tooltip="MarmotMaterialFiniteStrain"]
"1" -> "2" [dir=forward tooltip="public-inheritance"]
"2" -> "3" [dir=forward tooltip="public-inheritance"]
}](../../_images/graphviz-a74f16411f7b436a7131e27ba2373521ae7c7301.png)
Finite-strain hyperelastic-plastic J2 model with isotropic hardening.
Elastic law: hyperelastic, compressible Neo-Hookean (Pence–Gou potential, variant B).
- Material parameters
- State variables
Fp — plastic deformation gradient
alphaP — strain-like hardening variable
- Implementation variants (implementationType)
0: scalar return mapping (not yet implemented)
1: Full return mapping; all derivatives computed analytically
2: FDAF — Full return mapping; derivatives approximated via forward finite differences
3: FDAC — Full return mapping; derivatives approximated via central finite differences
4: CSDA — Full return mapping; derivatives via complex-step differentiation
Public Functions
-
FiniteStrainJ2Plasticity(const double *materialProperties, int nMaterialProperties, int materialLabel)
Construct the finite-strain J2 plasticity model.
-
virtual void computeStress(ConstitutiveResponse<3> &response, AlgorithmicModuli<3> &tangents, const Deformation<3> &deformation, const TimeIncrement &timeIncrement)
Compute the Kirchhoff stress and the algorithmic tangent for the current step. Performs an elastic trial; if yielding occurs, return mapping is performed.
Template parameter
<3>indicates 3D.- Parameters:
response – [inout]
tau- Kirchhoff stress tensor \(\boldsymbol{\tau}\).elasticEnergyDensity- elastic energy density \(\psi\).rho- density (from material parameters if provided).
tangents – [inout]
dTau_dF- algorithmic tangent \(\partial\boldsymbol{\tau}/\partial\boldsymbol{F}\).
deformation – [in]
F- deformation gradient \(\boldsymbol{F}\).
timeIncrement – [in]
t- old (pseudo-)time.dT-(pseudo-)time increment.
-
void computeStressWithScalarReturnMapping(ConstitutiveResponse<3> &response, AlgorithmicModuli<3> &tangents, const Deformation<3> &deformation, const TimeIncrement &timeIncrement)
Scalar-return mapping variant; would solve a 1D consistency equation for \(\Delta\lambda\).
Note
Not yet implemented.
-
void computeStressWithFullReturnMapping(ConstitutiveResponse<3> &response, AlgorithmicModuli<3> &tangents, const Deformation<3> &deformation, const TimeIncrement &timeIncrement)
Full return mapping; all derivatives computed analytically.
-
void computeStressFDAF(ConstitutiveResponse<3> &response, AlgorithmicModuli<3> &tangents, const Deformation<3> &deformation, const TimeIncrement &timeIncrement)
Full return mapping; derivatives via forward finite differences (FDAF). Approximates \(\partial \boldsymbol{R}/\partial \boldsymbol{X}\) and \(\partial^{2}\Psi/\partial \boldsymbol{C}^{\mathrm e}\,\partial \boldsymbol{C}^{\mathrm e}\).
-
void computeStressFDAC(ConstitutiveResponse<3> &response, AlgorithmicModuli<3> &tangents, const Deformation<3> &deformation, const TimeIncrement &timeIncrement)
Full return mapping; derivatives via central finite differences (FDAC). Approximates \(\partial \boldsymbol{R}/\partial \boldsymbol{X}\) and \(\partial^{2}\Psi/\partial \boldsymbol{C}^{\mathrm e}\,\partial \boldsymbol{C}^{\mathrm e}\).
-
void computeStressCSDA(ConstitutiveResponse<3> &response, AlgorithmicModuli<3> &tangents, const Deformation<3> &deformation, const TimeIncrement &timeIncrement)
Full return mapping; derivatives via complex-step differentiation approximation (CSDA). Evaluates \(\partial \boldsymbol{R}/\partial \boldsymbol{X}\) and \(\partial^{2}\Psi/\partial \boldsymbol{C}^{\mathrm e}\,\partial \boldsymbol{C}^{\mathrm e}\) using complex-step.
-
inline virtual int getNumberOfRequiredStateVars()
Get the number of required state variables.
- Returns:
10 (9 for
Fpand 1 foralphaP).
-
inline virtual double getDensity()
Return the material density if provided in material parameters.
-
virtual void assignStateVars(double *stateVars, int nStateVars)
Bind external storage for internal variables (
Fp,alphaP).Creates/updates the internal state manager that maps (
Fp,alphaP) onto this buffer.- Parameters:
stateVars – Pointer to a contiguous array provided by the caller.
nStateVars – Number of entries in that array.
- Throws:
std::invalid_argument – If nStateVars < getNumberOfRequiredStateVars().
-
virtual StateView getStateView(const std::string &result)
Expose a named view into the state vector.
- Parameters:
result – One of:
Fp(length 9),alphaP(length 1).- Returns:
A StateView object that provides access to the requested state variable; returns an empty view if
resultis unknown.
-
virtual void initializeYourself()
Initialize state (sets \(\boldsymbol F^{\mathrm p} = \boldsymbol I\))
-
inline std::tuple<double, Tensor33d, double> yieldFunction(const Tensor33d &Fe, const double betaP)
Yield function \(f(\boldsymbol F^{\mathrm e},\beta_{\mathrm p})\) and first derivatives.
- Parameters:
Fe – Elastic deformation gradient \(\boldsymbol F^{\mathrm e}\).
betaP – Stress-like hardening variable \(\beta_{\mathrm p}\).
- Returns:
\(\{\,f,\;\partial f/\partial \boldsymbol F^{\mathrm e},\;\partial f/\partial \beta_{\mathrm p}\,\}\).
-
template<typename T>
inline T yieldFunctionFromStress(const Tensor33t<T> &mandelStress, const T betaP) Yield function \(f(\boldsymbol M,\beta_p)\) w.r.t. Mandel stress (templated).
- Template Parameters:
T – Scalar type (double, hyper-dual, complex).
- Parameters:
mandelStress – Mandel stress \(\boldsymbol M\).
betaP – Stress-like hardening variable \(\beta_{\mathrm p}\).
- Returns:
\(f\).
-
inline std::tuple<double, Tensor33d, Tensor3333d, double> yieldFunctionFromStress(const Tensor33d &mandelStress, const double betaP)
Yield function w.r.t. Mandel stress (double version).
- Parameters:
mandelStress – Mandel stress \(\boldsymbol M\).
betaP – Stress-like hardening variable \(\beta_{\mathrm p}\).
- Returns:
\(\{\,f,\;\partial f/\partial \boldsymbol M,\;\partial^2 f/\partial \boldsymbol M\partial \boldsymbol M,\;\partial f/\partial \beta_{\mathrm p}\,\}\).
-
template<typename T>
inline std::tuple<T, Tensor33t<T>, T> yieldFunctionFromStressFirstOrderDerived(const Tensor33t<T> &mandelStress, const T betaP) First-order derivatives of \(f(\boldsymbol M,\beta_{\mathrm p})\) w.r.t. \(\boldsymbol M \).
- Template Parameters:
T – Scalar type (double, hyper-dual, complex).
- Parameters:
mandelStress – Mandel stress \(\boldsymbol M\).
betaP – Stress-like hardening variable \(\beta_{\mathrm p}\).
- Returns:
\(\{\,f,\;\partial f/\partial \boldsymbol M,\;\partial f/\partial \beta_{\mathrm p}\,\}\).
-
inline bool isYielding(const Tensor33d &Fe, const double betaP)
Check yield condition \(f(\boldsymbol F^{\mathrm e},\beta_{\mathrm p})>0\).
- Parameters:
Fe – Elastic deformation gradient \(\boldsymbol F^{\mathrm e}\).
betaP – Stress-like hardening variable \(\beta_{\mathrm p}\).
- Returns:
True if yielding occurs, i.e. \(f>0\).
-
inline std::tuple<Tensor33d, Tensor3333d> computeMandelStress(const Tensor33d &Fe)
Compute Mandel stress \(\boldsymbol M\) and its derivative w.r.t. \(\boldsymbol F^{\mathrm e}\).
- Parameters:
Fe – Elastic deformation gradient \(\boldsymbol F^{\mathrm e}\).
- Returns:
\(\{\,\boldsymbol M,\;\partial \boldsymbol M/\partial \boldsymbol F^{\mathrm e}\,\}\).
-
template<typename T>
inline Tensor33t<T> computeMandelStressOnly(const Tensor33t<T> &Fe) Compute Mandel stress only (templated scalar type).
- Template Parameters:
T – Scalar type (double, hyper-dual, complex).
- Parameters:
Fe – Elastic deformation gradient \(\boldsymbol F^{\mathrm e}\).
- Returns:
Mandel stress \(\boldsymbol M\).
-
template<typename T>
inline T computeBetaPOnly(const T alphaP) Isotropic hardening law.
- Parameters:
alphaP – Strain-like hardening variable \(\alpha_{\mathrm p}\).
- Returns:
Stress-like hardening variable \(\beta_{\mathrm p}\).
-
inline std::tuple<double, double> computeBetaP(const double alphaP)
Isotropic hardening law and its derivative.
- Parameters:
alphaP – Strain-like hardening variable \(\alpha_{\mathrm p}\).
- Returns:
\(\{\,\beta_{\mathrm p},\;\partial\beta_{\mathrm p}/\partial\alpha_{\mathrm p}\,\}\)
-
inline std::tuple<double, double> g(const Tensor33d &Fe_trial, const Tensor33d &df_dS, double dLambda, double betaP)
Scalar consistency function \(g(\Delta\lambda)\) and its derivative for the scalar-return mapping.
-
template<typename T>
inline Tensor33t<T> computeReturnMappingDirection(Tensor33t<T> Fe, T betaP) Flow direction for return mapping, \(\partial f/\partial \boldsymbol S\) (templated).
- Template Parameters:
T – Scalar type (double, hyper-dual, complex).
- Parameters:
Fe – Elastic deformation gradient \(\boldsymbol F^{\mathrm e}\).
betaP – Stress-like hardening \(\beta_{\mathrm p}\).
- Returns:
\(\partial f/\partial \boldsymbol S\).
-
template<typename T>
inline std::tuple<Tensor33d, Tensor33d> computePlasticIncrement(Tensor33t<T> df_dS, T dLambda) Incremental plastic flow via exponential map.
- Template Parameters:
T – Scalar type (double, hyper-dual, complex).
- Parameters:
df_dS – Flow direction in Mandel-stress space \(\partial f/\partial \boldsymbol S\).
dLambda – Plastic multiplier increment \(\Delta\lambda\).
- Returns:
\(\{\,\Delta\boldsymbol F^{\mathrm p},\;\partial\Delta\boldsymbol F^{\mathrm p}/\partial\boldsymbol S:\partial f/\partial\boldsymbol S\,\}\).
-
template<typename T>
inline VectorXt<T> computeResidualVector(const VectorXt<T> &X, const Tensor33d &FeTrial, const double alphaPTrial) Residual vector \(\boldsymbol R(\boldsymbol X)\) for the full return mapping (templated).
Vector \(\boldsymbol X\) has 11 unknowns: 9 for \(\boldsymbol F^{\mathrm e}\), 1 for \(\alpha_{\mathrm p}\), 1 for \(\Delta\lambda\). Residual components:
\(R_{0..8}\): elastic consistency (updated \(\boldsymbol F^{\mathrm e}\) vs trial),
\(R_{9}\): \(\alpha_{\mathrm p}\) update,
\(R_{10}\): yield \(f\).
- Template Parameters:
T – Scalar type (double, hyper-dual, complex).
- Parameters:
X – Current iterate \([\mathrm{vec}(\boldsymbol F^{\mathrm e}),\;\alpha_{\mathrm p},\;\Delta\lambda]^\mathrm T\).
FeTrial – Trial elastic deformation gradient \(\boldsymbol F^\text{e,trial}\).
alphaPTrial – Trial strain-like hardening variable \(\alpha_{\mathrm p}^\text{trial}\).
- Returns:
\(\boldsymbol R\).
-
inline std::tuple<Eigen::VectorXd, Eigen::MatrixXd> computeResidualVectorAndTangent(const Eigen::VectorXd &X, const Tensor33d &FeTrial, const double alphaPTrial)
Residual \(\boldsymbol R \) and Jacobian \(\partial\boldsymbol R/\partial\boldsymbol X\) for the full return mapping (double precision).
- Parameters:
X – Current iterate \([\mathrm{vec}(\boldsymbol F^{\mathrm e}),\;\alpha_{\mathrm p},\;\Delta\lambda]^\mathrm T\).
FeTrial – Trial elastic gradient \(\boldsymbol F^\text{e,trial}\).
alphaPTrial – Trial strain-like hardening variable \(\alpha_{\mathrm p}^\text{trial}\).
- Returns:
\(\{\,\boldsymbol R,\;\partial\boldsymbol R/\partial\boldsymbol X\,\}\).
Public Members
-
const double K
Bulk modulus (read from
materialProperties[0])
-
const double G
Shear modulus (read from
materialProperties[1])
-
const double fy
Initial yield stress (read from
materialProperties[2])
-
const double fyInf
Saturated (asymptotic) yield stress (read from
materialProperties[3])
-
const double eta
Saturation rate parameter (read from
materialProperties[4])
-
const double H
Linear hardening modulus (read from
materialProperties[5])
-
const int implementationType
Algorithm variant selector (read from
materialProperties[6]).
-
const double density
Density (read from
materialProperties[7]) (if provided).
-
std::unique_ptr<FiniteStrainJ2PlasticityStateVarManager> stateVars
-
class FiniteStrainJ2PlasticityStateVarManager : public MarmotStateVarVectorManager
Inheritance diagram for Marmot::Materials::FiniteStrainJ2Plasticity::FiniteStrainJ2PlasticityStateVarManager:
![digraph {
graph [bgcolor="#00000000"]
node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2]
edge [color="#1414CE"]
"1" [label="Marmot::Materials::FiniteStrainJ2Plasticity::FiniteStrainJ2PlasticityStateVarManager" tooltip="Marmot::Materials::FiniteStrainJ2Plasticity::FiniteStrainJ2PlasticityStateVarManager" fillcolor="#BFBFBF"]
"2" [label="MarmotStateVarVectorManager" tooltip="MarmotStateVarVectorManager"]
"1" -> "2" [dir=forward tooltip="public-inheritance"]
}](../../_images/graphviz-9d4bb1b3e325d1e2c2c1f42c8c39be6ca7c4d3e9.png)
Collaboration diagram for Marmot::Materials::FiniteStrainJ2Plasticity::FiniteStrainJ2PlasticityStateVarManager:
![digraph {
graph [bgcolor="#00000000"]
node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2]
edge [color="#1414CE"]
"1" [label="Marmot::Materials::FiniteStrainJ2Plasticity::FiniteStrainJ2PlasticityStateVarManager" tooltip="Marmot::Materials::FiniteStrainJ2Plasticity::FiniteStrainJ2PlasticityStateVarManager" fillcolor="#BFBFBF"]
"2" [label="MarmotStateVarVectorManager" tooltip="MarmotStateVarVectorManager"]
"3" [label="MarmotStateVarVectorManager::StateVarVectorLayout" tooltip="MarmotStateVarVectorManager::StateVarVectorLayout"]
"1" -> "2" [dir=forward tooltip="public-inheritance"]
"2" -> "3" [dir=forward tooltip="usage"]
}](../../_images/graphviz-38beec2f855fc3a0b3442153ff6f0656fdec87b2.png)
State variable manager for
FpandalphaP; provides named views and layout.Public Functions
-
inline FiniteStrainJ2PlasticityStateVarManager(double *theStateVarVector)
Bind the manager to an external contiguous buffer holding (
Fp,alphaP).
Public Members
-
double &alphaP
Strain-like hardening variable \(\alpha_{\mathrm p}\).
Public Static Attributes
-
static const auto layout =
makeLayout( {{ .name = "Fp", .length = 9 },{ .name = "alphaP", .length = 1 },} ) Memory layout of internal variables: 9 for
Fpand 1 foralphaP(total 10).
-
inline FiniteStrainJ2PlasticityStateVarManager(double *theStateVarVector)