StandardElasticity
brickTFEL/Material/MohrCoulombYieldCriterion.hxx
fileStandardElastoViscoPlasticity
brickThis page describes how to implement a non-associated plastic behaviour based on the Mohr-Coulomb criterion. The algorithm mostly follows the work of [1] and relies on an apex smoothing introduced by [2].
This implementation has been introduced in OpenGeoSys, as illustrated below:
This page shows:
StandardElastoViscoplasticity
brickThe behaviour is described by a standard decomposition of the strain \(\underline{\varepsilon}^{\mathrm{to}}\) in an elastic and a plastic component, respectively denoted \(\underline{\varepsilon}^{\mathrm{el}}\) and \(\underline{\varepsilon}^{\mathrm{p}}\):
\[ \underline{\varepsilon}^{\mathrm{to}}=\underline{\varepsilon}^{\mathrm{el}}+\underline{\varepsilon}^{\mathrm{p}} \]
The stress \(\underline{\sigma}\) is related to the the elastic strain \(\underline{\varepsilon}^{\mathrm{el}}\) by a the orthotropic elastic stiffness \(\underline{\underline{\mathbf{D}}}\):
\[ \underline{\sigma}= \underline{\underline{\mathbf{D}}}\,\colon\,\underline{\varepsilon}^{\mathrm{el}} \]
The plastic part of the behaviour is described by the following yield surface: \[ F = 0 \]
F is defined as follow:
\[ F = \frac{I_1}{3} \sin \phi + \sqrt{J_2 K(\theta)^2 + a^2 \sin^2 \phi} - c \cos \phi \qquad{(1)}\] where \[ K(\theta) = \begin{cases} \cos \theta - \frac{1}{\sqrt{3}} \sin \phi \sin \theta & |\theta| < \theta_\mathrm{T} \\ A + B \sin 3\theta + C \sin^2 3 \theta & |\theta| \geq \theta_\mathrm{T} \end{cases} \qquad{(2)}\]
\[ \begin{aligned} A &= - \frac{1}{\sqrt{3}} \sin \phi \mathrm{sign \theta} \sin \theta_\mathrm{T} - B \mathrm{sign \theta} \sin 3\,\theta_\mathrm{T} - C \sin^2 3\theta_\mathrm{T} + \cos \theta_\mathrm{T}\\ % B &= \frac{\mathrm{sign \theta} \sin 6\theta_\mathrm{T} \left( \cos \theta_\mathrm{T} - \frac{1}{\sqrt{3}} \sin \phi \mathrm{sign \theta} \sin \theta_\mathrm{T} \right) - 6 \cos 6\theta_\mathrm{T} \left( \mathrm{sign \theta} \sin \theta_\mathrm{T} + \frac{1}{\sqrt{3}} \sin \phi \cos \theta_\mathrm{T} \right) }{18 \cos^3 3 \theta_\mathrm{T}}\\ % C &= \frac{-\cos 3 \theta_\mathrm{T} \left( \cos \theta_\mathrm{T} - \frac{1}{\sqrt{3}} \sin \phi \mathrm{sign \theta} \sin \theta_\mathrm{T} \right) - 3 \mathrm{sign \theta} \sin 3 \theta_\mathrm{T} \left( \mathrm{sign \theta} \sin \theta_\mathrm{T} + \frac{1}{\sqrt{3}} \sin \phi \cos \theta_\mathrm{T} \right) }{18 \cos^3 3 \theta_\mathrm{T}} \end{aligned} \qquad{(3)}\]
and
\[ I_1 = {\mathrm{tr}{\left(\underline{\sigma}\right)}} \qquad J_2 = \frac{1}{2} \underline{\sigma}^\mathrm{D} \,:\,\underline{\sigma}^\mathrm{D} \qquad J_3 = \det \underline{\sigma}^\mathrm{D} \qquad \theta = \frac{1}{3} \arcsin \left( -\frac{3\sqrt{3} J_3}{2 \sqrt{J_2^3}} \right) \]
\(\underline{\sigma}^\mathrm{D}\), \(c\) and \(\phi\) are respectively the deviatoric part of the tensor \(\underline{\sigma}\), the cohesion and friction angle.
The contribution of the smoothing is visualized in Figure 1.
The plastic potential differs from the yield surface in order to more accurately estimate dilatancy, but has an analogous structure:
\[ G_\mathrm{F} = \frac{I_1}{3} \sin \psi + \sqrt{J_2 K_G^2 + a^2 \tan^2 \phi \cos^2 \psi} - c \cos \psi \qquad{(4)}\]
where \(\psi\) is the dilatancy angle. \(K_G\), \(A_G\) and \(B_G\) follow from (2) and (3) by substituting the friction angle with the dilatancy angle :
\[ K_G(\theta) = \begin{cases} \cos \theta - \frac{1}{\sqrt{3}} \sin \psi \sin \theta & |\theta| < \theta_\mathrm{T} \\ A_G + B_G \sin 3\theta + C_G \sin^2 3 \theta & |\theta| \geq \theta_\mathrm{T} \end{cases} \]
\[ \begin{aligned} A_G &= - \frac{1}{\sqrt{3}} \sin \psi \mathrm{sign \theta} \sin \theta_\mathrm{T} - B \mathrm{sign \theta} \sin 3\theta_\mathrm{T} - C \sin^2 3\theta_\mathrm{T} + \cos \theta_\mathrm{T}\\ % B_G &= \frac{\mathrm{sign \theta} \sin 6\theta_\mathrm{T} \left( \cos \theta_\mathrm{T} - \frac{1}{\sqrt{3}} \sin \psi \mathrm{sign \theta} \sin \theta_\mathrm{T} \right) - 6 \cos 6\theta_\mathrm{T} \left( \mathrm{sign \theta} \sin \theta_\mathrm{T} + \frac{1}{\sqrt{3}} \sin \psi \cos \theta_\mathrm{T} \right) }{18 \cos^3 3 \theta_\mathrm{T}}\\ % C_G &= \frac{-\cos 3 \theta_\mathrm{T} \left( \cos \theta_\mathrm{T} - \frac{1}{\sqrt{3}} \sin \psi \mathrm{sign \theta} \sin \theta_\mathrm{T} \right) - 3 \mathrm{sign \theta} \sin 3 \theta_\mathrm{T} \left( \mathrm{sign \theta} \sin \theta_\mathrm{T} + \frac{1}{\sqrt{3}} \sin \psi \cos \theta_\mathrm{T} \right) }{18 \cos^3 3 \theta_\mathrm{T}} \end{aligned} \]
Equation ((4)) can be written as follow:
\[ G_\mathrm{F} = \frac{I_1}{3} \sin \psi + \sqrt{J_2 K_G^2 + a_G^2 \sin^2 \psi} - c \cos \psi \qquad{(5)}\]
with
\[ a_G = a \dfrac{\tan \phi}{\tan \psi} \]
this formulation allows the use of the
StandardElastoViscoplasticity
brick (see the last
part).
Plastic flow follows in a general manner for a \(I_1,\,J_2,\,\theta\)-type yield surface as
\[ \underline{n} = \frac{\partial G_{F}}{\partial I_1} \underline{I} + \left( \frac{\partial G_{F}}{\partial J_2} + \frac{\partial G_{F}}{\partial \theta} \frac{\partial \theta}{\partial J_2} \right) \underline{\sigma}^{D} + \frac{\partial G_{F}}{\partial \theta} \frac{\partial \theta}{\partial J_3} J_3 {(\underline{\sigma}^{D})}^{-1} \,:\,\underline{\underline{\mathbf{P}}}^{D} \]
or by use of the Caley-Hamilton theorem as
\[ \underline{n} = \frac{\partial G_{F}}{\partial I_1} \underline{I} + \left( \frac{\partial G_{F}}{\partial J_2} + \frac{\partial G_{F}}{\partial \theta} \frac{\partial \theta}{\partial J_2} \right) \underline{\sigma}^{D} + \frac{\partial G_{F}}{\partial \theta} \frac{\partial \theta}{\partial J_3} (\underline{\sigma}^{D})^{2}\,:\,\underline{\underline{\mathbf{P}}}^{D} \]
The previous constitutive equations will be integrated using a standard implicit scheme.
Assuming a plastic loading, the system of equations to be solved is: \[ \left\{ \begin{aligned} \Delta\,\underline{\varepsilon}^{\mathrm{el}}-\Delta\,\underline{\varepsilon}^{\mathrm{to}}+\Delta\,p\,{\left.\underline{n}\right|_{t+\theta\,\Delta\,t}} &= 0 \\ F &= 0 \\ \end{aligned} \right. \]
where \({\left.X\right|_{t+\theta\,\Delta\,t}}\) is the value of \(X\) at \(t+\theta\,\Delta\,t\), \(\theta\) being a numerical parameter.
In the following, the first (tensorial) equation is noted \(f_{\underline{\varepsilon}^{\mathrm{el}}}\) and the second (scalar) equation is noted \(f_{p}\).
In practice, it is physically sound to make satisfy exactly the yield condition at the end of the time step (otherwise, stress extrapolation can lead to stress state outside the yield surface and spurious oscillations can also be observed). This leads to the choice \(\theta=1\).
The jacobian \(J\) of the implicit system can be decomposed by blocks: \[ J= \begin{pmatrix} {\displaystyle \frac{\displaystyle \partial f_{\underline{\varepsilon}^{\mathrm{el}}}}{\displaystyle \partial \Delta\,\underline{\varepsilon}^{\mathrm{el}}}} & {\displaystyle \frac{\displaystyle \partial f_{\underline{\varepsilon}^{\mathrm{el}}}}{\displaystyle \partial \Delta\,p}} & \\\\ {\displaystyle \frac{\displaystyle \partial f_{p}}{\displaystyle \partial \Delta\,\underline{\varepsilon}^{\mathrm{el}}}} & {\displaystyle \frac{\displaystyle \partial f_{p}}{\displaystyle \partial \Delta\,p}} \\ \end{pmatrix} \qquad{(6)}\]
The expression of the previous terms is given by:
\[ \left\{ \begin{aligned} {\displaystyle \frac{\displaystyle \partial f_{\underline{\varepsilon}^{\mathrm{el}}}}{\displaystyle \partial \Delta\,\underline{\varepsilon}^{\mathrm{el}}}} &= \underline{\underline{\mathbf{I}}} + \Delta \lambda \frac{\partial \underline{n}}{\partial \Delta \underline{\epsilon}_{el}} \\ {\displaystyle \frac{\displaystyle \partial f_{\underline{\varepsilon}^{\mathrm{el}}}}{\displaystyle \partial \Delta\,p}} &= \underline{n} \\ {\displaystyle \frac{\displaystyle \partial f_{p}}{\displaystyle \partial \Delta\,\underline{\varepsilon}^{\mathrm{el}}}} &= E^{-1} \underline{n}_F \,:\,\underline{\underline{\mathbf{C}}} \\ {\displaystyle \frac{\displaystyle \partial f_{p}}{\displaystyle \partial \Delta\,p}} &= 0 \end{aligned} \right. \]
Assuming an elastic loading, the system of equations to be solved is trivially: \[ \left\{ \begin{aligned} \Delta\,\underline{\varepsilon}^{\mathrm{el}}-\Delta\,\underline{\varepsilon}^{\mathrm{to}}&= 0 \\ \Delta\,p &= 0 \\ \end{aligned} \right. \]
The jacobian associated with this system is the identity matrix.
While not mandatory (the @DSL
keyword can be place
anywhere in the file), its is convenient to start the implementation by
declaring the domain specific language to be used. For an integration by
a \(\theta\)-scheme, the
Implicit
domain specific language is choosen:
@DSL Implicit;
The @Behaviour
keyword is used to give the name of the
behaviour.
@Behaviour MohrCoulombAbboSloan;
The following instructions give some information about the author, the date
@Author Thomas Nagel;
@Date 05/02/2019;
The @Algorithm
keyword is used to give the name of the
algorithm.
@Algorithm NewtonRaphson;
StandardElasticity
brickThe implicit scheme used satisfies the requirements of the
StandardElasticity
brick as described here.
The StandardElasticity
brick which provides:
StandardElasticity
is introduced as follows:@Brick StandardElasticity;
The following instruction changes the default value of the stopping criterion \(\epsilon\) used by the Newton-Raphson method and the time integration scheme parameter \(\theta\) :
@Theta 1.0;
@Epsilon 1.e-14;
Thanks to the StandardElasticity
brick, all the
modelling hypotheses can be supported. The following statement, starting
with the @ModellingHypotheses
, enables all the modelling
hypotheses:
@ModellingHypotheses {".+"};
The @RequireStiffnessTensor
keyword requires the
stiffness tensor to be computed by the calling code. This generally
means that some extra material properties will be introduced and handled
by the interface before the behaviour integration.
@RequireStiffnessTensor<UnAltered>;
The @StateVariable
keyword introduces the
EquivalentPlasticStrain \(\lambda\).
@StateVariable real lam;
.setGlossaryName("EquivalentPlasticStrain"); lam
The @MaterialProperty
keyword introduces several
properties, here:
@MaterialProperty stress c;
.setEntryName("Cohesion");
c@MaterialProperty real phi;
.setEntryName("FrictionAngle");
phi@MaterialProperty real psi;
.setEntryName("DilatancyAngle");
psi@MaterialProperty real lodeT;
.setEntryName("TransitionAngle");
lodeT@MaterialProperty stress a;
.setEntryName("TensionCutOffParameter"); a
In MFront
, an integration variable is defined to store a
variable and use it in various code block.
Here several local variables are declared such as the bolean variable F: if true, plastic loading
@LocalVariable Stensor np;
@LocalVariable bool F; // if true, plastic loading
@LocalVariable real sin_psi;
@LocalVariable real sin_phi;
@LocalVariable real cos_phi;
@LocalVariable real cos_lodeT;
@LocalVariable real sin_lodeT;
@LocalVariable real tan_lodeT;
@LocalVariable real cos_3_lodeT;
@LocalVariable real sin_3_lodeT;
@LocalVariable real cos_6_lodeT;
@LocalVariable real sin_6_lodeT;
@LocalVariable real tan_3_lodeT;
@LocalVariable real a_G;
The @InitLocalVariables
code block is called before the
behaviour integration.
@InitLocalVariables
{
First, we define some variables :
constexpr auto sqrt3 = Cste<real>::sqrt3;
constexpr auto isqrt3 = Cste<real>::isqrt3;
// conversion to rad
*= pi / 180.;
phi *= pi / 180.;
psi *= pi / 180.;
lodeT = sin(psi);
sin_psi = cos(phi);
cos_phi = sin(phi);
sin_phi = sin(lodeT);
sin_lodeT = cos(lodeT);
cos_lodeT = tan(lodeT);
tan_lodeT = cos(3. * lodeT);
cos_3_lodeT = sin(3. * lodeT);
sin_3_lodeT = cos(6. * lodeT);
cos_6_lodeT = sin(6. * lodeT);
sin_6_lodeT = tan(3. * lodeT);
tan_3_lodeT = (a * tan(phi)) / tan(psi) a_G
Then the computeElasticPrediction
method (introducted
with the StandardElasticity
brick) is used to compute \(\sigma^{el}\)
// elastic prediction
const auto sig_el = computeElasticPrediction();
The three invariant \(I_{1}^{el}\), \(J_{2}^{el}\) and \(J_{3}^{el}\) corresponding to the elastic prediction are calculated:
const auto s_el = deviator(sig_el);
const auto I1_el = trace(sig_el);
const auto J2_el = max((s_el | s_el) / 2., local_zero_tolerance);
const auto J3_el = det(s_el);
The \(\theta^{el}\) angle is defined:
const auto arg = min(max(-3. * sqrt3 * J3_el / (2. * J2_el * sqrt(J2_el)),
-1. + local_zero_tolerance),
1. - local_zero_tolerance);
const auto lode_el = 1. / 3. * asin(arg);
K is initiliazed as \(K^{el}\) value:
auto K = 0.0;
if (abs(lode_el) < lodeT) {
= cos(lode_el) - isqrt3 * sin_phi * sin(lode_el);
K } else {
const auto sign =
(max(lode_el / max(abs(lode_el), local_zero_tolerance), -1.), 1.);
minconst auto term1 = cos_lodeT - isqrt3 * sin_phi * sin_lodeT;
const auto term2 = sign * sin_lodeT + isqrt3 * sin_phi * cos_lodeT;
const auto term3 = 18. * cos_3_lodeT * cos_3_lodeT * cos_3_lodeT;
const auto B = ( sign * sin_6_lodeT * term1 -
6. * cos_6_lodeT * term2 ) / term3;
const auto C = (- cos_3_lodeT * term1 - 3. * sign * sin_3_lodeT * term2) / term3;
const auto A =
-isqrt3 * sin_phi * sign * sin_lodeT - B * sign * sin_3_lodeT -
* sin_3_lodeT*sin_3_lodeT + cos_lodeT;
C = A + B * arg + C * arg*arg;
K }
To finish \(F^{el}\) is calculated
and the normal np
is initialized.
const auto sMC =
/ 3 * sin_phi + sqrt(J2_el * K * K + a * a * sin_phi * sin_phi);
I1_el = sMC - c * cos_phi > 0.;
F = Stensor(real(0)); np
The implementation of the implicit system and its derivative is done
in the @Integrator
code block:
@Integrator{
Some expressions are defined
constexpr auto sqrt3 = Cste<real>::sqrt3;
constexpr auto isqrt3 = Cste<real>::isqrt3;
constexpr auto id = Stensor::Id();
constexpr auto id4 = Stensor4::Id();
If there is no plasticity (elastic strain) there is no need for additional calculation because the various variables have already been initialized with elastic hypothesis. If there is plastic strain the rest is necessary.
if (F) {
\(K\) and \(\dfrac{\partial K}{\partial \theta}\) are computed:
const auto s = deviator(sig);
const auto I1 = trace(sig);
const auto J2 = max((s | s) / 2., local_zero_tolerance);
const auto J3 = real(det(s) < 0. ? min(det(s), -local_zero_tolerance)
: max(det(s), local_zero_tolerance));
const auto arg = min(max(-3. * sqrt3 * J3 / (2. * J2 * sqrt(J2)),
-1. + local_zero_tolerance),
1. - local_zero_tolerance);
const auto lode = 1. / 3. * asin(arg);
const auto cos_lode = cos(lode);
const auto sin_lode = sin(lode);
const auto cos_3_lode = cos(3. * lode);
const auto sin_6_lode = sin(6. * lode);
const auto cos_6_lode = cos(6. * lode);
const auto sin_3_lode = arg;
const auto tan_3_lode = tan(3. * lode);
auto K = 0.;
auto dK_dlode = 1.;
if (abs(lode) < lodeT) {
= cos_lode - isqrt3 * sin_phi * sin_lode;
K = -sin_lode - isqrt3 * sin_phi * cos_lode;
dK_dlode } else {
const auto sign =
(max(lode / max(abs(lode), local_zero_tolerance), -1.), 1.);
minconst auto term1 = cos_lodeT - isqrt3 * sin_phi * sin_lodeT;
const auto term2 = sign * sin_lodeT + isqrt3 * sin_phi * cos_lodeT;
const auto term3 = 18. * cos_3_lodeT * cos_3_lodeT * cos_3_lodeT;
const auto B = ( sign * sin_6_lodeT * term1 -
6. * cos_6_lodeT * term2 ) / term3;
const auto C = (- cos_3_lodeT * term1 - 3. * sign * sin_3_lodeT * term2) / term3;
const auto A =
-isqrt3 * sin_phi * sign * sin_lodeT - B * sign * sin_3_lodeT -
* sin_3_lodeT*sin_3_lodeT + cos_lodeT;
C = A + B * sin_3_lode + C * sin_3_lode*sin_3_lode;
K = 3. * B * cos_3_lode + 3. * C * sin_6_lode;
dK_dlode }
\(K_G\), \(\dfrac{\partial K_G}{\partial \theta}\) and \(\dfrac{\partial^2 K_G}{\partial^2 \theta}\) are computed :
auto KG = 0.0; // move into a function to avoid code duplication
auto dKG_dlode = 1.;
auto dKG_ddlode = 1.;
if (abs(lode) < lodeT)
{
= cos_lode - isqrt3 * sin_psi * sin_lode;
KG = -sin_lode - isqrt3 * sin_psi * cos_lode;
dKG_dlode = -cos_lode + isqrt3 * sin_psi * sin_lode;
dKG_ddlode }
else
{
const auto sign =
(max(lode / max(abs(lode), local_zero_tolerance), -1.), 1.);
minconst auto term1 = cos_lodeT - isqrt3 * sin_psi * sin_lodeT;
const auto term2 = sign * sin_lodeT + isqrt3 * sin_psi * cos_lodeT;
const auto term3 = 18. * cos_3_lodeT * cos_3_lodeT * cos_3_lodeT;
const auto B = ( sign * sin_6_lodeT * term1 -
6. * cos_6_lodeT * term2 ) / term3;
const auto C = (- cos_3_lodeT * term1 - 3. * sign * sin_3_lodeT * term2) / term3;
const auto A =
-isqrt3 * sin_psi * sign * sin_lodeT - B * sign * sin_3_lodeT -
* sin_3_lodeT*sin_3_lodeT + cos_lodeT;
C = A + B * sin_3_lode + C * sin_3_lode * sin_3_lode;
KG = 3. * B * cos_3_lode + 3. * C * sin_6_lode;
dKG_dlode = -9. * B * sin_3_lode + 18. * C * cos_6_lode;
dKG_ddlode }
The flow direction is computed :
// flow direction
const auto dev_s_squared = computeJ3Derivative(
); // replaces dev_s_squared = deviator(square(s));
sigconst auto dG_dI1 = sin_psi / 3.;
const auto root = max(sqrt(J2 * KG * KG + a_G * a_G * sin_psi * sin_psi),
);
local_zero_toleranceconst auto dG_dJ2 = KG / (2. * root) * (KG - tan_3_lode * dKG_dlode);
const auto dG_dJ3 = J2 * KG * tan_3_lode / (3. * J3 * root) * dKG_dlode;
const auto n = eval(dG_dI1 * id + dG_dJ2 * s + dG_dJ3 * dev_s_squared);
The yield function F
is computed:
// yield function
const auto rootF = max(sqrt(J2 * K * K + a * a * sin_phi * sin_phi), local_zero_tolerance);
const auto Fy1 = I1 * sin_phi / 3 + rootF;
const auto Fy = Fy1 - c * cos_phi;
Derivatives are calculated before computing the Jacobian matrix:
\(\dfrac{\partial F}{\partial I_1}\), \(\dfrac{\partial F}{\partial J_2}\), \(\dfrac{\partial F}{\partial J_3}\), \(n_F\), \(\dfrac{\partial G}{\partial \theta}\), \(\dfrac{\partial^2 G}{\partial^2 \theta}\) \(\dfrac{\partial^2 G}{\partial \theta \partial J_2}\), \(\dfrac{\partial^2 G}{\partial^2 J_2}\), \(\dfrac{\partial^2 G}{\partial^2 J_3}\), \(\dfrac{\partial ^2G}{\partial J_2 \partial J_3}\)
// yield function derivative for Jacobian
const auto dF_dI1 = sin_phi / 3.;
const auto dF_dJ2 = K / (2. * rootF) * (K - tan_3_lode * dK_dlode);
const auto dF_dJ3 = J2 * K * tan_3_lode / (3. * J3 * rootF) * dK_dlode;
const auto nF = eval(dF_dI1 * id + dF_dJ2 * s + dF_dJ3 * dev_s_squared);
// building dfeel_ddeel
const auto Pdev = id4 - (id ^ id) / 3;
const auto dG_dlode = KG * J2 / (root)*dKG_dlode;
const auto dG_ddlode =
/ root * (dKG_dlode * dKG_dlode * (1. - J2 * KG * KG / (root * root)) + KG * dKG_ddlode);
J2 const auto dG_ddlodeJ2 = KG / root * dKG_dlode * (1. - J2 * KG * KG / (2 * root * root));
const auto dG_ddJ2 =
-KG * KG * KG * KG / (4. * root * root * root) + dG_dlode * tan_3_lode / (2 * J2 * J2) -
/ (2 * J2) * (2 * dG_ddlodeJ2 - tan_3_lode / (2 * J2) * dG_ddlode -
tan_3_lode 3 / (2 * J2 * cos_3_lode * cos_3_lode) * dG_dlode);
const auto dG_ddJ3 = -tan_3_lode / (3 * J3 * J3) * dG_dlode +
/ (3 * J3) * (dG_ddlode * tan_3_lode / (3 * J3) +
tan_3_lode * 1. / (J3 * cos_3_lode * cos_3_lode));
dG_dlode const auto dG_ddJ2J3 = dG_ddlodeJ2 * tan_3_lode / (3 * J3) -
/ (2 * J2) * (dG_ddlode * tan_3_lode / (3 * J3) +
tan_3_lode * 1. / (J3 * cos_3_lode * cos_3_lode)); dG_dlode
\(f^{el}\), \(\dfrac{\partial f^{el}}{\partial \Delta\,\underline{\varepsilon}^{\mathrm{el}}}\), \(\dfrac{\partial f^{el}}{\partial \Delta\,\lambda}\) are computed:
// elasticity
+= dlam * n;
feel += theta * dlam * (dG_dJ2 * Pdev + dG_dJ3 * computeJ3SecondDerivative(sig) +
dfeel_ddeel * (s ^ s) + dG_ddJ3 * (dev_s_squared ^ dev_s_squared) +
dG_ddJ2 * ((dev_s_squared ^ s) + (s ^ dev_s_squared))) *
dG_ddJ2J3 ;
D= n; dfeel_ddlam
These equations are equivalent to:
\[ \begin{aligned} & f^{el} = \Delta\,\underline{\varepsilon}^{\mathrm{el}}- \Delta\,\underline{\varepsilon}^{\mathrm{to}}\ + \Delta\,\lambda\,\underline{n} = 0 \\ & \dfrac{\partial f^{el}}{\partial \Delta\,\underline{\varepsilon}^{\mathrm{el}}} = 1 + \Delta \lambda \frac{\partial \underline{n}}{\partial \Delta \underline{\varepsilon}^{\mathrm{el}}} = 1 + \Delta \lambda \frac{\partial \underline{n}}{\partial \underline{\sigma}} \frac{\partial \underline{\sigma}}{\partial \Delta \underline{\varepsilon}^{\mathrm{el}}} \, \\ & \dfrac{\partial f^{el}}{\partial \Delta\,\lambda} = \underline{n} \end{aligned} \]
because this implementation takes into account the fact that the
Implicit
DSL automatically initializes feel
to
the current estimation of \(\Delta\,\underline{\varepsilon}^{\mathrm{el}}\)
and the Jacobian to identity. Moreover, the
StandardElasticity
brick automatically subtracts \(\Delta\,\underline{\varepsilon}^{\mathrm{to}}\)
to feel
.
And to finish \(f^{\lambda}\), \(\dfrac{\partial f^{\lambda}}{\partial \Delta\,\lambda}\), \(\dfrac{\partial f^{\lambda}}{\partial \Delta\,\underline{\varepsilon}^{\mathrm{el}}}\) are computed:
// plasticity
= Fy / D(0, 0);
flam = strain(0);
dflam_ddlam = theta * (nF | D) / D(0, 0);
dflam_ddeel = n; np
The Jacobian can now be computed (see (6))
The yield function was approached along different stress paths, shown in Figure 2 within the \(\pi\)-plane. This shows that a) yield is correctly detected, and b) the stress-state is correctly pulled back onto the yield surface.
The second test (included as a benchmark) follows an analytical solution for a stress-free cavity in an infinite medium under a variable far-field stress. The solution computes stress- and displacement fields as well as the location of the plastified zone. It has been used in [1] where the test is described in more detail (with a partial extension towards non-associated flow).
TFEL/Material/MohrCoulombYieldCriterion.hxx
fileBecause F et G have an analogous structure, it’s possible to simplify
the MFront file and to do the calculation in the hxx
TFEL/Material/MohrCoulombYieldCriterion.hxx
instead of the
MFront file.
Severals simplifications are done:
sin_phi
,
sin_psi
, … are now defined in the hxx filecomputeMohrCoulombStressCriterion
functioncomputeMohrCoulombStressCriterionNormal
functioncomputeMohrCoulombStressCriterionSecondDerivative
functionExcept for some name changes (for example p instead lam for the
EquivalentPlasticStrain) and the functions previously introduced
(computeMohrCoulombStressCriterion
,computeMohrCoulombStressCriterionNormal
and computeMohrCoulombStressCriterionSecondDerivative
) the
rest of the MFront file is identical to that described above in the
section Implementation
).
The new MFront file is however much shorter and clearer:
@DSL Implicit;
@Behaviour MohrCoulombAbboSloan;
@Author Thomas Nagel;
@Date 05 / 02 / 2019;
@Algorithm NewtonRaphson;
@Brick StandardElasticity;
@Includes{
#include "TFEL/Material/MohrCoulombYieldCriterion.hxx"
}
@Theta 1.0; // time integration scheme
@Epsilon 1e-14; // tolerance of local stress integration algorithm
@ModellingHypotheses{".+"};
@RequireStiffnessTensor<UnAltered>;
@Parameter pi = 3.14159265359;
@StateVariable real p;
.setGlossaryName("EquivalentPlasticStrain");
p
@MaterialProperty stress c;
.setEntryName("Cohesion");
c@MaterialProperty real phi;
.setEntryName("FrictionAngle");
phi@MaterialProperty real psi;
.setEntryName("DilatancyAngle");
psi@MaterialProperty real lodeT;
.setEntryName("TransitionAngle");
lodeT@MaterialProperty stress a;
.setEntryName("TensionCutOffParameter");
a
@LocalVariable Stensor ngp;
@LocalVariable bool F; // if true, plastic loading
@LocalVariable MohrCoulombParameters<StressStensor> pf;
@LocalVariable MohrCoulombParameters<StressStensor> pg;
@LocalVariable real a_G;
@InitLocalVariables {
= (a * tan((phi*pi)/180)) / tan((psi*pi)/180);
a_G constexpr const auto u = MohrCoulombParameters<StressStensor>::DEGREE;
= makeMohrCoulombParameters<StressStensor, u>(c, phi, lodeT, a);
pf = makeMohrCoulombParameters<StressStensor, u>(c, psi, lodeT, a_G);
pg const auto sel = computeElasticPrediction();
const auto smc = computeMohrCoulombStressCriterion(pf, sel);
= smc > stress(0);
F = Stensor(real(0));
ngp }
@Integrator {
if (F) {
// in C++11:
auto Fy = stress{};
auto nf = Stensor{};
auto Fg = stress{};
auto ng = Stensor{};
auto dng = Stensor4{};
std::tie(Fy, nf) = computeMohrCoulombStressCriterionNormal(pf, sig);
std::tie(Fg, ng, dng) = computeMohrCoulombStressCriterionSecondDerivative(pg, sig);
// elasticity
+= dp * ng;
feel += theta * dp * dng * D;
dfeel_ddeel = ng;
dfeel_ddp // plasticity
= Fy / D(0, 0);
fp = strain(0);
dfp_ddp = theta * (nf | D) / D(0, 0);
dfp_ddeel = ng;
ngp }
}
StandardElastoViscoPlasticity
brickThe Mohr-Coulomb model has been introduced in the
StandardElastoViscoPlasticity
brick.
The MFront
file is then very short: the whole model is
contained in the brick and can be called by the keyword
MohrCoulomb
The MFront file is now:
@DSL Implicit;
@Behaviour MohrCoulomAbboSloan3;
@Author HELFER Thomas 202608;
@Date 10 / 12 / 2019;
@Description {
}
@Algorithm NewtonRaphson;
@Epsilon 1.e-14;
@Theta 1;
@Brick StandardElastoViscoPlasticity{
: "Hooke" {
stress_potential : 150.e3,
young_modulus : 0.3
poisson_ratio },
: "Plastic" {
inelastic_flow : "MohrCoulomb" {
criterion : 3.e1, // cohesion
c : 0.523598775598299, // friction angle or dilatancy angle
phi : 0.506145483078356, // transition angle as defined by Abbo and Sloan
lodeT : 1e1 // tension cuff-off parameter
a },
: "MohrCoulomb" {
flow_criterion : 3.e1, // cohesion
c : 0.174532925199433, // friction angle or dilatancy angle
phi : 0.506145483078356, // transition angle as defined by Abbo and Sloan
lodeT : 3e1 // tension cuff-off parameter
a },
: "Linear" {R0 : "0"}
isotropic_hardening }
};