Reference documentation for deal.II version 9.2.0
|
#include <deal.II/base/config.h>
#include <deal.II/base/numbers.h>
#include <deal.II/base/table_indices.h>
#include <deal.II/base/template_constraints.h>
#include <deal.II/base/tensor.h>
#include <algorithm>
#include <array>
#include <functional>
Go to the source code of this file.
Namespaces | |
internal | |
internal::SymmetricTensorImplementation | |
internal::SymmetricTensorAccessors | |
Enumerations | |
enum | SymmetricTensorEigenvectorMethod { SymmetricTensorEigenvectorMethod::hybrid, SymmetricTensorEigenvectorMethod::ql_implicit_shifts, SymmetricTensorEigenvectorMethod::jacobi } |
Functions | |
template<int dim, typename Number > | |
constexpr SymmetricTensor< 2, dim, Number > | unit_symmetric_tensor () |
template<int dim, typename Number > | |
constexpr SymmetricTensor< 4, dim, Number > | deviator_tensor () |
template<int dim, typename Number > | |
constexpr SymmetricTensor< 4, dim, Number > | identity_tensor () |
template<int dim, typename Number > | |
constexpr SymmetricTensor< 2, dim, Number > | invert (const SymmetricTensor< 2, dim, Number > &) |
template<int dim, typename Number > | |
constexpr SymmetricTensor< 4, dim, Number > | invert (const SymmetricTensor< 4, dim, Number > &) |
template<int dim2, typename Number > | |
constexpr Number | trace (const SymmetricTensor< 2, dim2, Number > &) |
template<int dim, typename Number > | |
constexpr SymmetricTensor< 2, dim, Number > | deviator (const SymmetricTensor< 2, dim, Number > &) |
template<int dim, typename Number > | |
constexpr Number | determinant (const SymmetricTensor< 2, dim, Number > &) |
constexpr TableIndices< 2 > | internal::SymmetricTensorAccessors::merge (const TableIndices< 2 > &previous_indices, const unsigned int new_index, const unsigned int position) |
constexpr TableIndices< 4 > | internal::SymmetricTensorAccessors::merge (const TableIndices< 4 > &previous_indices, const unsigned int new_index, const unsigned int position) |
template<int rank_, int dim, typename Number , typename OtherNumber > | |
constexpr SymmetricTensor< rank_, dim, typename ProductType< Number, OtherNumber >::type > | operator+ (const SymmetricTensor< rank_, dim, Number > &left, const SymmetricTensor< rank_, dim, OtherNumber > &right) |
template<int rank_, int dim, typename Number , typename OtherNumber > | |
constexpr SymmetricTensor< rank_, dim, typename ProductType< Number, OtherNumber >::type > | operator- (const SymmetricTensor< rank_, dim, Number > &left, const SymmetricTensor< rank_, dim, OtherNumber > &right) |
template<int rank_, int dim, typename Number , typename OtherNumber > | |
constexpr Tensor< rank_, dim, typename ProductType< Number, OtherNumber >::type > | operator+ (const SymmetricTensor< rank_, dim, Number > &left, const Tensor< rank_, dim, OtherNumber > &right) |
template<int rank_, int dim, typename Number , typename OtherNumber > | |
constexpr Tensor< rank_, dim, typename ProductType< Number, OtherNumber >::type > | operator+ (const Tensor< rank_, dim, Number > &left, const SymmetricTensor< rank_, dim, OtherNumber > &right) |
template<int rank_, int dim, typename Number , typename OtherNumber > | |
constexpr Tensor< rank_, dim, typename ProductType< Number, OtherNumber >::type > | operator- (const SymmetricTensor< rank_, dim, Number > &left, const Tensor< rank_, dim, OtherNumber > &right) |
template<int rank_, int dim, typename Number , typename OtherNumber > | |
constexpr Tensor< rank_, dim, typename ProductType< Number, OtherNumber >::type > | operator- (const Tensor< rank_, dim, Number > &left, const SymmetricTensor< rank_, dim, OtherNumber > &right) |
template<int dim, typename Number > | |
constexpr Number | third_invariant (const SymmetricTensor< 2, dim, Number > &t) |
template<int dim, typename Number > | |
constexpr Number | trace (const SymmetricTensor< 2, dim, Number > &d) |
template<int dim, typename Number > | |
constexpr Number | first_invariant (const SymmetricTensor< 2, dim, Number > &t) |
template<typename Number > | |
constexpr Number | second_invariant (const SymmetricTensor< 2, 1, Number > &) |
template<typename Number > | |
constexpr Number | second_invariant (const SymmetricTensor< 2, 2, Number > &t) |
template<typename Number > | |
constexpr Number | second_invariant (const SymmetricTensor< 2, 3, Number > &t) |
template<typename Number > | |
std::array< Number, 1 > | eigenvalues (const SymmetricTensor< 2, 1, Number > &T) |
template<typename Number > | |
std::array< Number, 2 > | eigenvalues (const SymmetricTensor< 2, 2, Number > &T) |
template<typename Number > | |
std::array< Number, 3 > | eigenvalues (const SymmetricTensor< 2, 3, Number > &T) |
template<int dim, typename Number > | |
void | internal::SymmetricTensorImplementation::tridiagonalize (const ::SymmetricTensor< 2, dim, Number > &A, ::Tensor< 2, dim, Number > &Q, std::array< Number, dim > &d, std::array< Number, dim - 1 > &e) |
template<int dim, typename Number > | |
std::array< std::pair< Number, Tensor< 1, dim, Number > >, dim > | internal::SymmetricTensorImplementation::ql_implicit_shifts (const ::SymmetricTensor< 2, dim, Number > &A) |
template<int dim, typename Number > | |
std::array< std::pair< Number, Tensor< 1, dim, Number > >, dim > | internal::SymmetricTensorImplementation::jacobi (::SymmetricTensor< 2, dim, Number > A) |
template<typename Number > | |
std::array< std::pair< Number, Tensor< 1, 2, Number > >, 2 > | internal::SymmetricTensorImplementation::hybrid (const ::SymmetricTensor< 2, 2, Number > &A) |
template<typename Number > | |
std::array< std::pair< Number, Tensor< 1, 3, Number > >, 3 > | internal::SymmetricTensorImplementation::hybrid (const ::SymmetricTensor< 2, 3, Number > &A) |
template<int dim, typename Number > | |
std::array< std::pair< Number, Tensor< 1, dim, Number > >, std::integral_constant< int, dim >::value > | eigenvectors (const SymmetricTensor< 2, dim, Number > &T, const SymmetricTensorEigenvectorMethod method=SymmetricTensorEigenvectorMethod::ql_implicit_shifts) |
template<int rank_, int dim, typename Number > | |
constexpr SymmetricTensor< rank_, dim, Number > | transpose (const SymmetricTensor< rank_, dim, Number > &t) |
template<int dim> | |
constexpr SymmetricTensor< 2, dim > | unit_symmetric_tensor () |
template<int dim> | |
constexpr SymmetricTensor< 4, dim > | deviator_tensor () |
template<int dim> | |
constexpr SymmetricTensor< 4, dim > | identity_tensor () |
template<int dim, typename Number > | |
constexpr SymmetricTensor< 4, dim, Number > | outer_product (const SymmetricTensor< 2, dim, Number > &t1, const SymmetricTensor< 2, dim, Number > &t2) |
template<int dim, typename Number > | |
constexpr SymmetricTensor< 2, dim, Number > | symmetrize (const Tensor< 2, dim, Number > &t) |
template<int rank_, int dim, typename Number > | |
constexpr SymmetricTensor< rank_, dim, Number > | operator* (const SymmetricTensor< rank_, dim, Number > &t, const Number &factor) |
template<int rank_, int dim, typename Number > | |
constexpr SymmetricTensor< rank_, dim, Number > | operator* (const Number &factor, const SymmetricTensor< rank_, dim, Number > &t) |
template<int rank_, int dim, typename Number , typename OtherNumber > | |
constexpr SymmetricTensor< rank_, dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > | operator* (const SymmetricTensor< rank_, dim, Number > &t, const OtherNumber &factor) |
template<int rank_, int dim, typename Number , typename OtherNumber > | |
constexpr SymmetricTensor< rank_, dim, typename ProductType< OtherNumber, typename EnableIfScalar< Number >::type >::type > | operator* (const Number &factor, const SymmetricTensor< rank_, dim, OtherNumber > &t) |
template<int rank_, int dim, typename Number , typename OtherNumber > | |
constexpr SymmetricTensor< rank_, dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > | operator/ (const SymmetricTensor< rank_, dim, Number > &t, const OtherNumber &factor) |
template<int rank_, int dim> | |
constexpr SymmetricTensor< rank_, dim > | operator* (const SymmetricTensor< rank_, dim > &t, const double factor) |
template<int rank_, int dim> | |
constexpr SymmetricTensor< rank_, dim > | operator* (const double factor, const SymmetricTensor< rank_, dim > &t) |
template<int rank_, int dim> | |
constexpr SymmetricTensor< rank_, dim > | operator/ (const SymmetricTensor< rank_, dim > &t, const double factor) |
template<int dim, typename Number , typename OtherNumber > | |
constexpr ProductType< Number, OtherNumber >::type | scalar_product (const SymmetricTensor< 2, dim, Number > &t1, const SymmetricTensor< 2, dim, OtherNumber > &t2) |
template<int dim, typename Number , typename OtherNumber > | |
constexpr ProductType< Number, OtherNumber >::type | scalar_product (const SymmetricTensor< 2, dim, Number > &t1, const Tensor< 2, dim, OtherNumber > &t2) |
template<int dim, typename Number , typename OtherNumber > | |
constexpr ProductType< Number, OtherNumber >::type | scalar_product (const Tensor< 2, dim, Number > &t1, const SymmetricTensor< 2, dim, OtherNumber > &t2) |
template<typename Number , typename OtherNumber > | |
constexpr void | double_contract (SymmetricTensor< 2, 1, typename ProductType< Number, OtherNumber >::type > &tmp, const SymmetricTensor< 4, 1, Number > &t, const SymmetricTensor< 2, 1, OtherNumber > &s) |
template<typename Number , typename OtherNumber > | |
constexpr void | double_contract (SymmetricTensor< 2, 1, typename ProductType< Number, OtherNumber >::type > &tmp, const SymmetricTensor< 2, 1, Number > &s, const SymmetricTensor< 4, 1, OtherNumber > &t) |
template<typename Number , typename OtherNumber > | |
constexpr void | double_contract (SymmetricTensor< 2, 2, typename ProductType< Number, OtherNumber >::type > &tmp, const SymmetricTensor< 4, 2, Number > &t, const SymmetricTensor< 2, 2, OtherNumber > &s) |
template<typename Number , typename OtherNumber > | |
constexpr void | double_contract (SymmetricTensor< 2, 2, typename ProductType< Number, OtherNumber >::type > &tmp, const SymmetricTensor< 2, 2, Number > &s, const SymmetricTensor< 4, 2, OtherNumber > &t) |
template<typename Number , typename OtherNumber > | |
constexpr void | double_contract (SymmetricTensor< 2, 3, typename ProductType< Number, OtherNumber >::type > &tmp, const SymmetricTensor< 4, 3, Number > &t, const SymmetricTensor< 2, 3, OtherNumber > &s) |
template<typename Number , typename OtherNumber > | |
constexpr void | double_contract (SymmetricTensor< 2, 3, typename ProductType< Number, OtherNumber >::type > &tmp, const SymmetricTensor< 2, 3, Number > &s, const SymmetricTensor< 4, 3, OtherNumber > &t) |
template<int dim, typename Number , typename OtherNumber > | |
constexpr Tensor< 1, dim, typename ProductType< Number, OtherNumber >::type > | operator* (const SymmetricTensor< 2, dim, Number > &src1, const Tensor< 1, dim, OtherNumber > &src2) |
template<int dim, typename Number , typename OtherNumber > | |
constexpr Tensor< 1, dim, typename ProductType< Number, OtherNumber >::type > | operator* (const Tensor< 1, dim, Number > &src1, const SymmetricTensor< 2, dim, OtherNumber > &src2) |
template<int rank_1, int rank_2, int dim, typename Number , typename OtherNumber > | |
constexpr Tensor< rank_1+rank_2 - 2, dim, typename ProductType< Number, OtherNumber >::type >::tensor_type | operator* (const Tensor< rank_1, dim, Number > &src1, const SymmetricTensor< rank_2, dim, OtherNumber > &src2) |
template<int rank_1, int rank_2, int dim, typename Number , typename OtherNumber > | |
constexpr Tensor< rank_1+rank_2 - 2, dim, typename ProductType< Number, OtherNumber >::type >::tensor_type | operator* (const SymmetricTensor< rank_1, dim, Number > &src1, const Tensor< rank_2, dim, OtherNumber > &src2) |
template<int dim, typename Number > | |
std::ostream & | operator<< (std::ostream &out, const SymmetricTensor< 2, dim, Number > &t) |
template<int dim, typename Number > | |
std::ostream & | operator<< (std::ostream &out, const SymmetricTensor< 4, dim, Number > &t) |
|
strong |
An enumeration for the algorithm to be employed when performing the computation of normalized eigenvectors and their corresponding eigenvalues by the eigenvalues() and eigenvectors() methods operating on SymmetricTensor objects.
The specialized algorithms utilized in computing the eigenvectors are presented in
Definition at line 3131 of file symmetric_tensor.h.
|
inlineconstexpr |
Return a unit symmetric tensor of rank 2, i.e., the \(\text{dim}\times\text{dim}\) identity matrix \(\mathbf I\).
Definition at line 3257 of file symmetric_tensor.h.
|
inlineconstexpr |
Return the tensor of rank 4 that, when multiplied by a symmetric rank 2 tensor \(\mathbf T\) returns the deviator \(\text{dev}\ \mathbf T\). It is the operator representation of the linear deviator operator \(\mathbb P\), also known as the volumetric projection tensor, calculated as:
\begin{align*} \mathbb{P} &=\mathbb{I} -\frac{1}{\text{dim}} \mathbf I \otimes \mathbf I \\ \mathcal{P}_{ijkl} &= \frac 12 \left(\delta_{ik} \delta_{jl} + \delta_{il} \delta_{jk} \right) - \frac{1}{\text{dim}} \delta_{ij} \delta_{kl} \end{align*}
For every tensor T
, there holds the identity deviator<dim,Number>(T) == deviator_tensor<dim,Number>() * T
, up to numerical round-off.
\[ \text{dev}\mathbf T = \mathbb P : \mathbf T \]
\[ \frac{\partial \text{dev}\mathbf{T}}{\partial \mathbf T} = \mathbb P. \]
Definition at line 3331 of file symmetric_tensor.h.
|
inlineconstexpr |
Return the fourth-order symmetric identity tensor \(\mathbb I\) which maps symmetric second-order tensors, such as \(\mathbf A\), to themselves.
\[ \mathbb I : \mathbf A = \mathbf A \]
Note that this tensor, even though it is the identity, has a somewhat funny form, and in particular does not only consist of zeros and ones. For example, for dim=2
, the identity tensor has all zero entries except for
\[ \mathcal{I}_{0000} = \mathcal{I}_{1111} = 1 \]
\[ \mathcal{I}_{0101} = \mathcal{I}_{0110} = \mathcal{I}_{1001} = \mathcal{I}_{1010} = \frac 12. \]
In index notation, we can write the general form
\[ \mathcal{I}_{ijkl} = \frac 12 \left( \delta_{ik} \delta_{jl} + \delta_{il} \delta_{jl} \right). \]
To see why this factor of \(1 / 2\) is necessary, consider computing \(\mathbf A= \mathbb I : \mathbf B\). For the element \(A_{01}\) we have \(A_{01} = \mathcal{I}_{0100} B_{00} + \mathcal{I}_{0111} B_{11} + \mathcal{I}_{0101} B_{01} + \mathcal{I}_{0110} B_{10}\). On the other hand, we need to have \(A_{01} = B_{01}\), and symmetry implies \(B_{01}=B_{10}\), leading to \(A_{01} = (\mathcal{I}_{0101} + \mathcal{I}_{0110}) B_{01}\), or, again by symmetry, \(\mathcal{I}_{0101} = \mathcal{I}_{0110} = \frac 12\). Similar considerations hold for the three-dimensional case.
This issue is also explained in the introduction to step-44.
Definition at line 3414 of file symmetric_tensor.h.
|
constexpr |
Invert a symmetric rank-2 tensor.
Definition at line 3467 of file symmetric_tensor.h.
|
constexpr |
Invert a symmetric rank-4 tensor. Since symmetric rank-4 tensors are mappings from and to symmetric rank-2 tensors, they can have an inverse.
If a tensor is not invertible, then the result is unspecified, but will likely contain the results of a division by zero or a very small number at the very least.
Definition at line 3488 of file symmetric_tensor.h.
|
inlineconstexpr |
|
inlineconstexpr |
Compute the deviator of a symmetric tensor, which is defined as \(\text{dev} \mathbf T = \mathbf T - \frac{1}{\text{dim}} \text{tr}\mathbf T \; \mathbf I\), where \(\mathbf I\) is the identity operator. This quantity equals the original tensor minus its contractive or dilative component and refers to the shear in, for example, elasticity.
Definition at line 3234 of file symmetric_tensor.h.
|
inlineconstexpr |
Compute the determinant of a rank 2 symmetric tensor. The determinant is also commonly referred to as the third invariant of rank-2 tensors.
For a one-dimensional tensor, the determinant equals the only element and is therefore equivalent to the trace.
For greater notational simplicity, there is also a third_invariant()
function that returns the determinant of a tensor.
Definition at line 2645 of file symmetric_tensor.h.
|
inlineconstexpr |
Addition of two symmetric tensors of equal rank. The result is another SymmetricTensor that has a number type that is compatible with the operation.
If possible (e.g. when Number
and OtherNumber
are of the same type, or if the result of Number() + OtherNumber()
is another Number
), you should use operator+=
instead since this does not require the creation of a temporary variable.
Definition at line 2525 of file symmetric_tensor.h.
|
inlineconstexpr |
Subtraction of two symmetric tensors of equal rank. The result is another SymmetricTensor that has a number type that is compatible with the operation.
If possible (e.g. when Number
and OtherNumber
are of the same type, or if the result of Number() - OtherNumber()
is another Number
), you should use operator-=
instead since this does not require the creation of a temporary variable.
Definition at line 2550 of file symmetric_tensor.h.
|
constexpr |
Addition of a SymmetricTensor and a general Tensor of equal rank. The result is a general Tensor that has a number type that is compatible with the operation.
Definition at line 2570 of file symmetric_tensor.h.
|
constexpr |
Addition of a general Tensor with a SymmetricTensor of equal rank. The result is a general Tensor that has a number type that is compatible with the operation.
Definition at line 2587 of file symmetric_tensor.h.
|
constexpr |
Subtraction of a general Tensor from a SymmetricTensor of equal rank. The result is a general Tensor that has a number type that is compatible with the operation.
Definition at line 2604 of file symmetric_tensor.h.
|
constexpr |
Subtraction of a SymmetricTensor from a general Tensor of equal rank. The result is a general Tensor that has a number type that is compatible with the operation.
Definition at line 2621 of file symmetric_tensor.h.
|
constexpr |
Compute the determinant of a rank 2 symmetric tensor. This function therefore computes the same value as the determinant()
functions and is only provided for greater notational simplicity (since there are also functions first_invariant() and second_invariant()).
\[ I_3 (\mathbf A) = III (\mathbf A) = \det (\mathbf A) \]
Definition at line 2686 of file symmetric_tensor.h.
|
inlineconstexpr |
Compute and return the trace of a tensor of rank 2, i.e. the sum of its diagonal entries. The trace is the first invariant of a rank-2 tensor.
\[ \text{tr} \mathbf A = \sum_i A_{ii} \]
Definition at line 2705 of file symmetric_tensor.h.
|
constexpr |
Compute the trace of a rank 2 symmetric tensor. This function therefore computes the same value as the trace()
functions and is only provided for greater notational simplicity (since there are also functions second_invariant() and third_invariant()).
\[ I_1 (\mathbf A) = I (\mathbf A) = \text{tr} \mathbf A = \sum_i A_{ii} \]
Definition at line 2728 of file symmetric_tensor.h.
|
constexpr |
Compute the second invariant of a tensor of rank 2. The second invariant of a tensor \(\mathbf A\) is defined as \(I_2 (\mathbf A) = II(\mathbf A) = \frac 12 \left[ (\text{tr} \mathbf A)^2 - \text{tr} (\mathbf{A}^2) \right]\).
For the kind of arguments to this function, i.e., a rank-2 tensor of size 1, the result is simply zero.
Definition at line 2748 of file symmetric_tensor.h.
|
constexpr |
Compute the second invariant of a tensor of rank 2. The second invariant of a tensor \(\mathbf A\) is defined as \(I_2 (\mathbf A) = II(\mathbf A) = \frac 12 \left[ (\text{tr} \mathbf A)^2 - \text{tr} (\mathbf{A}^2) \right]\).
For the kind of arguments to this function, i.e., a symmetric rank-2 tensor of size 2, the result is (counting indices starting at one) \(I_2(\mathbf A) = II(\mathbf A) = \frac 12 \left[ (A_{11} + A_{22})^2 - (A_{11}^2+2 A_{12}^2+ A_{22}^2) \right] = A_{11} A_{22} - A_{12}^2\). As expected, for the \(2\times 2\) symmetric tensors this function handles, this equals the determinant of the tensor. (This is so because for \(2\times 2\) symmetric tensors, there really are only two invariants, so the second and third invariant are the same; the determinant is the third invariant.)
Definition at line 2776 of file symmetric_tensor.h.
|
constexpr |
Compute the second invariant of a tensor of rank 2. The second invariant of a tensor \(\mathbf A\) is defined as \(I_2 (\mathbf A) = II(\mathbf A) = \frac 12 \left[ (\text{tr} \mathbf A)^2 - \text{tr} (\mathbf{A}^2) \right]\).
Definition at line 2794 of file symmetric_tensor.h.
std::array< Number, 1 > eigenvalues | ( | const SymmetricTensor< 2, 1, Number > & | T | ) |
Return the eigenvalues of a symmetric \(1 \times 1\) tensor. The (single) entry of the tensor is, of course, equal to the (single) eigenvalue.
std::array< Number, 2 > eigenvalues | ( | const SymmetricTensor< 2, 2, Number > & | T | ) |
Return the eigenvalues of a symmetric \(2\times 2\) tensor. The array of eigenvalues is sorted in descending order.
For \(2\times 2\) tensors, the eigenvalues of tensor \(\mathbf T\) are the roots of the characteristic polynomial \(0 = \lambda^2 - \lambda\;\text{tr}\mathbf{T} + \det \mathbf{T}\) as given by \(\lambda_1, \lambda_2 = \frac{1}{2} \left[ \text{tr} \mathbf{T} \pm \sqrt{(\text{tr} \mathbf{T})^2 - 4 \det \mathbf{T}} \right]\).
std::array< Number, 3 > eigenvalues | ( | const SymmetricTensor< 2, 3, Number > & | T | ) |
Return the eigenvalues of a symmetric \(3\times 3\) tensor. The array of eigenvalues is sorted in descending order.
For \(3\times 3\) tensors, the eigenvalues of tensor \(\mathbf T\) are the roots of the characteristic polynomial \(0 = \lambda^3 - \lambda^2\;\text{tr}\mathbf T - \frac{1}{2} \lambda \left[\text{tr}(\mathbf{T}^2) - (\text{tr}\mathbf T)^2\right] - \det \mathbf T\).
std::array< std::pair< Number, Tensor< 1, dim, Number > >, std::integral_constant< int, dim >::value > eigenvectors | ( | const SymmetricTensor< 2, dim, Number > & | T, |
const SymmetricTensorEigenvectorMethod | method = SymmetricTensorEigenvectorMethod::ql_implicit_shifts |
||
) |
Return the eigenvalues and eigenvectors of a real-valued rank-2 symmetric tensor \(\mathbf T\). The array of matched eigenvalue and eigenvector pairs is sorted in descending order (determined by the eigenvalues).
The specialized algorithms utilized in computing the eigenvectors are presented in
|
constexpr |
Return the transpose of the given symmetric tensor. Since we are working with symmetric objects, the transpose is of course the same as the original tensor. This function mainly exists for compatibility with the Tensor class.
Definition at line 3214 of file symmetric_tensor.h.
|
inlineconstexpr |
unit_symmetric_tensor<dim>() is the specialization of the function unit_symmetric_tensor<dim,Number>() which uses double
as the data type for the elements.
Definition at line 3293 of file symmetric_tensor.h.
|
inlineconstexpr |
This version of the deviator_tensor<dim>() function is a specialization of deviator_tensor<dim,Number>() that uses double
as the data type for the elements of the tensor.
Definition at line 3367 of file symmetric_tensor.h.
|
inlineconstexpr |
This version of the identity_tensor<dim>() function is the specialization of identity_tensor<dim,Number>() which uses double
as the data type for the elements of the tensor.
Definition at line 3448 of file symmetric_tensor.h.
|
inlineconstexpr |
Return the tensor of rank 4 that is the outer product of the two tensors given as arguments, i.e. the result \(\mathbb A = \mathbf{T}_1 \otimes \mathbf{T}_2\) satisfies \(\mathbb A : \mathbf B = (\mathbf{T}_2 : \mathbf B) \mathbf{T}_1\) for all symmetric tensors \(\mathbf B\). In index notation
\[ \mathcal{A}_{ijkl} = (T_1)_{ij} (T_2)_{kl} \]
For example, the deviator tensor \(\mathbb P = \mathbb I - \frac{1}{\text{dim}} \mathbf I \otimes \mathbf I\) can be computed as identity_tensor<dim>() - 1/d * outer_product (unit_symmetric_tensor<dim>(), unit_symmetric_tensor<dim>())
, since the (double) contraction with the unit tensor yields the trace of a symmetric tensor ( \(\mathbf I : \mathbf B = \text{tr} \mathbf B\)).
Definition at line 3520 of file symmetric_tensor.h.
|
inlineconstexpr |
Return the symmetrized version of a full rank-2 tensor, i.e. \(\text{sym}\mathbf A = \frac 12 \left(\mathbf A + \mathbf{A}^T\right)\), as a symmetric rank-2 tensor. This is the version for general dimensions.
Definition at line 3547 of file symmetric_tensor.h.
|
inlineconstexpr |
Multiplication of a symmetric tensor of general rank with a scalar from the right. This version of the operator is used if the scalar has the same data type as is used to store the elements of the symmetric tensor.
Definition at line 3572 of file symmetric_tensor.h.
|
constexpr |
Multiplication of a symmetric tensor of general rank with a scalar from the left. This version of the operator is used if the scalar has the same data type as is used to store the elements of the symmetric tensor.
Definition at line 3590 of file symmetric_tensor.h.
|
inlineconstexpr |
Multiplication of a symmetric tensor with a scalar number from the right.
The purpose of this operator is to enable only multiplication of a tensor by a scalar number (i.e., a floating point number, a complex floating point number, etc.). The function is written in a way that only allows the compiler to consider the function if the second argument is indeed a scalar number – in other words, OtherNumber
will not match, for example std::vector<double>
as the product of a tensor and a vector clearly would make no sense. The mechanism by which the compiler is prohibited of considering this operator for multiplication with non-scalar types are explained in the documentation of the EnableIfScalar class.
The return type of the function is chosen so that it matches the types of both the tensor and the scalar argument. For example, if you multiply a SymmetricTensor<2,dim,double>
by std::complex<double>
, then the result will be a SymmetricTensor<2,dim,std::complex<double>>
. In other words, the type with which the returned tensor stores its components equals the type you would get if you multiplied an individual component of the input tensor by the scalar factor.
Definition at line 3629 of file symmetric_tensor.h.
|
inlineconstexpr |
Multiplication of a symmetric tensor with a scalar number from the left. See the discussion with the operator with switched arguments for more information about template arguments and the return type.
Definition at line 3659 of file symmetric_tensor.h.
|
inlineconstexpr |
Division of a symmetric tensor of general rank by a scalar.
Definition at line 3679 of file symmetric_tensor.h.
|
inlineconstexpr |
Multiplication of a symmetric tensor of general rank with a scalar from the right.
Definition at line 3698 of file symmetric_tensor.h.
|
inlineconstexpr |
Multiplication of a symmetric tensor of general rank with a scalar from the left.
Definition at line 3715 of file symmetric_tensor.h.
|
inlineconstexpr |
Division of a symmetric tensor of general rank by a scalar.
Definition at line 3731 of file symmetric_tensor.h.
|
constexpr |
Compute the scalar product \(\mathbf A: \mathbf B=\sum_{i,j} A_{ij}B_{ij}\) between two tensors \(\mathbf A, \mathbf B\) of rank 2. In the current case where both arguments are symmetric tensors, this is equivalent to calling the expression A*B
which uses SymmetricTensor::operator*()
.
Definition at line 3749 of file symmetric_tensor.h.
|
inlineconstexpr |
Compute the scalar product \(\mathbf A: \mathbf B=\sum_{i,j} A_{ij}B_{ij}\) between two tensors \(\mathbf A, \mathbf B\) of rank 2. We don't use operator*
for this operation since the product between two tensors is usually assumed to be the contraction over the last index of the first tensor and the first index of the second tensor. For example, if B
is a Tensor, calling A*B
(instead of scalar_product(A,B)
) provides \((\mathbf A \cdot\mathbf B)_{ij}=\sum_k A_{ik}B_{kj}\).
Definition at line 3772 of file symmetric_tensor.h.
|
constexpr |
Compute the scalar product \(\mathbf A:\mathbf B=\sum_{i,j} A_{ij}B_{ij}\) between two tensors \(\mathbf A, \mathbf B\) of rank 2. We don't use operator*
for this operation since the product between two tensors is usually assumed to be the contraction over the last index of the first tensor and the first index of the second tensor. For example, if A
is a Tensor, calling A*B
(instead of scalar_product(A,B)
) provides \((\mathbf A \cdot\mathbf B)_{ij}=\sum_k A_{ik}B_{kj}\).
Definition at line 3799 of file symmetric_tensor.h.
|
inlineconstexpr |
Double contraction between a rank-4 and a rank-2 symmetric tensor, resulting in the symmetric tensor of rank 2 that is given as first argument to this function. This operation is the symmetric tensor analogon of a matrix-vector multiplication.
This function does the same as SymmetricTensor::operator*(). It should not be used, however, since the member operator has knowledge of the actual data storage format and is at least 2 orders of magnitude faster. This function mostly exists for compatibility purposes with the general Tensor class.
Definition at line 3822 of file symmetric_tensor.h.
|
inlineconstexpr |
Double contraction between a rank-4 and a rank-2 symmetric tensor, resulting in the symmetric tensor of rank 2 that is given as first argument to this function. This operation is the symmetric tensor analogon of a matrix-vector multiplication.
This function does the same as SymmetricTensor::operator*(). It should not be used, however, since the member operator has knowledge of the actual data storage format and is at least 2 orders of magnitude faster. This function mostly exists for compatibility purposes with the general Tensor class.
Definition at line 3848 of file symmetric_tensor.h.
|
inlineconstexpr |
Double contraction between a rank-4 and a rank-2 symmetric tensor, resulting in the symmetric tensor of rank 2 that is given as first argument to this function. This operation is the symmetric tensor analogon of a matrix-vector multiplication.
This function does the same as SymmetricTensor::operator*(). It should not be used, however, since the member operator has knowledge of the actual data storage format and is at least 2 orders of magnitude faster. This function mostly exists for compatibility purposes with the general Tensor class.
Definition at line 3874 of file symmetric_tensor.h.
|
inlineconstexpr |
Double contraction between a rank-4 and a rank-2 symmetric tensor, resulting in the symmetric tensor of rank 2 that is given as first argument to this function. This operation is the symmetric tensor analogon of a matrix-vector multiplication.
This function does the same as SymmetricTensor::operator*(). It should not be used, however, since the member operator has knowledge of the actual data storage format and is at least 2 orders of magnitude faster. This function mostly exists for compatibility purposes with the general Tensor class.
Definition at line 3905 of file symmetric_tensor.h.
|
inlineconstexpr |
Double contraction between a rank-4 and a rank-2 symmetric tensor, resulting in the symmetric tensor of rank 2 that is given as first argument to this function. This operation is the symmetric tensor analogon of a matrix-vector multiplication.
This function does the same as SymmetricTensor::operator*(). It should not be used, however, since the member operator has knowledge of the actual data storage format and is at least 2 orders of magnitude faster. This function mostly exists for compatibility purposes with the general Tensor class.
Definition at line 3936 of file symmetric_tensor.h.
|
inlineconstexpr |
Double contraction between a rank-4 and a rank-2 symmetric tensor, resulting in the symmetric tensor of rank 2 that is given as first argument to this function. This operation is the symmetric tensor analogon of a matrix-vector multiplication.
This function does the same as SymmetricTensor::operator*(). It should not be used, however, since the member operator has knowledge of the actual data storage format and is at least 2 orders of magnitude faster. This function mostly exists for compatibility purposes with the general Tensor class.
Definition at line 3968 of file symmetric_tensor.h.
|
constexpr |
Multiply a symmetric rank-2 tensor (i.e., a matrix) by a rank-1 tensor (i.e., a vector). The result is a rank-1 tensor (i.e., a vector).
Definition at line 3994 of file symmetric_tensor.h.
|
constexpr |
Multiply a rank-1 tensor (i.e., a vector) by a symmetric rank-2 tensor (i.e., a matrix). The result is a rank-1 tensor (i.e., a vector).
Definition at line 4014 of file symmetric_tensor.h.
|
constexpr |
The dot product (single contraction) for tensors: Return a tensor of rank \((\text{rank}_1 + \text{rank}_2 - 2)\) that is the contraction of the last index of a tensor src1
of rank rank_1
with the first index of a tensor src2
of rank rank_2:
\[ \text{result}_{i_1,\ldots,i_{r1},j_1,\ldots,j_{r2}} = \sum_{k} \text{left}_{i_1,\ldots,i_{r1}, k} \text{right}_{k, j_1,\ldots,j_{r2}} \]
Definition at line 4052 of file symmetric_tensor.h.
|
constexpr |
The dot product (single contraction) for tensors: Return a tensor of rank \((\text{rank}_1 + \text{rank}_2 - 2)\) that is the contraction of the last index of a tensor src1
of rank rank_1
with the first index of a tensor src2
of rank rank_2:
\[ \text{result}_{i_1,\ldots,i_{r1},j_1,\ldots,j_{r2}} = \sum_{k} \text{left}_{i_1,\ldots,i_{r1}, k} \text{right}_{k, j_1,\ldots,j_{r2}} \]
Definition at line 4089 of file symmetric_tensor.h.
|
inline |
Output operator for symmetric tensors of rank 2. Print the elements consecutively, with a space in between, two spaces between rank 1 subtensors, three between rank 2 and so on. No special amends are made to represents the symmetry in the output, for example by outputting only the unique entries.
Definition at line 4108 of file symmetric_tensor.h.
|
inline |
Output operator for symmetric tensors of rank 4. Print the elements consecutively, with a space in between, two spaces between rank 1 subtensors, three between rank 2 and so on. No special amends are made to represents the symmetry in the output, for example by outputting only the unique entries.
Definition at line 4135 of file symmetric_tensor.h.