DAMASK_EICMD/python/damask/mechanics.py

290 lines
7.1 KiB
Python
Raw Normal View History

"""
Finite-strain continuum mechanics.
All routines operate on numpy.ndarrays of shape (...,3,3).
"""
2021-11-01 03:20:41 +05:30
from typing import Sequence
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:
2019-10-30 22:35:44 +05:30
"""
Calculate left Cauchy-Green deformation tensor (Finger deformation tensor).
2019-10-30 22:35:44 +05:30
Parameters
----------
2020-06-10 00:57:08 +05:30
F : numpy.ndarray of shape (...,3,3)
Deformation gradient.
Returns
-------
B : numpy.ndarray of shape (...,3,3)
2020-11-19 19:08:54 +05:30
Left Cauchy-Green deformation tensor.
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:
2019-11-27 16:49:37 +05:30
"""
Calculate right Cauchy-Green deformation tensor.
2019-11-27 16:49:37 +05:30
Parameters
----------
F : numpy.ndarray of shape (...,3,3)
Deformation gradient.
Returns
-------
C : numpy.ndarray of shape (...,3,3)
2020-11-19 19:08:54 +05:30
Right Cauchy-Green deformation tensor.
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:
2019-10-30 22:35:44 +05:30
"""
Calculate the Mises equivalent of a strain tensor.
2019-10-30 22:35:44 +05:30
Parameters
----------
2020-06-10 00:57:08 +05:30
epsilon : numpy.ndarray of shape (...,3,3)
Symmetric strain tensor of which the von Mises equivalent is computed.
2019-10-19 16:40:46 +05:30
Returns
-------
epsilon_vM : numpy.ndarray of shape (...)
Von Mises equivalent strain of epsilon.
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:
2019-10-30 22:35:44 +05:30
"""
Calculate the Mises equivalent of a stress tensor.
2019-10-30 22:35:44 +05:30
Parameters
----------
2020-06-10 00:57:08 +05:30
sigma : numpy.ndarray of shape (...,3,3)
Symmetric stress tensor of which the von Mises equivalent is computed.
2019-10-19 16:40:46 +05:30
Returns
-------
sigma_vM : numpy.ndarray of shape (...)
Von Mises equivalent stress of sigma.
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
----------
2020-11-19 18:35:59 +05:30
T_sym : numpy.ndarray of shape (...,3,3)
Symmetric tensor of which the maximum shear is computed.
2019-10-19 16:40:46 +05:30
Returns
-------
2020-11-19 18:35:59 +05:30
gamma_max : numpy.ndarray of shape (...)
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:
"""
Calculate the rotational part of a tensor.
Parameters
----------
2020-06-10 00:57:08 +05:30
T : numpy.ndarray of shape (...,3,3)
Tensor of which the rotational part is computed.
Returns
-------
R : damask.Rotation of shape (...)
Rotational part of the vector.
"""
return _rotation.Rotation.from_matrix(_polar_decomposition(T,'R')[0])
2021-11-01 03:20:41 +05:30
def strain(F: _np.ndarray, t: str, m: float) -> _np.ndarray:
2019-10-30 22:35:44 +05:30
"""
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
----------
2020-06-10 00:57:08 +05:30
F : numpy.ndarray of shape (...,3,3)
Deformation gradient.
2020-02-15 18:40:16 +05:30
t : {V, U}
Type of the polar decomposition, V for left stretch tensor
2020-11-19 19:08:54 +05:30
and 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
-------
epsilon : numpy.ndarray of shape (...,3,3)
Strain of F.
2021-05-07 23:12:23 +05:30
References
----------
https://en.wikipedia.org/wiki/Finite_strain_theory
https://de.wikipedia.org/wiki/Verzerrungstensor
2019-10-30 22:35:44 +05:30
"""
2020-02-15 18:40:16 +05:30
if t == 'V':
2020-11-18 03:26:22 +05:30
w,n = _np.linalg.eigh(deformation_Cauchy_Green_left(F))
2020-02-15 18:40:16 +05:30
elif t == 'U':
2020-11-18 03:26:22 +05:30
w,n = _np.linalg.eigh(deformation_Cauchy_Green_right(F))
2019-10-25 17:00:20 +05:30
2020-02-15 18:40:16 +05:30
if m > 0.0:
2020-11-16 11:42:37 +05:30
eps = 1.0/(2.0*abs(m)) * (+ _np.einsum('...j,...kj,...lj',w**m,n,n) - _np.eye(3))
2020-02-15 18:40:16 +05:30
elif m < 0.0:
2020-11-16 11:42:37 +05:30
eps = 1.0/(2.0*abs(m)) * (- _np.einsum('...j,...kj,...lj',w**m,n,n) + _np.eye(3))
2020-02-15 18:40:16 +05:30
else:
2020-11-16 11:42:37 +05:30
eps = _np.einsum('...j,...kj,...lj',0.5*_np.log(w),n,n)
2019-10-25 17:00:20 +05:30
2020-06-10 00:57:08 +05:30
return eps
2020-02-15 18:40:16 +05:30
2021-11-01 03:20:41 +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
----------
P : numpy.ndarray of shape (...,3,3)
First Piola-Kirchhoff stress.
F : numpy.ndarray of shape (...,3,3)
Deformation gradient.
Returns
-------
sigma : numpy.ndarray of shape (...,3,3)
Cauchy stress.
"""
return _tensor.symmetric(_np.einsum('...,...ij,...kj',1.0/_np.linalg.det(F),P,F))
2021-11-01 03:20:41 +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
----------
P : numpy.ndarray of shape (...,3,3)
First Piola-Kirchhoff stress.
F : numpy.ndarray of shape (...,3,3)
Deformation gradient.
Returns
-------
S : numpy.ndarray of shape (...,3,3)
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:
"""
2020-11-19 19:08:54 +05:30
Calculate left stretch of a tensor.
Parameters
----------
T : numpy.ndarray of shape (...,3,3)
Tensor of which the left stretch is computed.
Returns
-------
V : numpy.ndarray of shape (...,3,3)
Left stretch tensor from Polar decomposition of T.
"""
return _polar_decomposition(T,'V')[0]
2021-11-01 03:20:41 +05:30
def stretch_right(T: _np.ndarray) -> _np.ndarray:
"""
2020-11-19 19:08:54 +05:30
Calculate right stretch of a tensor.
Parameters
----------
T : numpy.ndarray of shape (...,3,3)
Tensor of which the right stretch is computed.
Returns
-------
U : numpy.ndarray of shape (...,3,3)
Left stretch tensor from Polar decomposition of T.
"""
return _polar_decomposition(T,'U')[0]
2021-11-01 03:20:41 +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
----------
2020-06-10 00:57:08 +05:30
T : numpy.ndarray of shape (...,3,3)
Tensor of which the singular values are computed.
requested : iterable of str
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
"""
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 needs to be out of V, R, U')
2021-11-01 03:20:41 +05:30
return tuple(output)
2020-02-15 18:40:16 +05:30
2021-11-01 03:20:41 +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
----------
2020-06-10 00:57:08 +05:30
T_sym : numpy.ndarray of 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
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)))