The Kinematics Library
-
mat ER_to_F(const mat &E, const mat &R)
Provides the transformation gradient \(\mathbf{F}\), from the Green-Lagrange strain \(\mathbf{E}\) and the rotation \(\mathbf{R}\):
\[\mathbf{F} = \mathbf{R} \cdot \mathbf{U} \quad \mathbf{E} = \frac{1}{2} \left( \sqrt{\mathbf{U}^2} - \mathbf{I} \right)\]mat E = randu(3,3); mat R = eye(3,3); mat F = ER_to_F(E, R);
-
mat eR_to_F(const mat &e, const mat &R)
Provides the transformation gradient \(\mathbf{F}\), from the logarithmic strain \(\mathbf{e}\) and the rotation \(\mathbf{R}\):
\[\mathbf{F} = \mathbf{V} \cdot \mathbf{R} \quad \mathbf{e} = \textrm{ln} \mathbf{V}\]mat e = randu(3,3); mat R = eye(3,3); mat F = eR_to_F(e, R);
-
mat G_UdX(const mat &F)
Provides the gradient of the displacement (Lagrangian) from the transformation gradient \(\mathbf{F}\):
\[\nabla_X \mathbf{U} = \mathbf{F} - \mathbf{I}\]mat F = randu(3,3); mat GradU = G_UdX(F);
-
mat G_Udx(const mat &F)
Provides the gradient of the displacement (Eulerian) from the transformation gradient \(\mathbf{F}\):
\[\nabla_x \mathbf{U} = \mathbf{I} - \mathbf{F}^{-1}\]mat F = randu(3,3); mat gradU = G_UdX(F);
-
mat R_Cauchy_Green(const mat &F)
Provides the Right Cauchy-Green tensor \(\mathbf{C}\): from the transformation gradient \(\mathbf{F}\):
\[\mathbf{C} = \mathbf{F}^T \cdot \mathbf{F}\]mat F = randu(3,3); mat C = R_Cauchy_Green(F);
-
mat L_Cauchy_Green(const mat &F)
Provides the Left Cauchy-Green tensor \(\mathbf{B}\): from the transformation gradient \(\mathbf{F}\):
\[\mathbf{B} = \mathbf{F} \cdot \mathbf{F}^T\]mat F = randu(3,3); mat B = L_Cauchy_Green(F);
-
RU_decomposition(mat &R, mat &U, const mat &F)
Provides the RU decomposition of the transformation gradient \(\mathbf{F}\):
\[\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}\]mat F = randu(3,3); mat R = zeros(3,3); mat U = zeros(3,3); RU_decomposition(R, U, F);
-
VR_decomposition(mat &R, mat &V, const mat &F)
Provides the VR decomposition of the transformation gradient \(\mathbf{F}\):
mat F = randu(3,3);
mat R = zeros(3,3);
mat V = zeros(3,3);
VR_decomposition(R, V, F);
-
vec Inv_X(const mat &X)
Provides the invariants of a symmetric tensor \(\mathbf{X}\):
\[\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)\]mat F = randu(3,3); mat C = R_Cauchy_Green(F); vec I = Inv_X(F);
-
mat Cauchy(const mat &F)
Provides the Cauchy tensor \(\mathbf{b}\): from the transformation gradient \(\mathbf{F}\):
\[\mathbf{b} = \left( \mathbf{F} \cdot \mathbf{F}^T \right)^{-1}\]mat F = randu(3,3); mat b = Cauchy(F);
-
mat Green_Lagrange(const mat &F)
Provides the Green-Lagrange tensor \(\mathbf{E}\): from the transformation gradient \(\mathbf{F}\):
\[\mathbf{E} = \frac{1}{2} \left( \mathbf{F}^T \cdot \mathbf{F} - \mathbf{I} \right)\]mat F = randu(3,3); mat E = Green_Lagrange(F);
-
mat Euler_Almansi(const mat &F)
Provides the Euler_Almansi tensor \(\mathbf{e}\): from the transformation gradient \(\mathbf{F}\):
\[\mathbf{e} = \frac{1}{2} \left( \mathbf{I} - \left( \mathbf{F} \cdot \mathbf{F}^T \right)^T \right)\]mat F = randu(3,3); mat e = Euler_Almansi(F);
-
mat Log_strain(const mat &F)
Provides the logarithmic strain tensor \(\mathbf{h}\): from the transformation gradient \(\mathbf{F}\):
mat F = randu(3,3);
mat h = Log_strain(F);
-
mat finite_L(const mat &F0, const mat &F1, const double &DTime)
Provides the approximation of the Eulerian velocity tensor \(\mathbf{L}\): from the transformation gradient \(\mathbf{F}_0\): at time \(t_0\):, \(\mathbf{F}_1\): at time \(t_1\): and the time difference \(\Delta t = t_1 - t_0\):
mat F0 = randu(3,3);
mat F1 = randu(3,3);
mat DTime = 0.1;
mat L = finite_L(F0, F1, DTime);
-
mat finite_D(const mat &F0, const mat &F1, const double &DTime)
Provides the approximation of the Eulerian symmetric rate tensor \(\mathbf{D}\): from the transformation gradient \(\mathbf{F}_0\): at time \(t_0\):, \(\mathbf{F}_1\): at time \(t_1\): and the time difference \(\Delta t = t_1 - t_0\): This is commonly referred as the rate of deformation (this necessitates although a specific discussion)
mat F0 = randu(3,3);
mat F1 = randu(3,3);
mat DTime = 0.1;
mat D = finite_D(F0, F1, DTime);
-
mat finite_W(const mat &F0, const mat &F1, const double &DTime)
Provides the approximation of the Eulerian antisymmetric spin tensor \(\mathbf{W}\): from the transformation gradient \(\mathbf{F}_0\): at time \(t_0\):, \(\mathbf{F}_1\): at time \(t_1\): and the time difference \(\Delta t = t_1 - t_0\): . This correspond to the Jaumann corotationnal rate:
mat F0 = randu(3,3);
mat F1 = randu(3,3);
mat DTime = 0.1;
mat W = finite_W(F0, F1, DTime);
-
mat finite_Omega(const mat &F0, const mat &F1, const double &DTime)
Provides the approximation of the Eulerian rigid-body rotation spin tensor \(\mathbf{\Omega}\): from the transformation gradient \(\mathbf{F}_0\): at time \(t_0\):, \(\mathbf{F}_1\): at time \(t_1\): and the time difference \(\Delta t = t_1 - t_0\): . This correspond to the Green-Naghdi corotationnal rate:
Note that the rotation matrix is obtained from a RU decomposition of the transformation gradient \(\mathbf{F}\):
mat F0 = randu(3,3);
mat F1 = randu(3,3);
mat DTime = 0.1;
mat Omega = finite_Omega(F0, F1, DTime);
-
mat finite_DQ(const mat &Omega0, const mat &Omega1, const double &DTime)
Provides the Hughes-Winget approximation of a increment of rotation or transformation from the spin/velocity \(\mathbf{\Omega}_0\): at time \(t_0\):, \(\mathbf{\Omega}_1\): at time \(t_1\): and the time difference \(\Delta t = t_1 - t_0\):
mat Omega0 = randu(3,3);
mat Omega1 = randu(3,3);
mat DTime = 0.1;
mat DQ = finite_DQ(Omega0, Omega1, DTime);