Merge branch 'extra-documentation' into 'development'

improved docstrings

See merge request damask/DAMASK!722
This commit is contained in:
Daniel Otto de Mentock 2023-02-17 12:58:09 +00:00
commit 6fd30dd266
1 changed files with 84 additions and 7 deletions

View File

@ -14,7 +14,7 @@ from . import _rotation
def deformation_Cauchy_Green_left(F: _np.ndarray) -> _np.ndarray:
"""
r"""
Calculate left Cauchy-Green deformation tensor (Finger deformation tensor).
Parameters
@ -27,12 +27,18 @@ def deformation_Cauchy_Green_left(F: _np.ndarray) -> _np.ndarray:
B : numpy.ndarray, shape (...,3,3)
Left Cauchy-Green deformation tensor.
Notes
-----
.. math::
\vb{B} = \vb{F} \vb{F}^\text{T}
"""
return _np.matmul(F,_tensor.transpose(F))
def deformation_Cauchy_Green_right(F: _np.ndarray) -> _np.ndarray:
"""
r"""
Calculate right Cauchy-Green deformation tensor.
Parameters
@ -45,12 +51,18 @@ def deformation_Cauchy_Green_right(F: _np.ndarray) -> _np.ndarray:
C : numpy.ndarray, shape (...,3,3)
Right Cauchy-Green deformation tensor.
Notes
-----
.. math::
\vb{C} = \vb{F}^\text{T} \vb{F}
"""
return _np.matmul(_tensor.transpose(F),F)
def equivalent_strain_Mises(epsilon: _np.ndarray) -> _np.ndarray:
"""
r"""
Calculate the Mises equivalent of a strain tensor.
Parameters
@ -63,12 +75,23 @@ def equivalent_strain_Mises(epsilon: _np.ndarray) -> _np.ndarray:
epsilon_vM : numpy.ndarray, shape (...)
Von Mises equivalent strain of epsilon.
Notes
-----
The von Mises equivalent of a strain tensor is defined as:
.. math::
\epsilon_\text{vM} = \sqrt{2/3 \epsilon^\prime_{ij} \epsilon^\prime_{ij}}
where :math:`\vb*{\epsilon}^\prime` is the deviatoric part
of the strain tensor.
"""
return _equivalent_Mises(epsilon,2.0/3.0)
def equivalent_stress_Mises(sigma: _np.ndarray) -> _np.ndarray:
"""
r"""
Calculate the Mises equivalent of a stress tensor.
Parameters
@ -81,6 +104,17 @@ def equivalent_stress_Mises(sigma: _np.ndarray) -> _np.ndarray:
sigma_vM : numpy.ndarray, shape (...)
Von Mises equivalent stress of sigma.
Notes
-----
The von Mises equivalent of a stress tensor is defined as:
.. math::
\sigma_\text{vM} = \sqrt{3/2 \sigma^\prime_{ij} \sigma^\prime_{ij}}
where :math:`\vb*{\sigma}^\prime` is the deviatoric part
of the stress tensor.
"""
return _equivalent_Mises(sigma,3.0/2.0)
@ -105,7 +139,7 @@ def maximum_shear(T_sym: _np.ndarray) -> _np.ndarray:
def rotation(T: _np.ndarray) -> _rotation.Rotation:
"""
r"""
Calculate the rotational part of a tensor.
Parameters
@ -118,6 +152,17 @@ def rotation(T: _np.ndarray) -> _rotation.Rotation:
R : damask.Rotation, shape (...)
Rotational part of the vector.
Notes
-----
The rotational part is calculated from the polar decomposition:
.. math::
\vb{R} = \vb{T} \vb{U}^{-1} = \vb{V}^{-1} \vb{T}
where :math:`\vb{V}` and :math:`\vb{U}` are the left
and right stretch tensor, respectively.
"""
return _rotation.Rotation.from_matrix(_polar_decomposition(T,'R')[0])
@ -212,7 +257,7 @@ def stress_second_Piola_Kirchhoff(P: _np.ndarray,
def stretch_left(T: _np.ndarray) -> _np.ndarray:
"""
r"""
Calculate left stretch of a tensor.
Parameters
@ -225,12 +270,23 @@ def stretch_left(T: _np.ndarray) -> _np.ndarray:
V : numpy.ndarray, shape (...,3,3)
Left stretch tensor from Polar decomposition of T.
Notes
-----
The left stretch tensor is calculated from the
polar decomposition:
.. math::
\vb{V} = \vb{T} \vb{R}^\text{T}
where :math:`\vb{R}` is a rotation.
"""
return _polar_decomposition(T,'V')[0]
def stretch_right(T: _np.ndarray) -> _np.ndarray:
"""
r"""
Calculate right stretch of a tensor.
Parameters
@ -243,6 +299,17 @@ def stretch_right(T: _np.ndarray) -> _np.ndarray:
U : numpy.ndarray, shape (...,3,3)
Left stretch tensor from Polar decomposition of T.
Notes
-----
The right stretch tensor is calculated from the
polar decomposition:
.. math::
\vb{U} = \vb{R}^\text{T} \vb{T}
where :math:`\vb{R}` is a rotation.
"""
return _polar_decomposition(T,'U')[0]
@ -260,6 +327,11 @@ def _polar_decomposition(T: _np.ndarray,
Requested outputs: R for the rotation tensor,
V for left stretch tensor, and U for right stretch tensor.
Returns
-------
VRU : tuple of numpy.ndarray, shape (...,3,3)
Requested components of the singular value decomposition.
"""
u, _, vh = _np.linalg.svd(T)
R = _np.einsum('...ij,...jk',u,vh)
@ -290,6 +362,11 @@ def _equivalent_Mises(T_sym: _np.ndarray,
s : float
Scaling factor (2/3 for strain, 3/2 for stress).
Returns
-------
eq : numpy.ndarray, shape (...)
Scaled second invariant of the deviatoric part of T_sym.
"""
d = _tensor.deviatoric(T_sym)
return _np.sqrt(s*_np.sum(d**2.0,axis=(-1,-2)))