StandardElastoViscoPlasticity
brickThis page describes the StandardElastoViscoPlasticity
brick. This brick is used to describe a specific class of strain based
behaviours based on an additive split of the total strain \(\underline{\epsilon}^{\mathrm{to}}\) into
an elastic part \(\underline{\epsilon}^{\mathrm{el}}\) and
one or several inelastic strains describing plastic (time-independent)
flows and/or viscoplastic (time-dependent) flows: \[
\underline{\epsilon}^{\mathrm{to}}=\underline{\epsilon}^{\mathrm{el}}
+\sum_{i_{\mathrm{p}}=0}^{n_{\mathrm{p}}}\underline{\epsilon}^{\mathrm{p}}_{i_{\mathrm{p}}}
+\sum_{i_{\mathrm{vp}}=0}^{n_{\mathrm{vp}}}\underline{\epsilon}^{\mathrm{vp}}_{i_{\mathrm{vp}}}
\]
This equation defines the equation associated with the elastic strain \(\underline{\epsilon}^{\mathrm{el}}\).
The brick decomposes the behaviour into two components:
Porous viscoplasticity
This page only introduces stress criteria which are not coupled with the evolution of the porosity. Stress criteria coupled with the evolution of porosity is described on a dedicated page.
@Brick "StandardElastoViscoPlasticity" {
// Here the stress potential is given by the Hooke law. We define:
// - the elastic properties (Young modulus and Poisson ratio).
// Here the Young modulus is a function of the temperature.
// The Poisson ratio is constant.
// - the thermal expansion coefficient
// - the reference temperature for the thermal expansion
: "Hooke" {
stress_potential : "2.e5 - (1.e5*((T - 100.)/960.)**2)",
young_modulus : 0.3,
poisson_ratio : "1.e-5 + (1.e-5 * ((T - 100.)/960.) ** 4)",
thermal_expansion : 0
thermal_expansion_reference_temperature },
// Here we define only one viscplastic flow defined by the Norton law,
// which is based:
// - the von Mises stress criterion
// - one isotorpic hardening rule based on Voce formalism
// - one kinematic hardening rule following the Armstrong-Frederick law
: "Norton" {
inelastic_flow : "Mises",
criterion : "Voce" {R0 : 200, Rinf : 100, b : 20},
isotropic_hardening : "Armstrong-Frederick" {
kinematic_hardening : "1.e6 - 98500 * (T - 100) / 96",
C : "5000 - 5* (T - 100)"
D },
: "(4200. * (T + 20.) - 3. * (T + 20.0)**2)/4900.",
K : "7. - (T - 100.) / 160.",
n : 3
Ksf }
};
This stress potential implements the Hooke law, i.e. a linear relation between the elastic strain and the stress, as follows:
\[ \underline{\sigma}=\underline{\mathbf{D}}\,\colon\,\underline{\epsilon}^{\mathrm{el}} \] where \(\underline{\mathbf{D}}\) is the elastic stiffness tensor.
This stress potential applies to isotropic and orthotropic materials. This stress potential provides:
The Hooke stress potential is fully described here.
IsotropicDamage
stress potentialThis stress potential adds to the Hooke stress potential the description of an isotropioc damage. The relation \[ \underline{\sigma}=\left(1-d\right)\,\underline{\mathbf{D}}\,\colon\,\underline{\epsilon}^{\mathrm{el}} \] where \(\underline{\mathbf{D}}\) is the elastic stiffness tensor and \(d\) is the isotropic damage variable.
This stress potential inherits all the features and options provided by the Hooke stress potential. The Hooke stress potential is fully described here.
DDIF2
stress
potentialThe DDIF2
behaviour is used to describe the brittle
nature of nuclear fuel ceramics and is usually coupled with a
description of the viscoplasticity of those ceramics (See for example
[3]).
This stress potential adds to the Hooke stress potential the description of cracking through an additional strain. As such, it inherits all the features provided by the Hooke stress potential.
The Hooke stress potential is fully described here.
The DDIF2 stress potential is fully described here.
Plastic
inelastic flowThe plastic flow is defined by:
The plastic strain rate satisfies: \[ \underline{\dot{\epsilon}}^{\mathrm{p}}=\dot{\lambda}\,{{\displaystyle \frac{\displaystyle \partial g}{\displaystyle \partial \underline{\sigma}}}} \]
The plastic multiplier satifies the Kuhn-Tucker relation: \[ \left\{ \begin{aligned} \dot{\lambda}\,f{\left(\underline{\sigma},p\right)}&=0\\ \dot{\lambda}&\geq 0 \end{aligned} \right. \]
The flow is associated is \(f\) is equal to \(g\). In practice \(f\) is defined by a stress criterion \(\phi\), a set of kinematic hardening rules, and an isotropic hardening rule, as follows:
\[ f{\left(\underline{\sigma},p\right)}= \phi{\left(\underline{\sigma}-\sum_{i}\underline{X}_{i}\right)}-\sum_{i}R_{i}{\left(p\right)} \]
where \(p\) is the equivalent plastic strain. Here we have decomposed the limit of the elastic domain as a sum, denoted \(\sum_{i}R_{i}{\left(p\right)}\), to indicate that one may define it by combining various predefined forms of isotropic hardening rules (Voce, Linear, etc.) defined hereafter.
Plastic
flowDuring the Newton iterations, the current estimate of the equivalent stress \(\sigma_{\mathrm{eq}}\) may be significantly higher than the elastic limit \(R\). This may lead to a divergence of the Newton algorithm.
One may reject the Newton steps leading to such high values of the
elastic limit by specifying a relative threshold denoted \(\alpha\), i.e. if \(\sigma_{\mathrm{eq}}\) is greater than
\(\alpha\,\cdot\,R\). A typical value
for \(\alpha\) is \(1.5\). This relative threshold is specified
by the maximum_equivalent_stress_factor
option.
In some cases, rejecting steps may also lead to a divergence of the
Newton algorithm, so one may specify a relative threshold \(\beta\) on the iteration number which
desactivate this check, i.e. the check is performed only if the current
iteration number is below \(\beta\,\cdot\,i_{\max{}}\) where \(i_{\max{}}\) is the maximum number of
iterations allowed for the Newton algorithm. A typical value for \(\beta\) is \(0.4\). This relative threshold is specified
by the equivalent_stress_check_maximum_iteration_factor
option.
@DSL Implicit;
@Behaviour PerfectPlasticity;
@Author Thomas Helfer;
@Date 17 / 08 / 2020;
@Description{};
@Epsilon 1.e-14;
@Theta 1;
@Brick StandardElastoViscoPlasticity{
: "Hooke" {young_modulus : 200e9, poisson_ratio : 0.3},
stress_potential : "Plastic" {
inelastic_flow : "Mises",
criterion : "Linear" {R0 : 150e6},
isotropic_hardening : 1.5,
maximum_equivalent_stress_factor : 0.4
equivalent_stress_check_maximum_iteration_factor}
};
Norton
inelastic flowThe plastic flow is defined by:
HyperbolicSine
inelastic flowThe viscoplastic flow is defined by:
HarmonicSumOfNortonHoffViscoplasticFlows
inelastic
flowThe equivalent viscoplastic strain rate \(\dot{p}\) is defined as:
\[ \dfrac{1}{\dot{p}}=\sum_{i=1}^{N}\dfrac{1}{\dot{p}_{i}} \]
where \(\dot{p}_{i}\) has an expression similar to the the Norton-Hoff viscoplastic flow:
\[ \dot{p}_{i}=A_{i}\,{\left(\dfrac{\left<\phi{\left(\underline{\sigma}-\sum_{j}\underline{X}_{j}\right)}-\sum_{j}R_{j}{\left(p\right)}\right>}{K_{i}}\right)}^{n_{i}} \]
in which appear:
@Brick StandardElastoViscoPlasticity{
: "Hooke" {young_modulus : 150e9, poisson_ratio : 0.3},
stress_potential : "HarmonicSumOfNortonHoffViscoplasticFlows" {
inelastic_flow : "Mises",
criterion : {8e-67, 8e-67},
A : {1,1},
K : {8.2,8.2}
n }
};
The exponential nature of the hyperbolic sinus function may lead to divergence of the Newton method. To avoid this, one may specify a relative threshold denoted \(K_{sf}\): if the stress estimate is greater than \(K_{sf}\,K\), the step is rejected.
UserDefinedViscoplasticity
inelastic flowThe UserDefinedViscoplasticity
inelastic flow allows the
user to specify the viscoplastic strain rate vp
as a
function of f
and p
where:
f
is the positive part of the \(\phi{\left(\underline{\sigma}-\sum_{i}\underline{X}_{i}\right)}-\sum_{i}R_{i}{\left(p\right)}\)
where \(\phi\) is the stress
criterion.p
is the equivalent viscoplastic strain.This function shall be given by a string option named
vp
. This function must depend on f
. Dependance
to p
is optional.
The function may also depend on other variables. Let A
be such a variable. The UserDefinedViscoplasticity
flow
will look if an option named A
has been given to the
flow:
MFront
file.If required, the derivatives of vp
with respect to
f
and p
can be provided through the options
dvp_df
and dvp_dp
. The derivatives
dvp_df
and dvp_dp
can depend on two additional
variables, vp
and seps
, which denotes the
viscoplastic strain rate and a stress threshold.
If those derivatives are not provided, automatic differentiation will
be used. The user shall be warned that the automatic differentiation
provided by the tfel::math::Evaluator
class may result in
inefficient code.
@Parameter temperature Ta = 600;
@Parameter strain p0 = 1e-8;
@Brick StandardElastoViscoPlasticity{
: "Hooke" {young_modulus : 150e9, poisson_ratio : 0.3},
stress_potential : "UserDefinedViscoplasticity" {
inelastic_flow : "Mises",
criterion : 8.2,
E : "8e-67 * exp(- T / Ta)",
A : 0.32,
m : "A * (f ** E) / ((p + p0) ** m)",
vp : "E * vp / (max(f, seps))"
dvp_df // dvp_dp is evaluated by automatic differentiation (which is not recommended)
}
};
Some stress criteria (Hosford 1972, Barlat 2004, Mohr-Coulomb) shows sharp edges that may cause the failure of the standard Newton algorithm, due to oscillations in the prediction of the flow direction.
Rejecting Newton steps leading to a too large variation of the flow direction between the new estimate of the flow direction and the previous estimate is a cheap and efficient method to overcome this issue. This method can be viewed as a bisectional linesearch based on the Newton prediction: the Newton steps magnitude is divided by two if its results to a too large change in the flow direction.
More precisely, the change of the flow direction is estimated by the computation of the cosine of the angle between the two previous estimates:
\[ \cos{\left(\alpha_{\underline{n}}\right)}=\dfrac{\underline{n}\,\colon\,\underline{n}_{p}}{\lVert \underline{n}\rVert\,\lVert \underline{n}\rVert} \]
with \(\lVert \underline{n}\rVert=\sqrt{\underline{n}\,\colon\,\underline{n}}\).
The Newton step is rejected if the value of \(\cos{\left(\alpha_{\underline{n}}\right)}\) is lower than a user defined threshold. This threshold must be in the range \(\left[-1:1\right]\), but due to the slow variation of the cosine near \(0\), a typical value of this threshold is \(0.99\) which is equivalent to impose that the angle between two successive estimates is below \(8\mbox{}^{\circ}\).
The following section describes stress criteria available by default.
However, the StandardElastoViscoPlasticity
brick can also
be extended by the user:
The von Mises stress is defined by: \[ \sigma_{\mathrm{eq}}=\sqrt{\dfrac{3}{2}\,\underline{s}\,\colon\,\underline{s}}=\sqrt{3\,J_{2}} \] where: - \(\underline{s}\) is the deviatoric stress defined as follows: \[ \underline{s}=\underline{\sigma}-\dfrac{1}{3}\,{\mathrm{tr}{\left(\underline{\sigma}\right)}}\,\underline{I} \] - \(J_{2}\) is the second invariant of \(\underline{s}\).
In terms of the eigenvalues of the stress, denoted by \(\sigma_{1}\), \(\sigma_{2}\) and \(\sigma_{3}\), the von Mises stress can also be defined by: \[ \sigma_{\mathrm{eq}}=\sqrt{\dfrac{1}{2}{\left({\left|\sigma_{1}-\sigma_{2}\right|}^{2}+{\left|\sigma_{1}-\sigma_{3}\right|}^{2}+{\left|\sigma_{2}-\sigma_{3}\right|}^{2}\right)}} \]
This stress criterion does not have any option.
: "Mises" criterion
The Drucker 1949 stress is defined by: \[ \sigma_{\mathrm{eq}}=\sqrt{3}\sqrt[6]{J_{2}^3-c\,J_{3}^{2}} \] where:
: "Drucker 1949" {c : 1.285} criterion
The user must provide the \(c\) coefficient.
The Hosford equivalent stress is defined by (see [4]): \[ \sigma_{\mathrm{eq}}^{H}=\sqrt[a]{\dfrac{1}{2}{\left({\left|\sigma_{1}-\sigma_{2}\right|}^{a}+{\left|\sigma_{1}-\sigma_{3}\right|}^{a}+{\left|\sigma_{2}-\sigma_{3}\right|}^{a}\right)}} \] where \(\sigma_{1}\), \(\sigma_{2}\) and \(\sigma_{3}\) are the eigenvalues of the stress.
Therefore, when \(a\) goes to infinity, the Hosford stress reduces to the Tresca stress. When \(n = 2\) the Hosford stress reduces to the von Mises stress.
The Hosford exponent a
is mandatory.
Specifying the eigen solver using the eigen_solver
option is optional. This option can have the value default
or the value Jacobi
.
: "Hosford" {a : 6} criterion
The Hosford yield surface may have sharp edges which may lead to divergence of the Newton algorithm du to oscillations of the flow direction. Specifying a threshold for the angle between. See Section 2.2 for details.
The user must provide the Hosford exponent \(a\).
In order to describe yield differential effects, the isotropic Cazacu 2004 equivalent stress criterion is defined by (see [5]):
\[ \sigma_{\mathrm{eq}}=\sqrt[3]{J_{2}^{3/2} - c \, J_{3}} \]
where:
: "Isotropic Cazacu 2004" {c : -1.056} criterion
This Hill
criterion, also called Hill1948
criterion, is based on the equivalent stress \(\sigma_{\mathrm{eq}}^{H}\) defined as
follows: \[
\begin{aligned}
\sigma_{\mathrm{eq}}^{H}&=\sqrt{\underline{\sigma}\,\colon\,\underline{\mathbf{H}}\,\colon\,\underline{\sigma}}\\
&=\sqrt{F\,{\left(\sigma_{11}-\sigma_{22}\right)}^2+
G\,{\left(\sigma_{22}-\sigma_{33}\right)}^2+
H\,{\left(\sigma_{33}-\sigma_{11}\right)}^2+
2\,L\sigma_{12}^{2}+
2\,M\sigma_{13}^{2}+
2\,N\sigma_{23}^{2}}
\end{aligned}
\]
Warning This convention is given in the book of Lemaître et Chaboche and seems to differ from the one described in most other books.
This stress criterion has \(6\)
mandatory options: F
, G
, H
,
L
, M
, N
. Each of these options
must be interpreted as material property.
Orthotropic axis convention If an orthotropic axis convention is defined (See the
@OrthotropicBehaviour
keyword’ documentation), the coefficients of the Hill tensor can be exchanged for some modelling hypotheses. The coefficientsF
,G
,H
,L
,M
,N
must always correspond to the three dimensional case.
: "Hill" {F : 0.371, G : 0.629, H : 4.052, L : 1.5, M : 1.5, N : 1.5}, criterion
Within the framework of the theory of representation, generalizations to orthotropic conditions of the invariants of the deviatoric stress have been proposed by Cazacu and Barlat (see [6]):
Those invariants may be used to generalize isotropic yield criteria based on \(J_{2}\) and \(J_{3}\) invariants to orthotropy. The Cazacu 2001 equivalent stress criterion is defined as the orthotropic counterpart of the Drucker 1949 yield criterion, as follows (see [6]):
\[ \sigma_{\mathrm{eq}}=\sqrt{3}\sqrt[6]{\left(J_{2}^{O}\right)^3-c\,\left(J_{3}^{O}\right)^{2}} \]
This criterion requires the following options:
a
, as an array of \(6\) material properties.b
, as an array of \(11\) material properties.c
, as a material property.: "Cazacu 2001" {
criterion : {0.586, 1.05, 0.823, 0.96, 1, 1},
a : {1.44, 0.061, -1.302, -0.281, -0.375, 1, 1, 1, 1, 0.445, 1},
b : 1.285
c },
Proper support of orthotropic axes conventions has not been implemented yet for the computation of the \(J_{2}^{O}\) and \(J_{3}^{O}\). Thus, the following restrictions apply:
Tridimensional
modelling hypothesis is supported.qPlate
orthotropic axis convention is used, only
the Tridimensional
and PlaneStress
modelling
hypotheses are supported.Pipe
orthotropic axis convention is used, only
theI Tridimensional
, Axisymmetrical
,
AxisymmetricalGeneralisedPlainStrain
, and
AxisymmetricalGeneralisedPlainStres
modelling hypotheses
are supported.Using the invariants \(J_{2}^{O}\) and \(J_{3}^{O}\) previously defined, Cazacu and Barlat proposed the following criterion (See [5]):
\[ \sigma_{\mathrm{eq}}=\sqrt[3]{\left(J_{2}^{O}\right)^{3/2} - c\,J_{3}^{O}} \]
This criterion requires the following options:
a
, as an array of \(6\) material properties.b
, as an array of \(11\) material properties.c
, as a material property.: "Orthotropic Cazacu 2004" {
criterion : {0.586, 1.05, 0.823, 0.96, 1, 1},
a : {1.44, 0.061, -1.302, -0.281, -0.375, 1, 1, 1, 1, 0.445, 1},
b : 1.285
c },
Proper support of orthotropic axes conventions has not been implemented yet for the computation of the \(J_{2}^{O}\) and \(J_{3}^{O}\). Thus, the following restrictions apply:
Tridimensional
modelling hypothesis is supported.qPlate
orthotropic axis convention is used, only
the Tridimensional
and PlaneStress
modelling
hypotheses are supported.Pipe
orthotropic axis convention is used, only
the Tridimensional
, Axisymmetrical
,
AxisymmetricalGeneralisedPlainStrain
, and
AxisymmetricalGeneralisedPlainStres
modelling hypotheses
are supported.The Barlat equivalent stress is defined as follows (See [7]): \[ \sigma_{\mathrm{eq}}^{B}= \sqrt[a]{ \frac{1}{4}\left( \sum_{i=0}^{3} \sum_{j=0}^{3} {\left|s'_{i}-s''_{j}\right|}^{a} \right) } \]
where \(s'_{i}\) and \(s''_{i}\) are the eigenvalues of two transformed stresses \(\underline{s}'\) and \(\underline{s}''\) by two linear transformation \(\underline{\mathbf{L}}'\) and \(\underline{\mathbf{L}}''\): \[ \left\{ \begin{aligned} \underline{s}' &= \underline{\mathbf{L'}} \,\colon\,\underline{\sigma}\\ \underline{s}'' &= \underline{\mathbf{L''}}\,\colon\,\underline{\sigma}\\ \end{aligned} \right. \]
The linear transformations \(\underline{\mathbf{L}}'\) and \(\underline{\mathbf{L}}''\) are defined by \(9\) coefficients (each) which describe the material orthotropy. There are defined through auxiliary linear transformations \(\underline{\mathbf{C}}'\) and \(\underline{\mathbf{C}}''\) as follows: \[ \begin{aligned} \underline{\mathbf{L}}' &=\underline{\mathbf{C}}'\,\colon\,\underline{\mathbf{M}} \\ \underline{\mathbf{L}}''&=\underline{\mathbf{C}}''\,\colon\,\underline{\mathbf{M}} \end{aligned} \] where \(\underline{\mathbf{M}}\) is the transformation of the stress to its deviator: \[ \underline{\mathbf{M}}=\underline{\mathbf{I}}-\dfrac{1}{3}\underline{I}\,\otimes\,\underline{I} \]
The linear transformations \(\underline{\mathbf{C}}'\) and \(\underline{\mathbf{C}}''\) of the deviator stress are defined as follows: \[ \underline{\mathbf{C}}'= \begin{pmatrix} 0 & -c'_{12} & -c'_{13} & 0 & 0 & 0 \\ -c'_{21} & 0 & -c'_{23} & 0 & 0 & 0 \\ -c'_{31} & -c'_{32} & 0 & 0 & 0 & 0 \\ 0 & 0 & 0 & c'_{44} & 0 & 0 \\ 0 & 0 & 0 & 0 & c'_{55} & 0 \\ 0 & 0 & 0 & 0 & 0 & c'_{66} \\ \end{pmatrix} \quad \text{and} \quad \underline{\mathbf{C}}''= \begin{pmatrix} 0 & -c''_{12} & -c''_{13} & 0 & 0 & 0 \\ -c''_{21} & 0 & -c''_{23} & 0 & 0 & 0 \\ -c''_{31} & -c''_{32} & 0 & 0 & 0 & 0 \\ 0 & 0 & 0 & c''_{44} & 0 & 0 \\ 0 & 0 & 0 & 0 & c''_{55} & 0 \\ 0 & 0 & 0 & 0 & 0 & c''_{66} \\ \end{pmatrix} \]
When all the coefficients \(c'_{ji}\) and \(c''_{ji}\) are equal to \(1\), the Barlat equivalent stress reduces to the Hosford equivalent stress.
This stress criterion has \(3\) mandatory options:
l1
;l2
;Orthotropic axis convention If an orthotropic axis convention is defined (See the
@OrthotropicBehaviour
keyword’ documentation), the coefficients of the linear transformationscan be exchanged for some modelling hypotheses. The coefficients given by the user must always correspond to the three dimensional case.
Specifying the eigen solver using the eigen_solver
option is optional. This option can have the value default
or the value Jacobi
.
: "Barlat" {
criterion : 8,
a : {-0.069888, 0.079143, 0.936408, 0.524741, 1.00306, 1.36318, 0.954322,
l1 1.06906, 1.02377},
: {0.981171, 0.575316, 0.476741, 1.14501, 0.866827, -0.079294, 1.40462,
l2 1.1471, 1.05166}
}
The Barlat 2004 yield surface may have sharp edges which may lead to divergence of the Newton algorithm du to oscillations of the flow direction. Specifying a threshold for the angle between. See Section 2.2 for details.
Note
The follwing hardening rules can be combined to define more complex hardening rules. For example, the following code adds to Voce hardening:
: "Voce" {R0 : 600e6, Rinf : 900e6, b : 1}, isotropic_hardening : "Voce" {R0 : 0, Rinf : 300e6, b : 10}, isotropic_hardening
The previous code is equalivent to the following hardening rule:
\[ R{\left(p\right)}=R_{0}^{0}+{\left(R_{\infty}^{0}-R_{0}^{0}\right)}\,{\left(1-\exp{\left(-b^{0}\,p\right)}\right)}+R_{\infty}^{1}\,{\left(1-\exp{\left(-b^{1}\,p\right)}\right)} \]
with:
- \(R_{0}^{0}=600\,.\,10^{6}\,Pa\)
- \(R_{\infty}^{0}=900\,.\,10^{6}\,Pa\)
- \(R_{\infty}^{1}=300\,.\,10^{6}\,Pa\)
- \(b^{0}=1\)
- \(b^{1}=10\)
Linear
isotropic hardening ruleThe Linear
isotropic hardening rule is defined by: \[
R{\left(p\right)}=R_{0}+H\,p
\]
The Linear
isotropic hardening rule expects one of the
two following material properties:
R0
: the yield strengthH
: the hardening slopeNote
If one of the previous material property is not defined, the generated code is optimised and there will be no parameter asscoiated with it. To avoid this, you must define the material property and assign it to a zero value.
The following code can be added in a block defining an inelastic flow:
: "Linear" {R0 : 120e6, H : 438e6}, isotropic_hardening
Swift
isotropic hardening ruleThe Swift
isotropic hardening rule is defined by: \[
R{\left(p\right)}=R_{0}\,{\left(\dfrac{p+p_{0}}{p_{0}}\right)}^{n}
\]
The Swift
isotropic hardening rule expects three
material properties:
R0
: the yield strengthp0
n
The following code can be added in a block defining an inelastic flow:
: "Swift" {R0 : 120e6, p0 : 1e-8, n : 5.e-2} isotropic_hardening
Power
isotropic hardening rule (since TFEL 3.4)The Power
isotropic hardening rule is defined by: \[
R{\left(p\right)}=R_{0}\,{\left(p+p_{0}\right)}^{n}
\]
The Power
isotropic hardening rule expects at least the
following material properties:
R0
: the coefficient of the power lawn
: the exponent of the power lawThe p0
material property is optional and
generally is considered a numerical parameter to avoir an initial
infinite derivative of the power law when the exponent is lower than
\(1\).
The following code can be added in a block defining an inelastic flow:
: "Linear" {R0 : 50e6},
isotropic_hardening : "Power" {R0 : 120e6, p0 : 1e-8, n : 5.e-2} isotropic_hardening
Voce
isotropic hardening ruleThe Voce
isotropic hardening rule is defined by: \[
R{\left(p\right)}=R_{\infty}+{\left(R_{0}-R_{\infty}\right)}\,exp{\left(-b\,p\right)}
\]
The Voce
isotropic hardening rule expects three material
properties:
R0
: the yield strengthRinf
: the utimate strengthb
The following code can be added in a block defining an inelastic flow:
: "Voce" {R0 : 200, Rinf : 100, b : 20} isotropic_hardening
The UserDefined
isotropic hardening rule allows the user
to specify the radius of the yield surface as a function of the
equivalent plastic strain p
.
This function shall be given by a string option named R
and must depend on p
. The function may also depend on other
variables. Let A
be such a variable. The
UserDefined
isotropic hardening rule will look if an option
named A
has been given:
MFront
file.If required, the derivative of R
with respect to
p
can be provided through the option dR_dp
.
The derivative dR_dp
can depend on the variable
R
.
If this derivative is not provided, automatic differentiation will be
used. The user shall be warned that the automatic differentiation
provided by the tfel::math::Evaluator
class may result in
inefficient code.
@Parameter stress R0 = 200e6;
@Parameter stress Hy = 40e6;
@Parameter real b = 100;
@Brick StandardElastoViscoPlasticity{
: "Hooke" {young_modulus : 150e9, poisson_ratio : 0.3},
stress_potential : "Plastic" {
inelastic_flow : "Mises",
criterion : "UserDefined" {
isotropic_hardening : "R0 + Hy * (1 - exp(-b * p))", // Yield radius
R : "b * (R0 + Hy - R)"
dR_dp }
}
};
The Data
isotropic hardening rule allows the user to
define an isotropic hardening rule using a curve defined by a set of
pairs of equivalent strain and equivalent stress.
This isotropic hardening rule can be parametrised using three entries:
values
: which must a dictionnary giving the value of
the yield surface radius as a function of the equivalent plastic
strain.interpolation
: which allows to select the interpolation
type. Possible values are linear
(default choice) and
cubic_spline
.extrapolation
: which allows to select the extrapolation
type. Possible values are bound_to_last_value
(or
constant
) and extrapolation
(default
choice).@Brick StandardElastoViscoPlasticity{
: "Hooke" {young_modulus : 150e9, poisson_ratio : 0.3},
stress_potential : "Plastic" {
inelastic_flow : "Mises",
criterion : "Data" {
isotropic_hardening : {0 : 150e6, 1e-3 : 200e6, 2e-3 : 400e6},
values : "linear"
interpolation }
}
};
The StrainRateSensitive
isotropic hardening rule is
defined by: \[
R{\left(p, \dot{p}\right)} =
R_{0}{\left(p\right)}\,R_{rs}{\left(\dot{p}\right)},
\]
where:
This isotropic hardening rule can be parametrised using the following options:
rate_independent_isotropic_hardening
: this option
introduces a contribution to \(R_{0}{\left(p\right)}\). This option can be
repeated multiple times.rate_sensitivity_factor
: this option introduces the
rate sensivity factor \(R_{rs}{\left(\dot{p}\right)}\). The list of
available rate sensivity factors is given in section 2.5.@Brick StandardElastoViscoPlasticity{
: "Hooke" {young_modulus : 210e9, poisson_ratio : 0.3},
stress_potential : "Plastic" {
inelastic_flow : "Mises",
criterion : "StrainRateSensitive" {
isotropic_hardening :
rate_independent_isotropic_hardening "Swift" {R0 : "alpha * Ks * (p0s ** ns)", p0 : "p0s", n : "ns"},
: "Voce" {
rate_independent_isotropic_hardening : "(1 - alpha) * Q1 * (1 - Q2)",
R0 : "(1 - alpha) * Q1",
Rinf : "Q3"
b },
:
rate_sensitivity_factor "CowperSymonds" {dp0 : "dp0cs", n : "ncs", Rs_eps : 1e-8}
}
}
};
\[ R_{rs}{\left(\dot{p}\right)}=1+A\,\left(\frac{\dot{p}}{\dot{p}_{0}}\right)^{n} \]
When the exponent \(n\) is lower than one, the following regularised version, based on the user defined value \(R_{\epsilon}\) is available:
\[ R_{rs}{\left(\dot{p}\right)}=1+ \left\{ \begin{aligned} R_{\epsilon}\,\left(\frac{\dot{p}}{\dot{p}_{\epsilon}}\right) &&\text{if}&&\dot{p}<\dot{p}_{\epsilon}\\ A\,\left(\frac{\dot{p}}{\dot{p}_{0}}\right)^{n} &&\text{if}&&\dot{p}\geq\dot{p}_{\epsilon} \end{aligned} \right. \]
\[ R_{rs}{\left(\dot{p}\right)}=1+ \left\{ \begin{aligned} 0 &&\text{if}&&\dot{p}<\dot{p}_{0}\\ A\,\log\left(\frac{\dot{p}}{\dot{p}_{0}}\right)&&\text{if}&&\dot{p}\geq\dot{p}_{0} \end{aligned} \right. \]
Prager
kinematic hardening ruleThe following code can be added in a block defining an inelastic flow:
: "Prager" {C : 33e6}, kinematic_hardening
Armstrong-Frederick
kinematic hardening ruleThe Armstrong-Frederick
kinematic hardening rule can be
described as follows (see [8]): \[
\left\{
\begin{aligned}
\underline{X}&=\dfrac{2}{3}\,C\,\underline{a} \\
\underline{\dot{a}}&=\dot{p}\,\underline{n}-D\,\dot{p}\,\underline{a}
\\
\end{aligned}
\right.
\]
The following code can be added in a block defining an inelastic flow:
: "Armstrong-Frederick" {C : 1.5e9, D : 5} kinematic_hardening
Burlet-Cailletaud
kinematic hardening ruleThe Burlet-Cailletaud
kinematic hardening rule is
defined as follows (see [9]):
\[ \left\{ \begin{aligned} \underline{X}&=\dfrac{2}{3}\,C\,\underline{a} \\ \underline{\dot{a}}&=\dot{p}\,\underline{n} -\eta\,D\,\dot{p}\,\underline{a} -{\left(1-\eta\right)}\,D\,\dfrac{2}{3}\,\dot{p}\,{\left(\underline{a}\,\colon\,\underline{n}\right)}\,\underline{n} \\ \end{aligned} \right. \]
The following code can be added in a block defining an inelastic flow:
: "Burlet-Cailletaud" {C : 250e7, D : 100, eta : 0} kinematic_hardening
Chaboche 2012
kinematic hardening ruleThe Chaboche 2012
kinematic hardening rule is defined as
follows (see [10]):
\[ \underline{\dot{a}} =\underline{\dot{\varepsilon}}^{p}-\frac{3\,D}{2\,C}\,\Phi\left(p\right)\, \Psi^{\left(\underline{X}\right)}\left(\underline{X}\right)\,\dot{p}\,\underline{X} =\underline{\dot{\varepsilon}}^{p}- D\,\Phi\left(p\right)\,\Psi\left(\underline{a}\right)\dot{p}\,\underline{a} \]
with:
The following code can be added in a block defining an inelastic flow:
: "Chaboche 2012" {
kinematic_hardening : 250e7,
C : 100,
D : 2,
m : 0.6,
w }
The Delobelle-Robinet-Schaffler (DRS) kinematic hardening rule has been introduced to describe orthotropic viscoplasticity of Zircaloy alloys [11, 12]. It describes both dynamic and static recovery by the following evolution law: \[ \underline{\dot{a}}= \dot{p}\,\underline{\mathbf{E}}_{c}\,\colon\,\underline{n} -D\,\dot{p}\,\underline{\mathbf{R}}_{d}\,\colon\,\underline{a} -f\,{\left(\frac{a_{\mathrm{eq}}}{a_{0}}\right)}^{m}\,{{\displaystyle \frac{\displaystyle \partial a_{\mathrm{eq}}}{\displaystyle \partial \underline{a}}}} \] with \(a_{\mathrm{eq}}=\sqrt{\underline{a}\,\colon\,\underline{\mathbf{R}}_{s}\,\colon\,\underline{a}}\) and \({{\displaystyle \frac{\displaystyle \partial a_{\mathrm{eq}}}{\displaystyle \partial \underline{a}}}}=\frac{\underline{\mathbf{R}}_{s}\,\colon\,\underline{a}}{a_{\mathrm{eq}}}\)
The three fourth order tensors \(\underline{\mathbf{E}}_{c}\), \(\underline{\mathbf{R}}_{d}\) and \(\underline{\mathbf{R}}_{s}\) are assumed to have the same structure as the Hill tensors and are defined by \(6\) components (see this page for details).
The f
and a0
parameters are optional and
defaults to \(1\).
: "DRS" {
kinematic_hardening : 150.e9, // kinematic moduli
C : 1e2, // back-strain callback coefficient
D : 10,
f : 5,
m : {0.33, 0.33, 0.33, 1, 1, 1},
Ec : {0.33, 0.63, 0.33, 1, 1, 1},
Rs : {0.33, 0.33, 0.33, 1, 1, 1} //
Rd },