criteria.hpp
Overview
Yield and failure criteria functions including von Mises, Hill, Drucker, and anisotropic criteria. Provides equivalent stress computations and their derivatives for use in plasticity algorithms.
API Reference
-
double Drucker_stress(const arma::vec &v, const double &b, const double &n)
Provides the Drucker equivalent stress, given its vector representation.
Returns the Drucker equivalent stress \( \mathbf{\sigma}^{P} \), considering
\[ \sigma^{P} = \sigma^{VM} \left(\frac{1 + b \cdot J_3 \left(\mathbf{\sigma} \right)}{\left(J_2 \left(\mathbf{\sigma} \right) \right)^{3/2} } \right)^{m} \]considering the input stress \( \mathbf{\sigma} \), \( \mathbf{\sigma}^{VM} \) is the Von Mises computed equivalent stress, and \( b \) and \( m \) are parameter that define the equivalent stress. Note that if n > 10, the Drucker criteria is sufficiently close to the Mises that the Mises norm is used to avoid numerical instabilities from high-power computationsvec sigma = randu(6); double b = 1.2; double m = 0.5; double sigma_Drucker = Drucker_stress(sigma, b, n);
- Parameters:
v – The stress vector
b – Parameter that defines the equivalent stress
n – Parameter that defines the equivalent stress
- Returns:
The Drucker equivalent stress (double)
-
arma::vec dDrucker_stress(const arma::vec &v, const double &b, const double &n)
Provides derivative of the Drucker equivalent stress, given its vector representation.
Returns the derivative of the Drucker equivalent stress with respect to stress. It main use is to define evolution equations for strain based on an associated rule of a convex yield surface
\[ \frac{\partial \sigma^{P}}{\partial \mathbf{\sigma}} = \sqrt{3} \left( \frac{1. + b * J_3}{J_2^{3/2}} \right)^{1/n-1} \, \frac{1}{2} \frac{\sqrt{J_2}}{\mathbf{\sigma}'} + b \, m \left( 6 J_2^2 \right) \left( 6 J_2 \mathbf{\sigma}' \cdot \mathbf{\sigma}' - 4 J_2^2 \mathbf{I} + \frac{3}{m-9} J_3 \mathbf{\sigma}' \right) \]considering the input stress \( \mathbf{\sigma} \), \( \mathbf{\sigma}^{VM} \) is the Von Mises computed equivalent stress, and \( b \) and \( m \) are parameter that define the equivalent stress. Note that if n > 10, the Drucker criteria is sufficiently close to the Mises that the derivative of the Mises norm is used to avoid numerical instabilities from high-power computationsvec sigma = randu(6); double b = 1.2; double m = 0.5; vec dsigma_Druckerdsigma = dDrucker_stress(sigma, b, n);
- Parameters:
v – The stress vector
b – Parameter that defines the equivalent stress
n – Parameter that defines the equivalent stress
- Returns:
The derivative of the Drucker equivalent stress (arma::vec)
-
double Tresca_stress(const arma::vec &v)
Provides the Tresca equivalent stress, given its vector representation.
Returns the Tresca equivalent stress \( \mathbf{\sigma}^{T} \), considering
\[ \sigma^{T} = \sigma_{I} - \sigma_{III}, \]considering the input stress \( \mathbf{\sigma} \). Note that the principal stress are classified such that \( \sigma_{I} \geq \sigma_{II} \geq \sigma_{III} \).vec sigma = randu(6); double sigma_Tresca = Tresca_stress(sigma);
- Parameters:
v – The stress vector
- Returns:
The Tresca equivalent stress (double)
-
arma::vec dTresca_stress(const arma::vec &v)
Provides derivative of the Tresca equivalent stress, given its vector representation.
Returns the derivative of the Tresca equivalent stress with respect to stress. It main use is to define evolution equations for strain based on an associated rule of a Tresca convex but not smooth yield surface
Warning
Note that so far that the correct derivative it is not implemented! Only stress flow \( \eta_{stress}=\frac{3/2\sigma_{dev}}{\sigma_{Mises}} \) is returned
vec sigma = randu(6); double b = 1.2; double m = 0.5; vec dsigma_Trescadsigma = dTresca_stress(sigma, b, n);
- Parameters:
v – The stress vector
- Returns:
The derivative of the Tresca equivalent stress (arma::vec)
-
double Drucker_ani_stress(const arma::vec &v, const arma::vec ¶ms, const double &b, const double &n)
Provides an anisotropic Drucker-type equivalent stress combining Drucker criteria with DFA anisotropy.
Returns an anisotropic equivalent stress that combines the Drucker yield criterion with the Deshpande-Fleck-Ashby (DFA) anisotropic formulation:
\[ \sigma^{P}_{ani} = \sigma^{DFA} \left( \frac{1 + b \cdot J_3(\mathbf{\sigma})}{(J_2(\mathbf{\sigma}))^{3/2}} \right)^{1/n} \]where \( \sigma^{DFA} \) is the DFA anisotropic equivalent stress (see DFA_stress()). This formulation is particularly useful for modeling porous shape memory alloys (SMA). Note that if \( n > 10 \), the criterion reduces to the DFA stress to avoid numerical instabilities.vec sigma = randu(6); vec params = {0.5, 0.6, 0.7, 1.5, 1.5, 1.6, 1.2}; double b = 1.2; double n = 2.0; double sigma_ani = Drucker_anisotrope_stress(sigma, params, b, n);
- Parameters:
v – The stress vector in Voigt notation (6 components)
params – The DFA parameters vector (F, G, H, L, M, N, K) - see P_DFA() for details
b – Parameter controlling the influence of the third invariant \( J_3 \)
n – Exponent parameter that defines the equivalent stress shape
- Returns:
The anisotropic Drucker equivalent stress (double)
-
arma::vec dDrucker_ani_stress(const arma::vec &v, const arma::vec ¶ms, const double &b, const double &n)
Provides the derivative of the anisotropic Drucker-type equivalent stress.
Returns the derivative of the anisotropic Drucker-DFA equivalent stress with respect to stress. This derivative is computed using the chain rule combining the derivatives of the DFA stress, \( J_2 \), and \( J_3 \) invariants. It is primarily used to define evolution equations for strain based on an associated flow rule. Note that if \( n > 10 \), the derivative reduces to the standard stress flow direction \( \eta_{stress} \) to avoid numerical instabilities.
vec sigma = randu(6); vec params = {0.5, 0.6, 0.7, 1.5, 1.5, 1.6, 1.2}; double b = 1.2; double n = 2.0; vec dsigma_ani = dDrucker_anisotrope_stress(sigma, params, b, n);
- Parameters:
v – The stress vector in Voigt notation (6 components)
params – The DFA parameters vector (F, G, H, L, M, N, K) - see P_DFA() for details
b – Parameter controlling the influence of the third invariant \( J_3 \)
n – Exponent parameter that defines the equivalent stress shape
- Returns:
The derivative of the anisotropic Drucker equivalent stress (arma::vec)
-
arma::mat P_Ani(const arma::vec &P_params)
Returns an anisotropic configurational tensor \( P \) in the Voigt format (6x6 matrix)
The vector of parameters must be constituted of 9 values, respectively: \( P_{11}, P_{22}, P_{33}, P_{12}, P_{13}, P_{23}, P_{44}=P_{1212}, P_{55}=P_{1313}, P_{66}=P_{2323} \)
\[\begin{split} P_{ani} = \left( \begin{array}{cccccc} P_{11} & P_{12} & P_{13} & 0 & 0 & 0 \\ P_{12} & P_{22} & P_{23} & 0 & 0 & 0 \\ P_{13} & P_{23} & P_{33} & 0 & 0 & 0 \\ 0 & 0 & 0 & 2 \, P_{44} & 0 & 0 \\ 0 & 0 & 0 & 0 & 2 \, P_{55} & 0 \\ 0 & 0 & 0 & 0 & 0 & 2 \, P_{66} \end{array} \right) \end{split}\]The equivalent anisotropic stress is written as: \( \sigma^{eq}_{ani} = \sqrt{ \mathbf{\sigma} : P : \mathbf{\sigma} } \)It reduces to:
\[\begin{split} \begin{align} \sigma^{ani} & = \left( P_{11}\,\sigma_{11}^2 + P_{22}\,\sigma_{22}^2 + P_{33}\,\sigma_{33}^2 \right. \\ & + 2\,P_{12}\,\sigma_{11}\,\sigma_{22} + 2\,P_{13}\,\sigma_{11}\,\sigma_{33} + 2\,P_{23}\,\sigma_{22} \sigma_{33} \\ & + \left. 2\,P_{44}\,\sigma_{12}^2 + 2\,P_{55}\,\sigma_{13}^2 + 2\,P_{66}\,\sigma_{23}^2 \right)^{1/2} \end{align} \end{split}\]For the isotropic Mises equivalent stress, \( P \) reduces to:\[\begin{split} P_{Mises} = \left( \begin{array}{cccccc} 1 & -1/2 & -1/2 & 0 & 0 & 0 \\ -1/2 & 1 & -1/2 & 0 & 0 & 0 \\ -1/2 & -1/2 & 1 & 0 & 0 & 0 \\ 0 & 0 & 0 & 3 & 0 & 0 \\ 0 & 0 & 0 & 0 & 3 & 0 \\ 0 & 0 & 0 & 0 & 0 & 3 \end{array} \right) \end{split}\]vec P_params = {1., 1.2, 1.3, -0.2, -0.2, -0.33, 1., 1., 1.4}; mat P = P_Ani(P_params);
- Parameters:
P_params – The vector of parameters (9 components)
- Returns:
The anisotropic configurational tensor (arma::mat)
-
arma::mat P_Hill(const arma::vec &P_params)
Provides an anisotropic configurational tensor considering the quadratic Hill yield criterion (Hill, 1948) in the Voigt format (6x6 matrix), given its vector representation.
The vector of parameters must be constituted of 6 values, respectively: \( F,G,H,L,M,N \), such that
\[\begin{split} P_{Hill48} = \left( \begin{array}{cccccc} G + H & -H & -G & 0 & 0 & 0 \\ -H & F + H & -F & 0 & 0 & 0 \\ -G & -F & F + G & 0 & 0 & 0 \\ 0 & 0 & 0 & 2 \, N & 0 & 0 \\ 0 & 0 & 0 & 0 & 2 \, M & 0 \\ 0 & 0 & 0 & 0 & 0 & 2 \, L \end{array} \right) \end{split}\]considering the input stress \( \mathbf{\sigma} \). Note that the equivalent anisotropic Hill 1948) tensor is written as : \( \sigma^{eq}_{ani} = \sqrt{ \mathbf{\sigma} : P : \mathbf{\sigma} } \) which reduces to :\[\begin{split} \begin{align} \sigma^{H48} & = \left( H\, \left( \sigma_{11} - \sigma_{22} \right)^2 + G\, \left( \sigma_{11} - \sigma_{33} \right)^2 + F\, \left( \sigma_{22} - \sigma_{33} \right)^2 \right. \\ & + \left. 2\,N\,\sigma_{12}^2 + 2\,M\,\sigma_{13}^2 + 2\,L\,\sigma_{23}^2 \right)^{1/2} \end{align} \end{split}\]Considering the Mises equivalent strain\[\begin{split} \begin{align} \sigma^{M} & = \frac{1}{2} \left( \left( \sigma_{11} - \sigma_{22} \right)^2 + \left( \sigma_{11} - \sigma_{33} \right)^2 + \left( \sigma_{22} - \sigma_{33} \right)^2 \right. \\ & + \left. 6\,\sigma_{12}^2 + 6\,\sigma_{13}^2 + 6\,\sigma_{23}^2 \right)^{1/2} \end{align} \end{split}\]In that case, \( H \) reduces to:\[\begin{split} P_{Hill48} = \left( \begin{array}{cccccc} 1 & -1/2 & -1/2 & 0 & 0 & 0 \\ -1/2 & 1 & -1/2 & 0 & 0 & 0 \\ -1/2 & -1/2 & -1/2 & 0 & 0 & 0 \\ 0 & 0 & 0 & 3 & 0 & 0 \\ 0 & 0 & 0 & 0 & 3 & 0 \\ 0 & 0 & 0 & 0 & 0 & 3 \end{array} \right) \end{split}\]So that \( F = H = G = 1/2 \) and \( L = M = N = 3/2 \)vec P_params = {0.5,0.6,0.7,3.,3.,3.2}; mat P = P_Hill(P_params);
- Parameters:
P_params – The vector of parameters
- Returns:
The anisotropic Hill48 configurational tensor (arma::mat)
-
arma::mat P_DFA(const arma::vec &P_params)
Provides an anisotropic configurational tensor considering the quadratic Deshpande–Fleck–Ashby yield criterion (Deshpande, Fleck & Ashby, 2001) in the Voigt format (6x6 matrix), given its vector representation.
The vector of parameters must be constituted of 7 values, respectively: \( F,G,H,L,M,N,K \), such that
\[\begin{split} P_{DFA} = \left( \begin{array}{cccccc} G + H + K/9 & -H + K/9 & -G + K/9 & 0 & 0 & 0 \\ -H + K/9 & F + H + K/9 & -F + K/9 & 0 & 0 & 0 \\ -G + K/9 & -F + K/9 & F + G + K/9 & 0 & 0 & 0 \\ 0 & 0 & 0 & 2 \, N & 0 & 0 \\ 0 & 0 & 0 & 0 & 2 \, M & 0 \\ 0 & 0 & 0 & 0 & 0 & 2 \, L \end{array} \right) \end{split}\]considering the input stress \( \mathbf{\sigma} \). Note that the equivalent anisotropic tensor is written as : \( \sigma^{eq}_{ani} = \sqrt{ \mathbf{\sigma} : P : \mathbf{\sigma} } \) which reduces to :\[\begin{split} \begin{align} \sigma^{DFA} & = \left( H\, \left( \sigma_{11} - \sigma_{22} \right)^2 + G\, \left( \sigma_{11} - \sigma_{33} \right)^2 + F\, \left( \sigma_{22} - \sigma_{33} \right)^2 \right. \\ & + \left. 2\,N\,\sigma_{12}^2 + 2\,M\,\sigma_{13}^2 + 2\,L\,\sigma_{23}^2 + \frac{K}{9} \left( \sigma_{11} + \sigma_{22} + \sigma_{33} \right)^2 \right)^{1/2} \end{align} \end{split}\]with \( p = \frac{1}{3} \textrm{tr} (\sigma) \):\[\begin{split} \begin{align} \sigma^{DFA} & = \left( H\, \left( \sigma_{11} - \sigma_{22} \right)^2 + G\, \left( \sigma_{11} - \sigma_{33} \right)^2 + F\, \left( \sigma_{22} - \sigma_{33} \right)^2 \right. \\ & + \left. 2\,N\,\sigma_{12}^2 + 2\,M\,\sigma_{13}^2 + 2\,L\,\sigma_{23}^2 + K\,p^2 \right)^{1/2} \end{align} \end{split}\]Considering the full anisotric formulation:\[\begin{split} \begin{align} \sigma^{ani} & = \left( P_{11}\,\sigma_{11}^2 + P_{22}\,\sigma_{22}^2 + P_{33}\,\sigma_{33}^2 \right. \\ & + 2\,P_{12}\,\sigma_{11}\,\sigma_{22} + 2\,P_{13}\,\sigma_{11}\,\sigma_{33} + 2\,P_{23}\,\sigma_{22} \sigma_{33} \\ & + \left. 2\,P_{44}\,\sigma_{12}^2 + 2\,P_{55}\,\sigma_{13}^2 + 2\,P_{66}\,\sigma_{23}^2 \right)^{1/2} \end{align} \end{split}\]the above matrix is identified.vec P_params = {0.5,0.6,0.7,1.5,1.5,1.6,1.}; mat P = P_Hill(P_params);
- Parameters:
P_params – The vector of parameters
- Returns:
The anisotropic DFA configurational tensor (arma::mat)
-
double Eq_stress_P(const arma::vec &sigma, const arma::mat &H)
Provides the anisotropic equivalent stress, given the stress in a vector format and given a configurational tensor P.
Returns anisotropic equivalent stress \( \sigma^{eq}_{ani} = \sqrt{ \mathbf{\sigma} : P : \mathbf{\sigma} } \).
vec P_params = {0.5,0.6,0.7,1.5,1.5,1.6}; mat P = P_Hill(P_params); vec sigma = randu(6); double sigma_ani = Eq_stress_P(sigma,P_Hill);
- Parameters:
sigma – The stress vector
H – The configurational tensor
- Returns:
The anisotropic equivalent stress (double)
-
arma::vec dEq_stress_P(const arma::vec &sigma, const arma::mat &H)
Provides the derivative of the anisotropic equivalent stress, given the stress in a vector format and given a configurational tensor P.
Returns the derivative of the anisotropic equivalent stress \( \frac{\partial \mathbf{\sigma}^{ani}}{\partial \mathbf{\sigma}} \), considering
\[ \sigma^{eq}_{ani} = \frac{\mathbf{H} : \mathbf{\sigma}}{ \sqrt{ \mathbf{\sigma} : P : \mathbf{\sigma} } } \]vec P_params = {0.5,0.6,0.7,1.5,1.5,1.6}; mat P = P_Hill(P_params); vec sigma = randu(6); vec dsigma_ani = dEq_stress_P(sigma,P_Hill);
- Parameters:
sigma – The stress vector
H – The configurational tensor
- Returns:
The derivative of the anisotropic equivalent stress (mat::vec)
-
double Hill_stress(const arma::vec &sigma, const arma::vec &P_params)
Provides the Hill48 anisotropic equivalent stress, given the stress in a vector format and a vector of parameters (F,G,H,L,M,N)
Returns the Hill48 anisotropic equivalent stress \( \sigma^{eq}_{ani} = \sqrt{ \mathbf{\sigma} : P : \mathbf{\sigma} } \). see the function P_Hill() for more details to obtain the tensor H from the set of parameters (F,G,H,L,M,N).
vec P_params = {0.5,0.6,0.7,1.5,1.5,1.6}; vec sigma = randu(6); double sigma_Hill = Hill_stress(sigma, P_params);
- Parameters:
sigma – The stress vector
P_params – The vector of parameters
- Returns:
The Hill48 anisotropic equivalent stress (double)
-
arma::vec dHill_stress(const arma::vec &sigma, const arma::vec &P_params)
Provides the derivative of the Hill48 anisotropic equivalent stress, given the stress in a vector format and a vector of parameters (F,G,H,L,M,N)
Returns the derivative of the anisotropic equivalent stress \( \frac{\partial \mathbf{\sigma}^{ani}}{\partial \mathbf{\sigma}} \), considering
\[ \sigma^{eq}_{ani} = \frac{\mathbf{H} : \mathbf{\sigma}}{ \sqrt{ \mathbf{\sigma} : P : \mathbf{\sigma} } } \]see the function P_Hill() for more details to obtain the tensor H from the set of parameters (F,G,H,L,M,N).vec P_params = {0.5,0.6,0.7,1.5,1.5,1.6}; vec sigma = randu(6); vec dsigma_Hill = dHill_stress(sigma, P_params);
- Parameters:
sigma – The stress vector
P_params – The vector of parameters
- Returns:
The derivative of the Hill48 anisotropic equivalent stress (arma::vec)
-
double DFA_stress(const arma::vec &sigma, const arma::vec &P_params)
Provides the Deshpande–Fleck–Ashby (2001) anisotropic equivalent stress, given the stress in a vector format and a vector of parameters (F,G,H,L,M,N,K)
Returns the DFA anisotropic equivalent stress \( \sigma^{eq}_{ani} = \sqrt{ \mathbf{\sigma} : P : \mathbf{\sigma} } \). see the function P_DFA() for more details to obtain the tensor H from the set of parameters (F,G,H,L,M,N,K).
vec P_params = {0.5,0.6,0.7,1.5,1.5,1.6,1.2}; vec sigma = randu(6); double sigma_DFA = DFA_stress(sigma, P_params);
- Parameters:
sigma – The stress vector
P_params – The vector of parameters
- Returns:
The DFA anisotropic equivalent stress (double)
-
arma::vec dDFA_stress(const arma::vec &sigma, const arma::vec &P_params)
Provides the derivative of the Deshpande–Fleck–Ashby (2001) anisotropic equivalent stress, given the stress in a vector format and a vector of parameters (F,G,H,L,M,N,K)
Returns the derivative of the DFA anisotropic equivalent stress \( \frac{\partial \mathbf{\sigma}^{ani}}{\partial \mathbf{\sigma}} \), considering
\[ \sigma^{eq}_{ani} = \frac{\mathbf{H} : \mathbf{\sigma}}{ \sqrt{ \mathbf{\sigma} : P : \mathbf{\sigma} } } \]see the function P_DFA() for more details to obtain the tensor H from the set of parameters (F,G,H,L,M,N,K).vec P_params = {0.5,0.6,0.7,1.5,1.5,1.6,1.2}; vec sigma = randu(6); vec dsigma_DFA = dDFA_stress(sigma, P_params);
- Parameters:
sigma – The stress vector
P_params – The vector of parameters
- Returns:
The derivative of the DFA anisotropic equivalent stress (arma::vec)
-
double Ani_stress(const arma::vec &sigma, const arma::vec &P_params)
Provides the anisotropic equivalent stress, given the stress in a vector format and a vector of parameters \( ( P_{11},P_{22},P_{33},P_{12},P_{13},P_{23},P_{44},P_{55},P_{66} ) \).
Returns the anisotropic equivalent stress \( \sigma^{eq}_{ani} = \sqrt{ \mathbf{\sigma} : P : \mathbf{\sigma} } \). see the function P_Ani() for more details to obtain the tensor P from the set of parameters.
vec P_params = {1.,1.2,1.3,-0.2,-0.2,-0.33,1.,1.,1.4}; vec sigma = randu(6); double sigma_Ani = Ani_stress(sigma, P_params);
- Parameters:
sigma – The stress vector
P_params – The vector of parameters
- Returns:
The anisotropic equivalent stress (double)
-
arma::vec dAni_stress(const arma::vec &sigma, const arma::vec &P_params)
Provides the derivative of the Hill48 anisotropic equivalent stress, given the stress in a vector format and a vector of parameters (F,G,H,L,M,N)
Returns the derivative of the anisotropic equivalent stress \( \frac{\partial \mathbf{\sigma}^{ani}}{\partial \mathbf{\sigma}} \), considering
\[ \sigma^{eq}_{ani} = \frac{\mathbf{H} : \mathbf{\sigma}}{ \sqrt{ \mathbf{\sigma} : P : \mathbf{\sigma} } } \]see the function P_Ani() for more details to obtain the tensor P from the set of parameters.vec P_params = {1.,1.2,1.3,-0.2,-0.2,-0.33,1.,1.,1.4}; vec sigma = randu(6); vec dsigma_Ani = dAni_stress(sigma, P_params);
- Parameters:
sigma – The stress vector
P_params – The vector of parameters
- Returns:
The derivative of the Hill48 anisotropic equivalent stress (arma::vec)
-
double Eq_stress(const arma::vec &sigma, const std::string &type_eq, const arma::vec &P_params = arma::zeros(1))
Provides an equivalent stress, given the stress in a vector format, a string indicating the type of equivalent stress and a vector of parameters.
Returns the anisotropic equivalent stress, depending on the equivalen type provided : “Mises”, “Tresca”, “Drucker”, “Hill”, “Ani” See the detailed function for each equivalent stress to determine the number and type of parameters to provide, in the form of an arma::vec (so for the Drucker criteria, P_params = {b,n} for instance)
vec P_params = {1.,1.2,1.3,-0.2,-0.2,-0.33,1.,1.,1.4}; vec sigma = randu(6); double sigma_Ani = Ani_stress(sigma, P_params);
- Parameters:
sigma – The stress vector
type_eq – The type of equivalent stress (“Mises”, “Tresca”, “Drucker”, “Hill”, “Ani”)
P_params – The vector of parameters
- Returns:
The equivalent stress (double)
-
arma::vec dEq_stress(const arma::vec &sigma, const std::string &type_eq, const arma::vec &P_params = arma::zeros(1))
Provides an derivative of the equivalent stress, given the stress in a vector format, a string indicating the type of equivalent stress and a vector of parameters.
Returns the anisotropic equivalent stress, depending on the equivalen type provided : “Mises”, “Tresca”, “Drucker”, “Hill”, “Ani” See the detailed function for each equivalent stress to determine the number and type of parameters to provide, in the form of an arma::vec (so for the Drucker criteria, P_params = {b,n} for instance)
vec P_params = {1.,1.2,1.3,-0.2,-0.2,-0.33,1.,1.,1.4}; vec sigma = randu(6); double sigma_Ani = Ani_stress(sigma, P_params);
- Parameters:
sigma – The stress vector
type_eq – The type of equivalent stress (“Mises”, “Tresca”, “Drucker”, “Hill”, “Ani”)
P_params – The vector of parameters
- Returns:
The deriative of the equivalent stress (arma::vec)