kinematics.hpp
Overview
Functions for finite strain kinematics, including deformation gradients, stretch tensors (left and right), rotation tensors, and various strain measures (Green-Lagrange, Euler-Almansi, logarithmic).
API Reference
-
arma::mat ER_to_F(const arma::mat &E, const arma::mat &R)
Provides the transformation gradient, from the Green-Lagrange strain and the rotation.
The transformation gradient \(\mathbf{F}\) is related to the Green-Lagrange strain \(\mathbf{E}\) and the rotation \(\mathbf{R}\) by the following equations:
\[ \mathbf{F} = \mathbf{R} \cdot \mathbf{U} \quad \mathbf{E} = \frac{1}{2} \left( \sqrt{\mathbf{U}^2} - \mathbf{I} \right) \]where \(\mathbf{U}\) is the right stretch tensor 3x3 matrix and I is the 3x3 identity matrix.Example:
mat E = randu(3,3); mat R = eye(3,3); mat F = ER_to_F(E, R);
- Parameters:
E – 3x3 matrix representing the Green-Lagrange strain tensor
R – 3x3 matrix representing the rotation tensor
- Returns:
a 3x3 matrix representing the transformation gradient \(\mathbf{F}\)
-
arma::mat eR_to_F(const arma::mat &e, const arma::mat &R)
Provides the transformation gradient, from the logarithmic strain and the rotation.
The transformation gradient \(\mathbf{F}\) is related to the logarithmic strain \(\mathbf{e}\) and the rotation \(\mathbf{R}\) by the following equations:
\[ \mathbf{F} = \mathbf{V} \cdot \mathbf{R} \quad \mathbf{e} = \textrm{ln} \mathbf{V} \]where \(\mathbf{V}\) the right stretch tensor 3x3 matrix
Example:
mat e = randu(3,3); mat R = eye(3,3); mat F = eR_to_F(e, R);
- Parameters:
e – 3x3 matrix representing the logarithmic strain tensor
R – 3x3 matrix representing the rotation tensor
- Returns:
a 3x3 matrix representing the transformation gradient
-
arma::mat G_UdX(const arma::mat &F)
Provides the gradient of the displacement (Lagrangian) from the transformation gradient.
The gradient of the displacement (Lagrangian) is related to the transformation gradient \(\mathbf{F}\) by the following equation:
\[ \nabla_X \mathbf{U} = \mathbf{F} - \mathbf{I} \]where \(\mathbf{I}\) is the 3x3 identity matrix.Example:
mat F = randu(3,3); mat GradU = G_UdX(F);
- Parameters:
F – 3x3 matrix representing the transformation gradient
- Returns:
a 3x3 matrix representing the gradient of the displacement (Lagrangian)
-
arma::mat G_Udx(const arma::mat &F)
Provides the gradient of the displacement (Eulerian) from the transformation gradient.
The gradient of the displacement (Eulerian) is related to the transformation gradient \(\mathbf{F}\) by the following equation:
\[ \nabla_x \mathbf{U} = \mathbf{I} - \mathbf{F}^{-1} \]where \(\mathbf{I}\) is the 3x3 identity matrix.Example:
mat F = randu(3,3); mat gradU = G_Udx(F);
- Parameters:
F – 3x3 matrix representing the transformation gradient
- Returns:
3x3 matrix representing the gradient of the displacement (Eulerian)
-
arma::mat R_Cauchy_Green(const arma::mat &F)
Provides the Right Cauchy-Green tensor from the transformation gradient.
The Right Cauchy-Green tensor \(\mathbf{C}\) is related to the transformation gradient \(\mathbf{F}\) by the following equation:
\[ \mathbf{C} = \mathbf{F}^T \cdot \mathbf{F} \]Example:
mat F = randu(3,3); mat C = R_Cauchy_Green(F);
- Parameters:
F – 3x3 matrix representing the transformation gradient
- Returns:
3x3 matrix representing the Right Cauchy-Green tensor
-
arma::mat L_Cauchy_Green(const arma::mat &F)
Provides the Left Cauchy-Green tensor from the transformation gradient.
The Left Cauchy-Green tensor \(\mathbf{B}\) is related to the transformation gradient \(\mathbf{F}\) by the following equation:
\[ \mathbf{B} = \mathbf{F} \cdot \mathbf{F}^T \]Example:
mat F = randu(3,3); mat B = L_Cauchy_Green(F);
- Parameters:
F – 3x3 matrix representing the transformation gradient
- Returns:
3x3 matrix representing the Left Cauchy-Green tensor
-
void RU_decomposition(arma::mat &R, arma::mat &U, const arma::mat &F)
Provides the RU decomposition of the transformation gradient.
The RU decomposition of the transformation gradient \(\mathbf{F}\) is related to the matrices \(\mathbf{R}\) and \(\mathbf{U}\) by the following equations:
\[ \mathbf{F} = \mathbf{R} \cdot \mathbf{U} \quad \mathbf{U} = \sqrt{\mathbf{F}^T \cdot \mathbf{F}} \quad \mathbf{R} = \mathbf{F} \cdot \mathbf{U}^{-1} \]Example:
mat F = randu(3,3); mat R = zeros(3,3); mat U = zeros(3,3); RU_decomposition(R, U, F);
- Parameters:
R – 3x3 matrix that is the rotation matrix \(\mathbf{R}\)
U – 3x3 matrix that is the right stretch tensor \(\mathbf{U}\)
F – 3x3 matrix representing the transformation gradient
-
void VR_decomposition(arma::mat &V, arma::mat &R, const arma::mat &F)
Provides the VR decomposition of the transformation gradient.
The VR decomposition of the transformation gradient \(\mathbf{F}\) is related to the matrices \(\mathbf{R}\) and \(\mathbf{V}\) by the following equations:
\[ \mathbf{F} = \mathbf{V} \cdot \mathbf{R} \quad \mathbf{V} = \sqrt{\mathbf{F} \cdot \mathbf{F}^T} \quad \mathbf{R} = \mathbf{V}^{-1} \cdot \mathbf{F} \]Example:
mat F = randu(3,3); mat V = zeros(3,3); mat R = zeros(3,3); VR_decomposition(V, R, F);
- Parameters:
R – 3x3 matrix that is the rotation matrix \(\mathbf{R}\)
V – 3x3 matrix that is the left stretch tensor \(\mathbf{V}\)
F – 3x3 matrix representing the transformation gradient
-
arma::vec Inv_X(const arma::mat &X)
Provides the invariants of a symmetric tensor.
The invariants of the symmetric tensor \(\mathbf{X}\) are defined as (note that this definition is not unique):
\[ \mathbf{I}_1 = \textrm{trace} \left( X \right) \quad \mathbf{I}_2 = \frac{1}{2} \left( \textrm{trace} \left( X \right)^2 - \textrm{trace} \left( X^2 \right) \right) \quad \mathbf{I}_3 = \textrm{det} \left( X \right) \]Example:
mat F = randu(3,3); mat C = R_Cauchy_Green(F); vec I = Inv_X(F);
- Parameters:
X – 3x3 matrix representing the symmetric tensor
- Returns:
Vector containing the three invariants of the tensor
-
arma::mat Cauchy(const arma::mat &F)
Provides the Cauchy tensor from the transformation gradient.
The Cauchy tensor \(\mathbf{b}\) is related to the transformation gradient \(\mathbf{F}\) by the following equation:
\[ \mathbf{b} = \left( \mathbf{F} \cdot \mathbf{F}^T \right)^{-1} \]Note that this tensor \(\mathbf{b}\) is the inverse of the Left Cauchy-Green tensor \(\mathbf{B}\).
Example:
mat F = randu(3,3); mat b = Cauchy(F);
- Parameters:
F – 3x3 matrix representing the transformation gradient
- Returns:
3x3 matrix representing the Cauchy tensor
-
arma::mat Green_Lagrange(const arma::mat &F)
Provides the Green-Lagrange tensor from the transformation gradient.
The Green-Lagrange tensor \(\mathbf{E}\) is related to the transformation gradient \(\mathbf{F}\) by the following equation:
\[ \mathbf{E} = \frac{1}{2} \left( \mathbf{F}^T \cdot \mathbf{F} - \mathbf{I} \right) \]Example:
mat F = randu(3,3); mat E = Green_Lagrange(F);
- Parameters:
F – 3x3 matrix representing the transformation gradient \(\mathbf{F}\)
- Returns:
3x3 matrix representing the Green-Lagrange tensor \(\mathbf{E}\)
-
arma::mat Euler_Almansi(const arma::mat &F)
Provides the Euler-Almansi tensor from the transformation gradient.
The Euler-Almansi tensor \(\mathbf{A}\) is related to the transformation gradient \(\mathbf{F}\) by the following equation:
\[ \mathbf{A} = \frac{1}{2} \left( \mathbf{I} - \left( \mathbf{F} \cdot \mathbf{F}^T \right)^T \right) \]Example:
mat F = randu(3,3); mat A = Euler_Almansi(F);
- Parameters:
F – 3x3 matrix representing the transformation gradient \(\mathbf{F}\)
- Returns:
3x3 matrix representing the Euler-Almansi tensor \(\mathbf{A}\)
-
arma::mat Log_strain(const arma::mat &F)
Provides the logarithmic strain tensor from the transformation gradient.
The logarithmic strain tensor \(\mathbf{e}\) is related to the transformation gradient \(\mathbf{F}\) by the following equation:
\[ \mathbf{e} = \textrm{ln} \left( V \right) = \frac{1}{2} \textrm{ln} \left( V^2 \right) = \frac{1}{2} \textrm{ln} \left( \mathbf{F} \cdot \mathbf{F}^T \right) \]Note that the log function is here the matrix logarithm of symmetric/hermitian positive definite matrix V^2. TSee the armadillo documentation for more details of the log matrix implementation.
Example:
mat F = randu(3,3); mat e = Log_strain(F);
- Parameters:
F – 3x3 matrix representing the transformation gradient \(\mathbf{F}\)
- Returns:
3x3 matrix representing the logarithmic strain tensor \(\mathbf{e}\)
-
arma::mat finite_L(const arma::mat &F0, const arma::mat &F1, const double &DTime)
Provides the approximation of the Eulerian velocity tensor from the transformation gradient at two different times.
The Eulerian velocity tensor \(\mathbf{L}\) is related to the transformation gradients \(\mathbf{F_0}\) at time \(t_0\), \(\mathbf{F_1}\) at time \(t_1\), and the time difference \(\Delta t = t_1 - t_0\) by the following equation:
\[ \mathbf{L} = \frac{1}{\Delta t} \left( \mathbf{F}_1 - \mathbf{F}_0 \right) \cdot \mathbf{F}_1^{-1} \]Example:
mat F0 = randu(3,3); mat F1 = randu(3,3); double DTime = 0.1; mat L = finite_L(F0, F1, DTime);
- Parameters:
F0 – 3x3 matrix representing the transformation gradient \(\mathbf{F_0}\) at time \(t_0\)
F1 – 3x3 matrix representing the transformation gradient \(\mathbf{F_1}\) at time \(t_1\)
DTime – double representing the time difference Delta \(\Delta t = t_1 - t_0\)
- Returns:
3x3 matrix representing the approximation of the Eulerian velocity tensor
-
arma::mat finite_D(const arma::mat &F0, const arma::mat &F1, const double &DTime)
Provides the approximation of the Eulerian symmetric rate tensor from the transformation gradient at time \(t_0\), the transformation gradient at time \(t_1\) and the time difference \(\Delta t = t_1 - t_0\).
The Eulerian symmetric rate tensor \(\mathbf{D}\) is related to the transformation gradient \(\mathbf{F}_0\) at time \(t_0\), the transformation gradient \(\mathbf{F}_1\) at time \(t_1\) and the time difference \(\Delta t = t_1 - t_0\) by the following equation:
\[ \mathbf{D} = \frac{1}{2} \left( \mathbf{L} + \mathbf{L}^T \right) \]where \(\mathbf{L}\) is the Eulerian velocity tensor. \(\mathbf{D}\) is commonly referred as the rate of deformation (this necessitates although a specific discussion).
Example:
mat F0 = randu(3,3); mat F1 = randu(3,3); double DTime = 0.1; mat D = finite_D(F0, F1, DTime);
- Parameters:
F0 – 3x3 matrix representing the transformation gradient \(\mathbf{F}_0\) at time \(t_0\)
F1 – 3x3 matrix representing the transformation gradient \(\mathbf{F}_1\) at time \(t_1\)
DTime – time difference \(\Delta t = t_1 - t_0\)
- Returns:
3x3 matrix representing the approximation of the Eulerian symmetric rate tensor \(\mathbf{D}\)
-
arma::mat finite_W(const arma::mat &F0, const arma::mat &F1, const double &DTime)
Provides the approximation of the Eulerian antisymmetric spin tensor from the transformation gradient at time \(t_0\), the transformation gradient at time \(t_1\) and the time difference \(\Delta t = t_1 - t_0\).
The Eulerian antisymmetric spin tensor \(\mathbf{W}\) is related to the transformation gradient \(\mathbf{F}_0\) at time \(t_0\), the transformation gradient \(\mathbf{F}_1\) at time \(t_1\) and the time difference \(\Delta t = t_1 - t_0\) by the following equation:
\[ \mathbf{W} = \frac{1}{2} \left( \mathbf{L} - \mathbf{L}^T \right) \]where \(\mathbf{L}\) is the Eulerian velocity tensor. \(\mathbf{W}\) corresponds to the Jaumann corotationnal rate.Example:
mat F0 = randu(3,3); mat F1 = randu(3,3); mat DTime = 0.1; mat W = finite_W(F0, F1, DTime);
- Parameters:
F0 – 3x3 matrix representing the transformation gradient \(\mathbf{F_0}\) at time \(t_0\)
F1 – 3x3 matrix representing the transformation gradient \(\mathbf{F_1}\) at time \(t_1\)
DTime – time difference \(\Delta t = t_1 - t_0\)
- Returns:
3x3 matrix representing the approximation of the Eulerian antisymmetric spin tensor \(\mathbf{W}\)
-
arma::mat finite_Omega(const arma::mat &F0, const arma::mat &F1, const double &DTime)
Provides the approximation of the Eulerian rigid-body rotation spin tensor from the transformation gradient at time \(t_0\), the transformation gradient at time \(t_1\) and the time difference \(\Delta t = t_1 - t_0\).
The Eulerian rigid-body rotation spin tensor \(\mathbf{\Omega}\) is related to the transformation gradient \(\mathbf{F}_0\) at time \(t_0\), the transformation gradient \(\mathbf{F}_1\) at time \(t_1\) and the time difference \(\Delta t = t_1 - t_0\) by the following equation:
\[ \mathbf{\Omega} = \frac{1}{\Delta t} \left( \mathbf{R}_1 - \mathbf{R}_0 \right) \cdot \mathbf{R}_1^{T} \]where \(\mathbf{R}_0\) and \(\mathbf{R}_1\) are the rotation matrices obtained from RU decompositions of the transformation gradients \(\mathbf{F}_0\) and \(\mathbf{F}_1\), respectively. This corresponds to the Green-Naghdi corotationnal rate.
Example:
mat F0 = randu(3,3); mat F1 = randu(3,3); mat DTime = 0.1; mat Omega = finite_Omega(F0, F1, DTime);
- Parameters:
F0 – 3x3 matrix representing the transformation gradient \(\mathbf{F_0}\) at time \(t_0\)
F1 – 3x3 matrix representing the transformation gradient \(\mathbf{F_1}\) at time \(t_1\)
DTime – time difference \(\Delta t = t_1 - t_0\)
- Returns:
3x3 matrix representing the approximation of the Eulerian rigid-body rotation spin tensor (Green-Naghdi corotational rate).
-
arma::mat finite_DQ(const arma::mat &Omega0, const arma::mat &Omega1, const double &DTime)
Provides the Hughes-Winget approximation of a increment of rotation or transformation from the spin/velocity at time \(t_0\), the spin/velocity at time \(t_1\) and the time difference \(\Delta t = t_1 - t_0\).
The Hughes-Winget approximation of a increment of rotation or transformation \(\mathbf{\Delta Q}\) is related to the spin/velocity \(\mathbf{\Omega}_0\) at time \(t_0\), the spin/velocity \(\mathbf{\Omega}_1\) at time \(t_1\) and the time difference \(\Delta t = t_1 - t_0\) by the following equation:
\[ \mathbf{\Delta Q} = \left( \mathbf{I} + \frac{1}{2} \Delta t \, \mathbf{\Omega}_0 \right) \cdot \left( \mathbf{I} - \frac{1}{2} \Delta t \, \mathbf{\Omega}_1 \right)^{-1} \]where \(\mathbf{I}\) is the 3x3 identity matrix.
Example:
mat Omega0 = randu(3,3); mat Omega1 = randu(3,3); mat DTime = 0.1; mat DQ = finite_DQ(Omega0, Omega1, DTime);
- Parameters:
Omega0 – spin/velocity at time \(t_0\)
Omega1 – spin/velocity at time \(t_1\)
DTime – time difference \(\Delta t = t_1 - t_0\)
- Returns:
3x3 matrix representing the Hughes-Winget approximation of a increment of rotation or transformation