DAMASK_EICMD/python/damask/mechanics.py

374 lines
8.9 KiB
Python
Raw Normal View History

"""
Finite-strain continuum mechanics.
All routines operate on numpy.ndarrays of shape (...,3,3).
"""
2023-02-21 20:57:06 +05:30
from typing import Sequence as _Sequence#, Literal as _Literal
import numpy as _np
2021-11-01 03:20:41 +05:30
from . import tensor as _tensor
from . import _rotation
2021-11-01 03:20:41 +05:30
def deformation_Cauchy_Green_left(F: _np.ndarray) -> _np.ndarray:
2023-02-17 00:00:28 +05:30
r"""
Calculate left Cauchy-Green deformation tensor (Finger deformation tensor).
2019-10-30 22:35:44 +05:30
Parameters
----------
2021-12-06 12:08:40 +05:30
F : numpy.ndarray, shape (...,3,3)
Deformation gradient.
Returns
-------
2021-12-06 12:08:40 +05:30
B : numpy.ndarray, shape (...,3,3)
2020-11-19 19:08:54 +05:30
Left Cauchy-Green deformation tensor.
2019-10-30 22:35:44 +05:30
2023-02-17 00:00:28 +05:30
Notes
-----
.. math::
2023-02-17 01:42:10 +05:30
\vb{B} = \vb{F} \vb{F}^\text{T}
2023-02-17 00:00:28 +05:30
2019-10-30 22:35:44 +05:30
"""
2020-11-19 18:35:59 +05:30
return _np.matmul(F,_tensor.transpose(F))
2019-11-27 16:49:37 +05:30
2021-11-01 03:20:41 +05:30
def deformation_Cauchy_Green_right(F: _np.ndarray) -> _np.ndarray:
2023-02-17 00:00:28 +05:30
r"""
Calculate right Cauchy-Green deformation tensor.
2019-11-27 16:49:37 +05:30
Parameters
----------
2021-12-06 12:08:40 +05:30
F : numpy.ndarray, shape (...,3,3)
Deformation gradient.
Returns
-------
2021-12-06 12:08:40 +05:30
C : numpy.ndarray, shape (...,3,3)
2020-11-19 19:08:54 +05:30
Right Cauchy-Green deformation tensor.
2019-11-27 16:49:37 +05:30
2023-02-17 00:00:28 +05:30
Notes
-----
.. math::
2023-02-17 01:42:10 +05:30
\vb{C} = \vb{F}^\text{T} \vb{F}
2023-02-17 00:00:28 +05:30
2019-11-27 16:49:37 +05:30
"""
2020-11-19 18:35:59 +05:30
return _np.matmul(_tensor.transpose(F),F)
2021-11-01 03:20:41 +05:30
def equivalent_strain_Mises(epsilon: _np.ndarray) -> _np.ndarray:
2023-02-14 21:40:22 +05:30
r"""
Calculate the Mises equivalent of a strain tensor.
2019-10-30 22:35:44 +05:30
Parameters
----------
2021-12-06 12:08:40 +05:30
epsilon : numpy.ndarray, shape (...,3,3)
Symmetric strain tensor of which the von Mises equivalent is computed.
2019-10-19 16:40:46 +05:30
Returns
-------
2021-12-06 12:08:40 +05:30
epsilon_vM : numpy.ndarray, shape (...)
Von Mises equivalent strain of epsilon.
2023-02-14 21:40:22 +05:30
Notes
-----
2023-02-17 00:00:28 +05:30
The von Mises equivalent of a strain tensor is defined as:
2023-02-14 21:40:22 +05:30
.. math::
2023-02-21 20:57:06 +05:30
\epsilon_\text{vM} = \sqrt{\frac{2}{3}\,\epsilon^\prime_{ij} \epsilon^\prime_{ij}}
2023-02-14 21:40:22 +05:30
2023-02-17 01:42:10 +05:30
where :math:`\vb*{\epsilon}^\prime` is the deviatoric part
2023-02-17 00:00:28 +05:30
of the strain tensor.
2019-10-30 22:35:44 +05:30
"""
2020-11-18 03:26:22 +05:30
return _equivalent_Mises(epsilon,2.0/3.0)
2021-11-01 03:20:41 +05:30
def equivalent_stress_Mises(sigma: _np.ndarray) -> _np.ndarray:
2023-02-14 21:40:22 +05:30
r"""
Calculate the Mises equivalent of a stress tensor.
2019-10-30 22:35:44 +05:30
Parameters
----------
2021-12-06 12:08:40 +05:30
sigma : numpy.ndarray, shape (...,3,3)
Symmetric stress tensor of which the von Mises equivalent is computed.
2019-10-19 16:40:46 +05:30
Returns
-------
2021-12-06 12:08:40 +05:30
sigma_vM : numpy.ndarray, shape (...)
Von Mises equivalent stress of sigma.
2023-02-14 21:40:22 +05:30
Notes
-----
2023-02-17 00:00:28 +05:30
The von Mises equivalent of a stress tensor is defined as:
2023-02-14 21:40:22 +05:30
.. math::
2023-02-21 20:57:06 +05:30
\sigma_\text{vM} = \sqrt{\frac{3}{2}\,\sigma^\prime_{ij} \sigma^\prime_{ij}}
2023-02-14 21:40:22 +05:30
2023-02-17 01:42:10 +05:30
where :math:`\vb*{\sigma}^\prime` is the deviatoric part
of the stress tensor.
2023-02-17 00:00:28 +05:30
2019-10-30 22:35:44 +05:30
"""
2020-11-18 03:26:22 +05:30
return _equivalent_Mises(sigma,3.0/2.0)
2021-11-01 03:20:41 +05:30
def maximum_shear(T_sym: _np.ndarray) -> _np.ndarray:
2019-10-30 22:35:44 +05:30
"""
2020-11-19 18:35:59 +05:30
Calculate the maximum shear component of a symmetric tensor.
2019-10-30 22:35:44 +05:30
Parameters
----------
2021-12-06 12:08:40 +05:30
T_sym : numpy.ndarray, shape (...,3,3)
2020-11-19 18:35:59 +05:30
Symmetric tensor of which the maximum shear is computed.
2019-10-19 16:40:46 +05:30
Returns
-------
2021-12-06 12:08:40 +05:30
gamma_max : numpy.ndarray, shape (...)
2020-11-19 18:35:59 +05:30
Maximum shear of T_sym.
2019-10-30 22:35:44 +05:30
"""
2020-11-19 18:35:59 +05:30
w = _tensor.eigenvalues(T_sym)
return (w[...,0] - w[...,2])*0.5
2020-03-03 03:41:05 +05:30
2021-11-01 03:20:41 +05:30
def rotation(T: _np.ndarray) -> _rotation.Rotation:
2023-02-14 21:40:22 +05:30
r"""
Calculate the rotational part of a tensor.
Parameters
----------
2021-12-06 12:08:40 +05:30
T : numpy.ndarray, shape (...,3,3)
Tensor of which the rotational part is computed.
Returns
-------
2021-12-06 12:08:40 +05:30
R : damask.Rotation, shape (...)
Rotational part of the vector.
2023-02-14 21:40:22 +05:30
Notes
-----
The rotational part is calculated from the polar decomposition:
.. math::
2023-02-16 22:29:19 +05:30
\vb{R} = \vb{T} \vb{U}^{-1} = \vb{V}^{-1} \vb{T}
2023-02-14 21:40:22 +05:30
2023-02-17 01:42:10 +05:30
where :math:`\vb{V}` and :math:`\vb{U}` are the left
2023-02-14 21:40:22 +05:30
and right stretch tensor, respectively.
"""
return _rotation.Rotation.from_matrix(_polar_decomposition(T,'R')[0])
2022-01-27 19:59:33 +05:30
def strain(F: _np.ndarray,
2023-02-21 20:57:06 +05:30
#t: _Literal['V', 'U'], should work, but rejected by SC
2022-01-27 19:59:33 +05:30
t: str,
m: float) -> _np.ndarray:
2023-02-21 20:57:06 +05:30
r"""
2020-11-16 05:42:23 +05:30
Calculate strain tensor (SethHill family).
2020-02-15 18:40:16 +05:30
2019-10-30 22:35:44 +05:30
Parameters
----------
2021-12-06 12:08:40 +05:30
F : numpy.ndarray, shape (...,3,3)
Deformation gradient.
t : {'V', 'U'}
Type of the polar decomposition, 'V' for left stretch tensor
or 'U' for right stretch tensor.
2020-02-15 18:40:16 +05:30
m : float
Order of the strain.
2019-10-25 17:00:20 +05:30
Returns
-------
2021-12-06 12:08:40 +05:30
epsilon : numpy.ndarray, shape (...,3,3)
Strain of F.
2023-02-21 20:57:06 +05:30
Notes
-----
The strain is defined as:
.. math::
\vb*{\epsilon}_V^{(m)} = \frac{1}{2m} (\vb{V}^{2m} - \vb{I}) \\\\
\vb*{\epsilon}_U^{(m)} = \frac{1}{2m} (\vb{U}^{2m} - \vb{I})
2021-05-07 23:12:23 +05:30
References
----------
2023-02-21 20:57:06 +05:30
| https://en.wikipedia.org/wiki/Finite_strain_theory
| https://de.wikipedia.org/wiki/Verzerrungstensor
2021-05-07 23:12:23 +05:30
2019-10-30 22:35:44 +05:30
"""
if t not in ['V', 'U']: raise ValueError('polar decomposition type not in {V, U}')
w,n = _np.linalg.eigh(deformation_Cauchy_Green_left(F) if t=='V' else deformation_Cauchy_Green_right(F))
return 0.5 * _np.einsum('...j,...kj,...lj',_np.log(w),n,n) if m == 0.0 \
else 0.5/m * (_np.einsum('...j,...kj,...lj', w**m,n,n) - _np.eye(3))
2020-02-15 18:40:16 +05:30
2022-01-27 19:59:33 +05:30
def stress_Cauchy(P: _np.ndarray,
F: _np.ndarray) -> _np.ndarray:
2020-11-19 18:35:59 +05:30
"""
Calculate the Cauchy stress (true stress).
Resulting tensor is symmetrized as the Cauchy stress needs to be symmetric.
Parameters
----------
2021-12-06 12:08:40 +05:30
P : numpy.ndarray, shape (...,3,3)
2020-11-19 18:35:59 +05:30
First Piola-Kirchhoff stress.
2021-12-06 12:08:40 +05:30
F : numpy.ndarray, shape (...,3,3)
2020-11-19 18:35:59 +05:30
Deformation gradient.
Returns
-------
2021-12-06 12:08:40 +05:30
sigma : numpy.ndarray, shape (...,3,3)
2020-11-19 18:35:59 +05:30
Cauchy stress.
"""
return _tensor.symmetric(_np.einsum('...,...ij,...kj',1.0/_np.linalg.det(F),P,F))
2022-01-27 19:59:33 +05:30
def stress_second_Piola_Kirchhoff(P: _np.ndarray,
F: _np.ndarray) -> _np.ndarray:
2020-11-19 18:35:59 +05:30
"""
Calculate the second Piola-Kirchhoff stress.
Resulting tensor is symmetrized as the second Piola-Kirchhoff stress
needs to be symmetric.
Parameters
----------
2021-12-06 12:08:40 +05:30
P : numpy.ndarray, shape (...,3,3)
2020-11-19 18:35:59 +05:30
First Piola-Kirchhoff stress.
2021-12-06 12:08:40 +05:30
F : numpy.ndarray, shape (...,3,3)
2020-11-19 18:35:59 +05:30
Deformation gradient.
Returns
-------
2021-12-06 12:08:40 +05:30
S : numpy.ndarray, shape (...,3,3)
2020-11-19 18:35:59 +05:30
Second Piola-Kirchhoff stress.
"""
return _tensor.symmetric(_np.einsum('...ij,...jk',_np.linalg.inv(F),P))
2021-11-01 03:20:41 +05:30
def stretch_left(T: _np.ndarray) -> _np.ndarray:
2023-02-14 21:40:22 +05:30
r"""
2020-11-19 19:08:54 +05:30
Calculate left stretch of a tensor.
Parameters
----------
2021-12-06 12:08:40 +05:30
T : numpy.ndarray, shape (...,3,3)
Tensor of which the left stretch is computed.
Returns
-------
2021-12-06 12:08:40 +05:30
V : numpy.ndarray, shape (...,3,3)
Left stretch tensor from Polar decomposition of T.
2023-02-14 21:40:22 +05:30
Notes
-----
The left stretch tensor is calculated from the
polar decomposition:
.. math::
2023-02-16 22:29:19 +05:30
\vb{V} = \vb{T} \vb{R}^\text{T}
2023-02-14 21:40:22 +05:30
2023-02-16 22:29:19 +05:30
where :math:`\vb{R}` is a rotation.
2023-02-14 21:40:22 +05:30
"""
return _polar_decomposition(T,'V')[0]
2021-11-01 03:20:41 +05:30
def stretch_right(T: _np.ndarray) -> _np.ndarray:
2023-02-14 21:40:22 +05:30
r"""
2020-11-19 19:08:54 +05:30
Calculate right stretch of a tensor.
Parameters
----------
2021-12-06 12:08:40 +05:30
T : numpy.ndarray, shape (...,3,3)
Tensor of which the right stretch is computed.
Returns
-------
2021-12-06 12:08:40 +05:30
U : numpy.ndarray, shape (...,3,3)
Left stretch tensor from Polar decomposition of T.
2023-02-14 21:40:22 +05:30
Notes
-----
The right stretch tensor is calculated from the
polar decomposition:
.. math::
2023-02-16 22:29:19 +05:30
\vb{U} = \vb{R}^\text{T} \vb{T}
2023-02-14 21:40:22 +05:30
2023-02-16 22:29:19 +05:30
where :math:`\vb{R}` is a rotation.
2023-02-14 21:40:22 +05:30
"""
return _polar_decomposition(T,'U')[0]
2022-01-27 19:59:33 +05:30
def _polar_decomposition(T: _np.ndarray,
requested: _Sequence[str]) -> tuple:
2019-10-30 22:35:44 +05:30
"""
Perform singular value decomposition.
2019-10-30 22:35:44 +05:30
Parameters
----------
2021-12-06 12:08:40 +05:30
T : numpy.ndarray, shape (...,3,3)
Tensor of which the singular values are computed.
2022-01-27 19:59:33 +05:30
requested : sequence of {'R', 'U', 'V'}
2023-02-21 20:57:06 +05:30
Requested outputs: 'R' for the rotation tensor,
'V' for left stretch tensor, and 'U' for right stretch tensor.
2019-10-30 22:35:44 +05:30
2023-02-14 17:24:51 +05:30
Returns
-------
VRU : tuple of numpy.ndarray, shape (...,3,3)
Requested components of the singular value decomposition.
2019-10-30 22:35:44 +05:30
"""
u, _, vh = _np.linalg.svd(T)
2020-11-16 11:42:37 +05:30
R = _np.einsum('...ij,...jk',u,vh)
2021-11-01 03:20:41 +05:30
output = []
2019-10-30 22:35:44 +05:30
if 'R' in requested:
2021-11-01 03:20:41 +05:30
output+=[R]
2019-10-30 22:35:44 +05:30
if 'V' in requested:
2021-11-01 03:20:41 +05:30
output+=[_np.einsum('...ij,...kj',T,R)]
2019-10-30 22:35:44 +05:30
if 'U' in requested:
2021-11-01 03:20:41 +05:30
output+=[_np.einsum('...ji,...jk',R,T)]
if len(output) == 0:
raise ValueError('output not in {V, R, U}')
2021-11-01 03:20:41 +05:30
return tuple(output)
2020-02-15 18:40:16 +05:30
2022-01-27 19:59:33 +05:30
def _equivalent_Mises(T_sym: _np.ndarray,
s: float) -> _np.ndarray:
2020-02-15 18:40:16 +05:30
"""
2020-11-16 11:42:37 +05:30
Base equation for Mises equivalent of a stress or strain tensor.
2020-02-15 18:40:16 +05:30
Parameters
----------
2021-12-06 12:08:40 +05:30
T_sym : numpy.ndarray, shape (...,3,3)
Symmetric tensor of which the von Mises equivalent is computed.
2020-02-15 18:40:16 +05:30
s : float
Scaling factor (2/3 for strain, 3/2 for stress).
2020-03-03 03:41:05 +05:30
2023-02-14 17:24:51 +05:30
Returns
-------
eq : numpy.ndarray, shape (...)
Scaled second invariant of the deviatoric part of T_sym.
2020-02-15 18:40:16 +05:30
"""
2020-11-19 19:08:54 +05:30
d = _tensor.deviatoric(T_sym)
2020-11-16 11:42:37 +05:30
return _np.sqrt(s*_np.sum(d**2.0,axis=(-1,-2)))