Math#
This module contains math functions.
Core
|
Return a sequence from batches of evenly spaced samples between pairs of milestone points. |
|
Rotation matrix with given rotation axis and dimension (2d or 3d). |
Field
|
Return the deformation gradient tensors of the n-th field. |
|
Return Lagrangian strain tensor or its principal values of the n-th field. |
|
Compute the Seth-Hill strains. |
|
Extract gradient or interpolated field values at quadrature points. |
|
Return values of a field or a tuple of fields. |
|
Calculate the norm of an array or the norms of a list of arrays. |
|
Interpolate method of a field. |
|
Return the gradient of a field or the gradient of a basis-array. |
Tensor
|
Return identity matrices with ones on the diagonal of the first two axes and zeros elsewhere. |
|
Return the symmetric parts of second-order tensors. |
|
Return the dyadic product of two first-order or two second-order tensors. |
|
Return the in-plane components of 2d or 3d second order tensors, where the planes are defined by their standard unit vectors. |
|
Return the inverses of second-order tensors. |
|
Return the determinants of second order tensors. |
|
Return the deviatoric parts of second-order tensors. |
|
Return the cofactors of second-order tensors. |
|
Compute the eigenvalues and right eigenvectors of a square array. |
|
Return the eigenvalues and eigenvectors of a complex Hermitian (conjugate symmetric) or a real symmetric matrix. |
|
Eigenvalues (and optional principal shear values) of a matrix A. |
|
Eigenvalues (and optional principal shear values) of a symmetric matrix A. |
Return the Equivalent von Mises values of symmetric second order-tensors. |
|
|
Transpose (mode=1) or major-transpose (mode=2) of matrix A. |
|
|
|
Return the sums along the diagonals of second order tensors. |
|
Return the ik-crossed dyadic product of two second-order tensors. |
|
Return the il-crossed dyadic product of two second-order tensors. |
|
Return the crossed dyadic product of two second-order tensors. |
|
Return the cross product of two vectors. |
|
Dot-product of A and B with inputs of n trailing axes.. |
|
Return the double-dot products of first-, second-, third- or fourth-order tensors. |
|
Return a three-dimensional second-order tensor in reduced symmetric (Voigt- notation) vector/matrix storage. |
|
|
|
** Equation system**
|
Solve a linear equation system for two-dimensional unknowns with optional elementwise-operating trailing axes. |
Detailed API Reference
- felupe.math.cdya(A, B, parallel=False, out=None, **kwargs)[source]#
Return the crossed dyadic product of two second-order tensors.
- Parameters:
A (ndarray of shape (M, M, ...)) – Array with second-order tensors.
B (ndarray of shape (M, M, ...)) – Array with second-order tensors.
parallel (bool, optional) – A flag to enable a threaded evaluation of the results (default is False).
**kwargs (dict, optional) – Optional keyword-arguments for
numpy.einsum()
, e.g.out=None
.
- Returns:
The array of ik-crossed dyadic products.
- Return type:
ndarray of shape (M, M, M, M, …)
Notes
The first two axes are the tensor dimensions and all remaining trailing axes are treated as batch dimensions. The definition of the crossed dyadic product is given in Eq. (1).
(1)#\[ \begin{align}\begin{aligned}\mathbb{C} &= \boldsymbol{A} \odot \boldsymbol{B} = \frac{1}{2} \left( \boldsymbol{A} \overset{ik}{\otimes} \boldsymbol{B} + \boldsymbol{A} \overset{il}{\otimes} \boldsymbol{B} \right)\\\mathbb{C}_{ijkl} &= \frac{1}{2} \left( A_{ik}~B_{jl} + A_{il}~B_{kj} \right)\end{aligned}\end{align} \]Examples
>>> import felupe as fem >>> import numpy as np >>> >>> A = fem.math.transpose(np.arange(9, dtype=float).reshape(1, 3, 3).T) >>> B = fem.math.transpose(np.arange(100, 109, dtype=float).reshape(1, 3, 3).T) >>> C = fem.math.cdya(A, B) >>> C.shape (3, 3, 3, 3, 1)
>>> C[..., 0].reshape(9, 9) array([[ 0. , 50. , 100. , 50. , 102. , 154. , 100. , 154. , 208. ], [ 0. , 50.5, 101. , 51.5, 104. , 156.5, 103. , 157.5, 212. ], [ 0. , 51. , 102. , 53. , 106. , 159. , 106. , 161. , 216. ], [300. , 351.5, 403. , 354.5, 408. , 461.5, 409. , 464.5, 520. ], [306. , 358. , 410. , 362. , 416. , 470. , 418. , 474. , 530. ], [312. , 364.5, 417. , 369.5, 424. , 478.5, 427. , 483.5, 540. ], [600. , 653. , 706. , 659. , 714. , 769. , 718. , 775. , 832. ], [612. , 665.5, 719. , 672.5, 728. , 783.5, 733. , 790.5, 848. ], [624. , 678. , 732. , 686. , 742. , 798. , 748. , 806. , 864. ]])
See also
felupe.math.dya
Dyadic product of two first-order or two second-order tensors.
felupe.math.cdya_ik
ik-crossed dyadic product of two second-order tensors.
felupe.math.cdya_il
il-crossed dyadic product of two second-order tensors.
- felupe.math.cdya_ik(A, B, parallel=False, **kwargs)[source]#
Return the ik-crossed dyadic product of two second-order tensors.
- Parameters:
A (ndarray of shape (N, N, ...)) – Array with second-order tensors.
B (ndarray of shape (M, M, ...)) – Array with second-order tensors.
parallel (bool, optional) – A flag to enable a threaded evaluation of the results (default is False).
**kwargs (dict, optional) – Optional keyword-arguments for
numpy.einsum()
, e.g.out=None
.
- Returns:
The array of ik-crossed dyadic products.
- Return type:
ndarray of shape (N, M, N, M, …)
Notes
The first two axes are the tensor dimensions and all remaining trailing axes are treated as batch dimensions. The definition of the ik-crossed dyadic product is given in Eq. (2).
(2)#\[ \begin{align}\begin{aligned}\mathbb{C} &= \boldsymbol{A} \overset{ik}{\otimes} \boldsymbol{B}\\\mathbb{C}_{ijkl} &= A_{ik}\ B_{jl}\end{aligned}\end{align} \]Examples
>>> import felupe as fem >>> import numpy as np >>> >>> A = fem.math.transpose(np.arange(9, dtype=float).reshape(1, 3, 3).T) >>> B = fem.math.transpose(np.arange(100, 109, dtype=float).reshape(1, 3, 3).T) >>> C = fem.math.cdya_ik(A, B) >>> C.shape (3, 3, 3, 3, 1)
>>> C[..., 0].reshape(9, 9) array([[ 0., 0., 0., 100., 101., 102., 200., 202., 204.], [ 0., 0., 0., 103., 104., 105., 206., 208., 210.], [ 0., 0., 0., 106., 107., 108., 212., 214., 216.], [300., 303., 306., 400., 404., 408., 500., 505., 510.], [309., 312., 315., 412., 416., 420., 515., 520., 525.], [318., 321., 324., 424., 428., 432., 530., 535., 540.], [600., 606., 612., 700., 707., 714., 800., 808., 816.], [618., 624., 630., 721., 728., 735., 824., 832., 840.], [636., 642., 648., 742., 749., 756., 848., 856., 864.]])
See also
felupe.math.dya
Dyadic product of two first-order or two second-order tensors.
felupe.math.cdya
Crossed dyadic product of two second-order tensors.
felupe.math.cdya_il
il-crossed dyadic product of two second-order tensors.
- felupe.math.cdya_il(A, B, parallel=False, **kwargs)[source]#
Return the il-crossed dyadic product of two second-order tensors.
- Parameters:
A (ndarray of shape (N, N, ...)) – Array with second-order tensors.
B (ndarray of shape (M, M, ...)) – Array with second-order tensors.
parallel (bool, optional) – A flag to enable a threaded evaluation of the results (default is False).
**kwargs (dict, optional) – Optional keyword-arguments for
numpy.einsum()
, e.g.out=None
.
- Returns:
The array of ik-crossed dyadic products.
- Return type:
ndarray of shape (N, M, N, M, …)
Notes
The first two axes are the tensor dimensions and all remaining trailing axes are treated as batch dimensions. The definition of the il-crossed dyadic product is given in Eq. (3).
(3)#\[ \begin{align}\begin{aligned}\mathbb{C} &= \boldsymbol{A} \overset{il}{\otimes} \boldsymbol{B}\\\mathbb{C}_{ijkl} &= A_{il}\ B_{kj}\end{aligned}\end{align} \]Examples
>>> import felupe as fem >>> import numpy as np >>> >>> A = fem.math.transpose(np.arange(9, dtype=float).reshape(1, 3, 3).T) >>> B = fem.math.transpose(np.arange(100, 109, dtype=float).reshape(1, 3, 3).T) >>> C = fem.math.cdya_il(A, B) >>> C.shape (3, 3, 3, 3, 1)
>>> C[..., 0].reshape(9, 9) array([[ 0., 100., 200., 0., 103., 206., 0., 106., 212.], [ 0., 101., 202., 0., 104., 208., 0., 107., 214.], [ 0., 102., 204., 0., 105., 210., 0., 108., 216.], [300., 400., 500., 309., 412., 515., 318., 424., 530.], [303., 404., 505., 312., 416., 520., 321., 428., 535.], [306., 408., 510., 315., 420., 525., 324., 432., 540.], [600., 700., 800., 618., 721., 824., 636., 742., 848.], [606., 707., 808., 624., 728., 832., 642., 749., 856.], [612., 714., 816., 630., 735., 840., 648., 756., 864.]])
See also
felupe.math.dya
Dyadic product of two first-order or two second-order tensors.
felupe.math.cdya_ik
ik-crossed dyadic product of two second-order tensors.
felupe.math.cdya
Crossed dyadic product of two second-order tensors.
- felupe.math.cof(A, sym=False, out=None)[source]#
Return the cofactors of second-order tensors.
- Parameters:
A (ndarray of shape (1, 1, ...), (2, 2, ...) or (3, 3, ...)) – The array of second-order tensors.
sym (bool, optional) – A flag to assume symmetric second-order tensors. Only the upper-triangle elements are taken into account. Default is False.
out (ndarray or None, optional) – If provided, the calculation is done into this array.
- Returns:
The cofactors of second-order tensors.
- Return type:
ndarray of shape (1, 1, …), (2, 2, …) or (3, 3, …)
Notes
The first two axes are the matrix dimensions and all remaining trailing axes are treated as batch dimensions.
The inverse of a three-dimensional second-order tensor is obtained by Eq. (4)
(4)#\[ \begin{align}\begin{aligned}\text{cof}(\boldsymbol{A}) &= \det (\boldsymbol{A}) \boldsymbol{A}^{-T}\\\text{cof}(\boldsymbol{A}) &= \begin{bmatrix} \boldsymbol{A}_2 \times \boldsymbol{A}_3 & \boldsymbol{A}_3 \times \boldsymbol{A}_1 & \boldsymbol{A}_1 \times \boldsymbol{A}_2 \end{bmatrix}\end{aligned}\end{align} \]with the column (grid) vectors \(\boldsymbol{A}_j\), see Eq. (5).
(5)#\[\boldsymbol{A} = \begin{bmatrix} \boldsymbol{A}_1 & \boldsymbol{A}_2 & \boldsymbol{A}_3 \end{bmatrix}\]Examples
>>> import felupe as fem >>> import numpy as np >>> >>> A = np.array([1.0, 1.3, 1.5, 1.3, 1.1, 1.4, 1.5, 1.4, 1.2]).reshape(3, 3, 1) >>> A[..., 0] array([[1. , 1.3, 1.5], [1.3, 1.1, 1.4], [1.5, 1.4, 1.2]])
>>> cofA = fem.math.cof(A) >>> cofA[..., 0] array([[-0.64, 0.54, 0.17], [ 0.54, -1.05, 0.55], [ 0.17, 0.55, -0.59]])
>>> (fem.math.det(A) * fem.math.transpose(fem.math.inv(A)))[..., 0] array([[-0.64, 0.54, 0.17], [ 0.54, -1.05, 0.55], [ 0.17, 0.55, -0.59]])
See also
felupe.math.det
Return the determinants of second order tensors.
felupe.math.inv
Return the inverses of second order tensors.
- felupe.math.cross(a, b)[source]#
Return the cross product of two vectors.
- Parameters:
a (ndarray of shape (N, ...)) – First array of vectors.
b (ndarray of shape (N, ...)) – Second array of vectors.
- Returns:
c – Vector cross-products.
- Return type:
ndarray of shape (N, …)
Notes
The first axis is the vector dimension and all remaining trailing axes are treated as batch dimensions.
\[\boldsymbol{c} = \boldsymbol{a} \times \boldsymbol{b}\]Note
This is
numpy.cross()
withaxisa = axisb = axisc = 0
.Examples
>>> import felupe as fem >>> import numpy as np >>> >>> a = np.arange(30, dtype=float).reshape(3, 2, 5) >>> b = np.arange(80, 20, -2, dtype=float).reshape(5, 2, 3).T >>> fem.math.cross(a, b) array([[[-800., -682., -564., -446., -328.], [-750., -632., -514., -396., -278.]], [[1600., 1364., 1128., 892., 656.], [1500., 1264., 1028., 792., 556.]], [[-800., -682., -564., -446., -328.], [-750., -632., -514., -396., -278.]]])
See also
numpy.cross
Return the cross product of two (arrays of) vectors.
- felupe.math.ddot(A, B, mode=(2, 2), parallel=False, **kwargs)[source]#
Return the double-dot products of first-, second-, third- or fourth-order tensors.
- Parameters:
A (ndarray of shape (M, M, ...)) – Array with first-, second-, third- or fourth-order tensors.
B (ndarray of shape (M, M, ...)) – Array with first-, second-, third- or fourth-order tensors.
mode (tuple of int, optional) – Mode of operation. Return the double-dot products of two second-order tensors with (2, 2), of a second-order and a fourth-order tensor with (2, 4) and so on. Default is (2, 2).
parallel (bool, optional) – A flag to enable a threaded evaluation of the results (default is False).
**kwargs (dict, optional) – Optional keyword-arguments for
numpy.einsum()
, e.g.out=None
.
- Returns:
Array with the double-dot products.
- Return type:
ndarray of shape (…)
Notes
The first two axes are the tensor dimensions and all remaining trailing axes are treated as batch dimensions.
The double-dot product is obtained by Eq. (6) for two second-order tensors,
(6)#\[ \begin{align}\begin{aligned}c &= \boldsymbol{A} : \boldsymbol{B}\\c &= A_{ij} : B_{ij}\end{aligned}\end{align} \]by Eq. (7) for two fourth-order tensors,
(7)#\[ \begin{align}\begin{aligned}\mathbb{C} &= \mathbb{A} : \mathbb{B}\\\mathbb{C}_{ijmn} &= \mathbb{A}_{ijkl} : \mathbb{B}_{klmn}\end{aligned}\end{align} \]by Eq. (8) for a second-order and a third-order tensor,
(8)#\[ \begin{align}\begin{aligned}\boldsymbol{c} = \boldsymbol{A} : \mathcal{B}\\\qquad c_{k} &= A_{ij} : \mathcal{B}_{ijk}\end{aligned}\end{align} \]by Eq. (9) for a third-order and a second-order tensor,
(9)#\[ \begin{align}\begin{aligned}\boldsymbol{c} = \mathcal{A} : \boldsymbol{A}\\\qquad c_{i} &= \mathcal{A}_{ijk} : B_{jk}\end{aligned}\end{align} \]by Eq. (10) for a second-order and a fourth-order tensor
(10)#\[ \begin{align}\begin{aligned}\boldsymbol{C} &= \boldsymbol{A} : \mathbb{B}\\C_{kl} &= A_{ij} : \mathbb{B}_{ijkl}\end{aligned}\end{align} \]and by Eq. (11) for a fourth-order and a second-order tensor.
(11)#\[ \begin{align}\begin{aligned}\boldsymbol{C} &= \mathbb{A} : \boldsymbol{A}\\C_{ij} &= \mathbb{A}_{ijkl} : B_{kl}\end{aligned}\end{align} \]Examples
>>> import felupe as fem >>> import numpy as np >>> >>> A = fem.math.transpose(np.arange(9, dtype=float).reshape(1, 3, 3).T) >>> B = A + 10
>>> fem.math.ddot(A, B)[..., 0] array(564.)
>>> fem.math.ddot(fem.math.dya(A, A), B, mode=(4, 2))[..., 0] array([[ 0., 564., 1128.], [1692., 2256., 2820.], [3384., 3948., 4512.]])
- felupe.math.deformation_gradient(field, n=0)[source]#
Return the deformation gradient tensors of the n-th field.
- felupe.math.det(A, out=None)[source]#
Return the determinants of second order tensors.
- Parameters:
A (ndarray of shape (M, M, ...)) – The array of second-order tensors.
out (ndarray or None, optional) – If provided, the calculation is done into this array.
- Returns:
Array with the determinants.
- Return type:
ndarray of shape (…)
Notes
The first two axes are the tensor dimensions and all remaining trailing axes are treated as batch dimensions.
The determinant of a second-order tensor is obtained by Eq. (12).
(12)#\[ \begin{align}\begin{aligned}\det \left( \boldsymbol{A} \right) &= A_{11}~A_{22}~A_{33} + A_{12}~A_{23}~A_{31} + A_{13}~A_{21}~A_{32}\\ &- A_{31}~A_{22}~A_{13} - A_{32}~A_{23}~A_{11} - A_{33}~A_{21}~A_{12}\end{aligned}\end{align} \]Examples
>>> import felupe as fem >>> import numpy as np >>> >>> A = fem.math.transpose(np.arange(9, dtype=float).reshape(1, 3, 3).T) / 10 >>> A += np.eye(3).reshape(3, 3, 1)
>>> A.shape (3, 3, 1)
>>> A[..., 0] array([[1. , 0.1, 0.2], [0.3, 1.4, 0.5], [0.6, 0.7, 1.8]])
>>> fem.math.det(A)[..., 0] array(2.02)
- felupe.math.dev(A, out=None)[source]#
Return the deviatoric parts of second-order tensors.
- Parameters:
A (ndarray of shape (M, M, ...)) – The array of second-order tensors.
out (ndarray or None, optional) – If provided, the calculation is done into this array.
- Returns:
The deviatoric parts of second-order tensors.
- Return type:
ndarray of shape (M, M, …)
Notes
The first two axes are the tensor dimensions and all remaining trailing axes are treated as batch dimensions.
The deviatoric part of a three-dimensional second-order tensor is obtained by Eq. (13).
(13)#\[\text{dev} \left( \boldsymbol{A} \right) = \boldsymbol{A} - \frac{\text{tr}(\boldsymbol{A})}{3} \boldsymbol{1}\]Examples
>>> import felupe as fem >>> import numpy as np >>> >>> A = fem.math.transpose(np.arange(9, dtype=float).reshape(1, 3, 3).T) >>> A[..., 0] array([[0., 1., 2.], [3., 4., 5.], [6., 7., 8.]])
>>> fem.math.sym(A)[..., 0] array([[0., 2., 4.], [2., 4., 6.], [4., 6., 8.]])
See also
numpy.eye
Return a 2-D array with ones on the diagonal and zeros elsewhere.
- felupe.math.dot(A, B, mode=(2, 2), parallel=False, **kwargs)[source]#
Dot-product of A and B with inputs of n trailing axes..
- felupe.math.dya(A, B, mode=2, parallel=False, **kwargs)[source]#
Return the dyadic product of two first-order or two second-order tensors.
- Parameters:
A (ndarray of shape (N, ...) or (N, N, ...)) – Array with first-order or second-order tensors.
B (ndarray of shape (M, ...) or (M, M, ...)) – Array with first-order or second-order tensors.
mode (int, optional) – Mode of operation. Return the dyadic products of two second-order tensors with 2 and the dyadic products of two first-order tensors with 1. Default is 2.
parallel (bool, optional) – A flag to enable a threaded evaluation of the results (default is False).
**kwargs (dict, optional) – Optional keyword-arguments for
numpy.multiply()
, e.g.out=None
.
- Returns:
The array of dyadic products.
- Return type:
ndarray of shape (N, M, …) or (N, N, M, M, …)
Notes
The first two axes are the tensor dimensions and all remaining trailing axes are treated as batch dimensions. The definition of the dyadic product is given in Eq. (14) for two second-order tensors
(14)#\[ \begin{align}\begin{aligned}\mathbb{C} &= \boldsymbol{A} \otimes \boldsymbol{B}\\\mathbb{C}_{ijkl} &= A_{ij}\ B_{kl}\end{aligned}\end{align} \]and in Eq. (15) for two first-order tensors.
(15)#\[ \begin{align}\begin{aligned}\boldsymbol{C} &= \boldsymbol{a} \otimes \boldsymbol{b}\\C_{ij} &= a_{i} \otimes b_{j}\end{aligned}\end{align} \]Examples
>>> import felupe as fem >>> import numpy as np >>> >>> A = fem.math.transpose(np.arange(9, dtype=float).reshape(1, 3, 3).T) >>> B = fem.math.transpose(np.arange(100, 109, dtype=float).reshape(1, 3, 3).T) >>> C = fem.math.dya(A, B) >>> C.shape (3, 3, 3, 3, 1)
>>> C[..., 0].reshape(9, 9) array([[ 0., 0., 0., 0., 0., 0., 0., 0., 0.], [100., 101., 102., 103., 104., 105., 106., 107., 108.], [200., 202., 204., 206., 208., 210., 212., 214., 216.], [300., 303., 306., 309., 312., 315., 318., 321., 324.], [400., 404., 408., 412., 416., 420., 424., 428., 432.], [500., 505., 510., 515., 520., 525., 530., 535., 540.], [600., 606., 612., 618., 624., 630., 636., 642., 648.], [700., 707., 714., 721., 728., 735., 742., 749., 756.], [800., 808., 816., 824., 832., 840., 848., 856., 864.]])
See also
felupe.math.cdya
Crossed-dyadic product of two second-order tensors.
felupe.math.cdya_ik
ik-crossed dyadic product of two second-order tensors.
felupe.math.cdya_il
il-crossed dyadic product of two second-order tensors.
- felupe.math.eig(a, eig=<function eig>, **kwargs)[source]#
Compute the eigenvalues and right eigenvectors of a square array.
- Parameters:
a (ndarray of shape (M, M, ...)) – Matrices for which the eigenvalues and right eigenvectors will be computed.
eig (callable, optional) – A callable for the eigenvalue and eigenvector evaluation compatible with
numpy.linalg.eig()
(default isnumpy.linalg.eig()
).**kwargs (dict, optional) – Optional keyword-arguments are passed to
eig(a, **kwargs)
.
- Returns:
A namedtuple with the following attributes
eigenvalues (ndarray of shape (M, …)) – The eigenvalues, each repeated according to its multiplicity. The eigenvalues are not necessarily ordered. The resulting array will be of complex type, unless the imaginary part is zero in which case it will be cast to a real type. When a is real the resulting eigenvalues will be real (0 imaginary part) or occur in conjugate pairs.
eigenvectors (ndarray of shape (M, M, …)) – The normalized (unit “length”) eigenvectors, such that the column
eigenvectors[:, i]
is the eigenvector corresponding to the eigenvalueeigenvalues[i]
.
Notes
Note
The first two axes are the tensor dimensions and all remaining trailing axes are treated as batch dimensions.
Examples
>>> import numpy as np >>> import felupe as fem >>> >>> x = np.array([1, 0, 0.3, 0, 1.3, 0, 0, 0, 0.7]).reshape(3, 3) >>> y = np.array([1, 0, 0, 0, 0.4, 0.1, 0, 0, 1.6]).reshape(3, 3) >>> F = np.stack([x, y], axis=2) >>> >>> C = fem.math.dot(fem.math.transpose(F), F) >>> w, v = fem.math.eig(C) >>> >>> w[0] array([1.15619667, 2.57066372])
The associated eigenvectors are extracted. The first column is the eigenvector for the first right Cauchy-Green deformation tensor and the second column for the second right Cauchy-Green deformation tensor.
>>> v[:, 0] array([[0.88697868, 0. ], [0. , 0.01659066], [0.46181038, 0.99986237]])
See also
felupe.math.eigh
Return the eigenvalues and eigenvectors of a complex Hermitian (conjugate symmetric) or a real symmetric matrix.
numpy.linalg.eig
Compute the eigenvalues and right eigenvectors of a square array.
numpy.linalg.eigh
Return the eigenvalues and eigenvectors of a complex Hermitian (conjugate symmetric) or a real symmetric matrix.
- felupe.math.eigh(a, UPLO='L')[source]#
Return the eigenvalues and eigenvectors of a complex Hermitian (conjugate symmetric) or a real symmetric matrix.
Returns two objects, a 1-D array containing the eigenvalues of a, and a 2-D square array or matrix (depending on the input type) of the corresponding eigenvectors (in columns).
- Parameters:
a (ndarray of shape (M, M, ...)) – Matrices for which the eigenvalues and right eigenvectors will be computed.
UPLO ({"L", "U"}, optional) – Specifies whether the calculation is done with the lower triangular part of a (‘L’, default) or the upper triangular part (‘U’). Irrespective of this value only the real parts of the diagonal will be considered in the computation to preserve the notion of a Hermitian matrix. It therefore follows that the imaginary part of the diagonal will always be treated as zero.
- Returns:
A namedtuple with the following attributes
eigenvalues ((M, …) ndarray) – The eigenvalues in ascending order, each repeated according to its multiplicity.
eigenvectors ((M, M, …) ndarray) – The normalized (unit “length”) eigenvectors, such that the column
eigenvectors[:, i]
is the eigenvector corresponding to the eigenvalueeigenvalues[i]
.
Notes
Note
The first two axes are the tensor dimensions and all remaining trailing axes are treated as batch dimensions.
Examples
>>> import numpy as np >>> import felupe as fem >>> >>> x = np.array([1, 0, 0.3, 0, 1.3, 0, 0, 0, 0.7]).reshape(3, 3) >>> y = np.array([1, 0, 0, 0, 0.4, 0.1, 0, 0, 1.6]).reshape(3, 3) >>> F = np.stack([x, y], axis=2) >>> >>> C = fem.math.dot(fem.math.transpose(F), F) >>> w, v = fem.math.eigh(C) >>> >>> w[-1] array([1.69 , 2.57066372])
The associated eigenvectors are extracted. The first column is the eigenvector, associated to the greated eigenvalue, for the first right Cauchy-Green deformation tensor and the second column for the second right Cauchy-Green deformation tensor.
>>> v[:, 0] array([[-4.61810381e-01, 0.00000000e+00], [ 1.11022302e-16, -9.99862366e-01], [ 8.86978676e-01, 1.65906569e-02]])
See also
felupe.math.eig
Compute the eigenvalues and right eigenvectors of a square array.
numpy.linalg.eig
Compute the eigenvalues and right eigenvectors of a square array.
numpy.linalg.eigh
Return the eigenvalues and eigenvectors of a complex Hermitian (conjugate symmetric) or a real symmetric matrix.
- felupe.math.eigvals(a, shear=False, eigvals=<function eigvals>, **kwargs)[source]#
Eigenvalues (and optional principal shear values) of a matrix A.
- felupe.math.eigvalsh(A, shear=False)[source]#
Eigenvalues (and optional principal shear values) of a symmetric matrix A.
- felupe.math.equivalent_von_mises(A)[source]#
Return the Equivalent von Mises values of symmetric second order-tensors.
- Parameters:
A (ndarray of shape (2, 2, ...) or (3, 3, ...)) – Symmetric second-order tensors for which the equivalent von Mises values will be computed.
- Returns:
The equivalent von Mises values.
- Return type:
ndarray of shape (…)
Notes
The first two axes are the tensor dimensions and all remaining trailing axes are treated as batch dimensions.
The equivalent von Mises value of a three-dimensional symmetric second order-tensor is given in Eq. (16).
(16)#\[\boldsymbol{A}_{v} = \sqrt{ \frac{3}{2} \text{dev}(\boldsymbol{A}) : \text{dev}(\boldsymbol{A}) }\]Examples
>>> import numpy as np >>> import felupe as fem >>> >>> stress = np.array([1, 1.1, 1.2, 1.1, 1.3, 1.4, 1.2, 1.4, 1.5]).reshape(3, 3, 1) >>> fem.math.equivalent_von_mises(stress) array([3.74432905])
>>> stress = np.diag([3, 1, 1]).reshape(3, 3, 1) >>> fem.math.equivalent_von_mises(stress) array([2.])
- felupe.math.extract(field, grad=True, sym=False, add_identity=True)[source]#
Extract gradient or interpolated field values at quadrature points.
- felupe.math.grad(x, **kwargs)[source]#
Return the gradient of a field or the gradient of a basis-array.
- felupe.math.identity(A=None, dim=None, shape=None, dtype=None)[source]#
Return identity matrices with ones on the diagonal of the first two axes and zeros elsewhere.
- Parameters:
A (ndarray of shape (N, M, ...) or None, optional) – The array of input matrices. If provided, the dimension and the shape of the trailing axes are taken from this array. If None,
dim
andshape
are required. Default is None.dim (int or None, optional) – The dimension of the matrix axes. If None, it is taken from
A
(default is None).shape (tuple of int or None, optional) – A tuple containing the shape of the trailing axes (batch dimensions). Default is None.
dtype (data-type or None, optional) – Data-type of the returned array. If None and
A
is not None, the data-type ofA
is used. If None andA
is None, the data-type of the returned array isfloat
.
- Returns:
The identity matrix.
- Return type:
ndarray of shape (N, M, *np.ones_like(…)) or (dim, dim, *np.ones_like(…))
Notes
The first two axes are the matrix dimensions and all remaining trailing axes are treated as batch dimensions. As the identity matrix is idependent of the input matrices
A
, the size of the trailing axes is reduced to one.The identity matrix is defined by Eq. (17)
(17)#\[\boldsymbol{1} = \boldsymbol{A} \boldsymbol{A}^{-1} = \boldsymbol{A}^{-1} \boldsymbol{A}\]and it is shown in Eq. (18).
(18)#\[\begin{split}\boldsymbol{1} = \begin{bmatrix} 1 & 0 & 0 \\ 0 & 1 & 0 \\ 0 & 0 & 1 \end{bmatrix}\end{split}\]Examples
>>> import felupe as fem >>> import numpy as np >>> >>> A = np.random.rand(3, 2, 8, 20) >>> I = fem.math.identity(A) >>> I.shape (3, 2, 1, 1)
With given dimension of the matrix axes the shape of the output is different.
>>> fem.math.identity(A, dim=2).shape (2, 2, 1, 1)
Note how the number of batch axes change if a
shape
is given.>>> fem.math.identity(A, shape=(4, 7, 3)).shape (3, 2, 1, 1, 1)
See also
numpy.eye
Return a 2-D array with ones on the diagonal and zeros elsewhere.
- felupe.math.inplane(A, vectors, **kwargs)[source]#
Return the in-plane components of 2d or 3d second order tensors, where the planes are defined by their standard unit vectors.
- Parameters:
A (ndarray of shape (N, N, ...)) – Second-order tensors.
vectors (list of ndarray of shape (N, ...)) – List of standard unit vectors of the planes.
**kwargs (dict, optional) – Optional keyword-arguments for
numpy.einsum()
.
Notes
Note
Out-of-plane tensor components are ignored.
The first two axes of the tensor and the first axis of each vector are the tensor/vector dimensions and all remaining trailing axes are treated as batch dimensions.
\[\boldsymbol{A}_{N-1} = \sum_{\alpha, \beta} \left( \boldsymbol{E}_\alpha^T \boldsymbol{A} \boldsymbol{E}_\beta\right) \boldsymbol{E}_\alpha \otimes \boldsymbol{E}_\beta\]- Returns:
In-plane components.
- Return type:
ndarray of shape (N - 1, N - 1, …)
Examples
>>> import felupe as fem >>> >>> A = np.array([[0, 3, 5], [3, 1, 4], [5, 4, 2]]) >>> vectors = [[1, 0, 0], [0, 1, 0]] >>> >>> A_inplane = fem.math.inplane(A, vectors) >>> A_inplane array([[0, 3], [3, 1]])
- felupe.math.inv(A, determinant=None, full_output=False, sym=False, out=None)[source]#
Return the inverses of second-order tensors.
- Parameters:
A (ndarray of shape (1, 1, ...), (2, 2, ...) or (3, 3, ...)) – The array of second-order tensors.
determinant (ndarray or None, optional) – The array with the pre-evaluated determinants of the second-order tensors ( default is None).
full_output (bool, optional) – A flag to return the array of inverses and the determinants (default is False).
sym (bool, optional) – A flag to assume symmetric second-order tensors. Only the upper-triangle elements are taken into account. Default is False.
out (ndarray or None, optional) – If provided, the calculation is done into this array.
- Returns:
The inverses of second-order tensors.
- Return type:
ndarray of shape (1, 1, …), (2, 2, …) or (3, 3, …)
Notes
The first two axes are the tensor dimensions and all remaining trailing axes are treated as batch dimensions.
The inverse of a three-dimensional second-order tensor is obtained by Eq. (19) and Eq. (20)
(19)#\[ \begin{align}\begin{aligned}\boldsymbol{A}^{-1} \boldsymbol{A} &= \boldsymbol{A} \boldsymbol{A}^{-1} = \boldsymbol{1}\\\boldsymbol{A}^{-1} &= \frac{1}{\det(\boldsymbol{A})} \text{cof}(\boldsymbol{A})^T\end{aligned}\end{align} \](20)#\[\begin{split}\boldsymbol{A}^{-1} = \frac{1}{\det(\boldsymbol{A})} \begin{bmatrix} \left( \boldsymbol{A}_2 \times \boldsymbol{A}_3 \right)^T \\ \left( \boldsymbol{A}_3 \times \boldsymbol{A}_1 \right)^T \\ \left( \boldsymbol{A}_1 \times \boldsymbol{A}_2 \right)^T \end{bmatrix}\end{split}\]with the column (grid) vectors \(\boldsymbol{A}_j\), see Eq. (21).
(21)#\[\boldsymbol{A} = \begin{bmatrix} \boldsymbol{A}_1 & \boldsymbol{A}_2 & \boldsymbol{A}_3 \end{bmatrix}\]Examples
>>> import felupe as fem >>> import numpy as np >>> >>> A = np.array([1.0, 1.3, 1.5, 1.3, 1.1, 1.4, 1.5, 1.4, 1.2]).reshape(3, 3, 1) >>> A[..., 0] array([[1. , 1.3, 1.5], [1.3, 1.1, 1.4], [1.5, 1.4, 1.2]])
>>> invA = fem.math.inv(A) >>> invA[..., 0] array([[-2.01892744, 1.70347003, 0.5362776 ], [ 1.70347003, -3.31230284, 1.73501577], [ 0.5362776 , 1.73501577, -1.86119874]])
>>> fem.math.dot(A, invA)[..., 0].round(3) array([[ 1., -0., 0.], [-0., 1., 0.], [-0., -0., 1.]])
See also
felupe.math.det
Return the determinants of second order tensors.
felupe.math.cof
Return the cofactors of second order tensors.
- felupe.math.linsteps(points, num=10, endpoint=True, axis=None, axes=None, values=0.0)[source]#
Return a sequence from batches of evenly spaced samples between pairs of milestone points.
- Parameters:
points (array_like) – The milestone values of the sequence.
num (int, optional) – Number of samples (without the end point) to generate. Must be non-negative (default is 10).
endpoint (bool, optional) – If True,
points[-1]
is the last sample. Otherwise, it is not included (default is True).axis (int or None, optional) – The axis in the result to store the samples. By default (None), the samples will be concatenated along the existing first axis. If
axis > 0
, the samples will be inserted at columnaxis
. Only positive integers are supported.axes (int or None, optional) – Number of output columns if
axis
is not None. Requiresaxes > axis
and will be set toaxes = axis + 1
by default (None).values (array_like, optional) – Default values if
axis
is not None (not used for columnaxis
). Default is 0.0.
- Returns:
samples – Concatenated batches of
num
equally spaced samples.- Return type:
ndarray
Examples
>>> import felupe as fem >>> >>> fem.math.linsteps([0, 0.5, 1.5, 3.5], num=2) array([0. , 0.25, 0.5 , 1. , 1.5 , 2.5 , 3.5 ])
Not including the end value of the sequence does not change the step size.
>>> fem.math.linsteps([0, 0.5, 1.5, 3.5], num=2, endpoint=False) array([0. , 0.25, 0.5 , 1. , 1.5 , 2.5 ])
>>> fem.math.linsteps([0, 0.5, 1.5, 3.5], num=2, axis=1, values=[-1, 0]) array([[-1. , 0. ], [-1. , 0.25], [-1. , 0.5 ], [-1. , 1. ], [-1. , 1.5 ], [-1. , 2.5 ], [-1. , 3.5 ]])
Output with four columns:
>>> fem.math.linsteps([0, 0.5, 1.5, 3.5], num=2, axis=1, axes=4) array([[0. , 0. , 0. , 0. ], [0. , 0.25, 0. , 0. ], [0. , 0.5 , 0. , 0. ], [0. , 1. , 0. , 0. ], [0. , 1.5 , 0. , 0. ], [0. , 2.5 , 0. , 0. ], [0. , 3.5 , 0. , 0. ]])
The ouput degenerates to the end value
>>> fem.math.linsteps([0, 0.5, 1.5, 3.5], num=0) array([3.5])
or an empty array.
>>> fem.math.linsteps([0, 0.5, 1.5, 3.5], num=0, endpoint=False) array([], dtype=float64)
See also
numpy.linspace
Return evenly spaced numbers over a specified interval.
- felupe.math.norm(array, axis=None)[source]#
Calculate the norm of an array or the norms of a list of arrays.
- felupe.math.rotation_matrix(alpha_deg, dim=3, axis=0)[source]#
Rotation matrix with given rotation axis and dimension (2d or 3d).
- Parameters:
- Returns:
rotation_matrix – Rotation matrix of dim 2 or 3 with given rotation axis.
- Return type:
ndarray
Notes
The two-dimensional rotation axis is denoted in Eq. (22).
(22)#\[\begin{split}\boldsymbol{R}(\alpha) = \begin{bmatrix} \cos(\alpha) & -\sin(\alpha) \\ \sin(\alpha) & \cos(\alpha) \end{bmatrix}\end{split}\]A three-dimensional rotation matrix is created by inserting zeros in the row and column at the given axis of rotation and one at the intersection, see Eq. (23). If the axis of rotation is the second axis, the two- dimensinal rotation matrix is transposed.
(23)#\[\begin{split}\boldsymbol{R}(\alpha) = \begin{bmatrix} \cos(\alpha) & -\sin(\alpha) & 0 \\ \sin(\alpha) & \cos(\alpha) & 0 \\ 0 & 0 & 1 \end{bmatrix}\end{split}\]Examples
>>> import numpy as np >>> import felupe as fem >>> >>> R = fem.math.rotation_matrix(alpha_deg=45, dim=2) >>> x = np.array([1., 0.]) >>> y = R @ x >>> y array([0.70710678, 0.70710678])
- felupe.math.solve_2d(A, b, solve=<function solve>, **kwargs)[source]#
Solve a linear equation system for two-dimensional unknowns with optional elementwise-operating trailing axes.
- Parameters:
A (ndarray of shape (M, N, M, N, ...)) – The left-hand side array of the equation system.
b (ndarray of shape (M, N, ...)) – The right-hand side array of the equation system.
n (int, optional) – The dimension of the unknowns (default is 1).
solve (callable, optional) – A function with a compatible signature to
numpy.linalg.solve()
(default isnumpy.linalg.solve()
).**kwargs (dict, optional) – Optional keyword arguments are passed to the callable
solve(A, b, **kwargs)
.
- Returns:
x – The solved array of unknowns.
- Return type:
ndarray of shape (M, N, …)
Notes
The first two axes of the rhs
b
and the first four axes of the lhsA
are the tensor dimensions and all remaining trailing axes are treated as batch dimensions. This function finds \(x_{kl}\) for Eq. (24).(24)#\[A_{ijkl} : x_{kl} = b_{ij}\]Examples
>>> import numpy as np >>> import felupe as fem >>> >>> np.random.seed(855436) >>> >>> A = np.random.rand(3, 3, 3, 3, 2, 4) >>> b = np.ones((3, 3, 2, 4)) >>> >>> x = fem.math.solve_2d(A, b)
>>> x[..., 0, 0] array([[ 1.1442917 , 2.14516919, 2.00237954], [-0.69463749, -2.46685827, -7.21630899], [ 4.44825615, 1.35899745, 1.08645703]])
>>> x.shape (3, 3, 2, 4)
- felupe.math.strain(field, C=None, fun=<function strain_stretch_1d>, tensor=True, asvoigt=False, n=0, **kwargs)[source]#
Return Lagrangian strain tensor or its principal values of the n-th field.
- Parameters:
field (FieldContainer or None) – A field container with a displacement field from which the deformation gradient tensor is extracted if
def_grad
is None.C (ndarray of shape (N, N, ...) or None, optional) – Optional array of right Cauchy-Green deformation tensors. If None, the right Cauchy-Green deformation tensor is obtained from the field. Default is None.
fun (callable, optional) – A callable for the one-dimensional strain-stretch relation. Function signature must be
lambda stretch, **kwargs: strain
(default is the log. strain,strain_stretch_1d()
withk=0
).tensor (bool, optional) – Assemble and return the strain tensors if True or return their principal values only if False. Default is True.
asvoigt (bool, optional) – Return the symmetric strain tensors in reduced vector (Voigt) storage. Default is False.
n (int, optional) – The index of the displacement field (default is 0).
**kwargs (dict, optional) – Optional keyword-arguments are passed to
fun(stretch, **kwargs)
.
- Returns:
The strain tensors or their principal values.
- Return type:
ndarray of shape (N, N, …) tensor, (N!, …) asvoigt or (N, …) princ. values
Notes
The right Cauchy-Green deformation tensor is defined by the dot-products of the column vectors of the deformation gradient tensor and enables a quadratic length measure, see Eq. (25).
(25)#\[\boldsymbol{C} = \boldsymbol{F}^T \boldsymbol{F}\]The principal stretches are obtained as the square-roots of the eigenvalues of the right Cauchy-Green deformation tensor, see Eq. (26).
(26)#\[\lambda^2_\alpha, \boldsymbol{N}_\alpha = \text{eig}\left( \boldsymbol{C} \right)\]The Lagrangian strain tensor is assembled with the one-dimensional strain-stretch relation from Eq. (27) and the eigenbases as dyadic products of the eigenvectors, see Eq. (28).
(27)#\[E_\alpha = \text{f}\left( \lambda_\alpha \right)\](28)#\[\boldsymbol{E} = \sum_\alpha E_\alpha \ \boldsymbol{N}_\alpha \otimes \boldsymbol{N}_\alpha\]See also
math.strain_stretch_1d
Compute the Seth-Hill strains.
- felupe.math.strain_stretch_1d(stretch, k=0)[source]#
Compute the Seth-Hill strains.
- Parameters:
stretch (ndarray of shape (M, ...)) – The array of stretches.
k (float, optional) – The strain-exponent of the Seth-Hill strain-stretch relation (default is 0). If 0, the logarithmic strains are computed.
- Returns:
The computed strains.
- Return type:
ndarray of shape (M, …)
Notes
\[E(\lambda, k) = \frac{1}{k} \left( \lambda^k - 1 \right)\]
- felupe.math.sym(A, out=None)[source]#
Return the symmetric parts of second-order tensors.
- Parameters:
A (ndarray of shape (M, M, ...)) – The array of second-order tensors.
out (ndarray or None, optional) – If provided, the calculation is done into this array.
- Returns:
Array with the symmetric parts of the second-order tensors.
- Return type:
ndarray of shape (M, M, …)
Notes
The first two axes are the tensor dimensions and all remaining trailing axes are treated as batch dimensions.
The symmetric part of a second-order tensor is obtained by Eq. (29).
(29)#\[\text{sym} \left( \boldsymbol{A} \right) = \frac{1}{2} \left( \boldsymbol{A} + \boldsymbol{A}^T \right)\]Examples
>>> import felupe as fem >>> import numpy as np >>> >>> A = fem.math.transpose(np.arange(18, dtype=float).reshape(2, 3, 3).T) >>> A[..., 0] array([[0., 1., 2.], [3., 4., 5.], [6., 7., 8.]])
>>> fem.math.sym(A)[..., 0] array([[0., 2., 4.], [2., 4., 6.], [4., 6., 8.]])
- felupe.math.tovoigt(A, strain=False)[source]#
Return a three-dimensional second-order tensor in reduced symmetric (Voigt- notation) vector/matrix storage.
- Parameters:
A (ndarray of shape (2, 2, ...) or (3, 3, ...)) – Array with three-dimensional second-order tensors.
strain (bool, optional) – A flag to double the off-diagonal (shear) values for strain-like tensors. Default is False.
- Returns:
A three-dimensional second-order tensor in reduced symmetric (Voigt) vector/ matrix storage.
- Return type:
ndarray of shape (3, …) or (6, …)
Notes
The first two axes are the tensor dimensions and all remaining trailing axes are treated as batch dimensions.
For a symmetric three-dimensional second-order tensor \(C_{ij} = C_{ji}\), the upper triangle entries are inserted into a 6x1 vector, starting from the main diagonal, followed by the consecutive next upper diagonals.
\[\begin{split}\boldsymbol{C} = \begin{bmatrix} C_{11} & C_{12} & C_{13} \\ C_{12} & C_{22} & C_{23} \\ C_{13} & C_{23} & C_{33} \end{bmatrix} \qquad \longrightarrow \boldsymbol{C} = \begin{bmatrix} C_{11} \\ C_{22} \\ C_{33} \\ C_{12} \\ C_{23} \\ C_{13} \end{bmatrix}\end{split}\]Examples
>>> import felupe as fem >>> import numpy as np >>> >>> C = np.array([1.0, 1.3, 1.5, 1.3, 1.1, 1.4, 1.5, 1.4, 1.2]).reshape(3, 3, 1) >>> C[..., 0] array([[1. , 1.3, 1.5], [1.3, 1.1, 1.4], [1.5, 1.4, 1.2]])
>>> fem.math.tovoigt(C)[..., 0] array([1. , 1.1, 1.2, 1.3, 1.4, 1.5])
- felupe.math.trace(A, out=None)[source]#
Return the sums along the diagonals of second order tensors.
- Parameters:
A (ndarray of shape (M, M, ...)) – The array of second-order tensors.
out (ndarray or None, optional) – If provided, the calculation is done into this array.
- Returns:
Array with the sums along the diagonals.
- Return type:
ndarray of shape (…)
Notes
The first two axes are the tensor dimensions and all remaining trailing axes are treated as batch dimensions.
The trace of a second-order tensor is obtained by Eq. (30).
(30)#\[ \begin{align}\begin{aligned}\text{tr} \left( \boldsymbol{A} \right) &= \boldsymbol{A} : \boldsymbol{1}\\A_{ii} &= A_{ij} : \delta_{ij}\end{aligned}\end{align} \]Examples
>>> import felupe as fem >>> import numpy as np >>> >>> A = fem.math.transpose(np.arange(9, dtype=float).reshape(1, 3, 3).T) >>> A.shape (3, 3, 1)
>>> A[..., 0] array([[0., 1., 2.], [3., 4., 5.], [6., 7., 8.]])
>>> fem.math.trace(A)[..., 0] array(12.)
See also
numpy.trace
Return the sum along diagonals of the array.