Reference documentation for deal.II version 9.3.3
|
#include <deal.II/base/config.h>
#include <deal.II/base/exceptions.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_accessors.h>
#include <deal.II/base/utilities.h>
#include <adolc/adouble.h>
#include <cmath>
#include <ostream>
#include <utility>
#include <vector>
Go to the source code of this file.
Classes | |
class | Tensor< 0, dim, Number > |
class | Tensor< rank_, dim, Number > |
Namespaces | |
namespace | internal |
namespace | internal::TensorImplementation |
Functions | |
template<int rank, int dim, typename Number , typename OtherNumber , typename std::enable_if< !std::is_integral< typename ProductType< Number, OtherNumber >::type >::value, int >::type = 0> | |
constexpr Tensor< rank, dim, typename ProductType< Number, OtherNumber >::type > | internal::TensorImplementation::division_operator (const Tensor< rank, dim, Number > &t, const OtherNumber &factor) |
Output functions for Tensor objects | |
template<int rank_, int dim, typename Number > | |
std::ostream & | operator<< (std::ostream &out, const Tensor< rank_, dim, Number > &p) |
template<int dim, typename Number > | |
std::ostream & | operator<< (std::ostream &out, const Tensor< 0, dim, Number > &p) |
Vector space operations on Tensor objects: | |
template<int dim, typename Number , typename Other > | |
constexpr ProductType< Other, Number >::type | operator* (const Other &object, const Tensor< 0, dim, Number > &t) |
template<int dim, typename Number , typename Other > | |
constexpr ProductType< Number, Other >::type | operator* (const Tensor< 0, dim, Number > &t, const Other &object) |
template<int dim, typename Number , typename OtherNumber > | |
constexpr ProductType< Number, OtherNumber >::type | operator* (const Tensor< 0, dim, Number > &src1, const Tensor< 0, dim, OtherNumber > &src2) |
template<int dim, typename Number , typename OtherNumber > | |
constexpr Tensor< 0, dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > | operator/ (const Tensor< 0, dim, Number > &t, const OtherNumber &factor) |
template<int dim, typename Number , typename OtherNumber > | |
constexpr Tensor< 0, dim, typename ProductType< Number, OtherNumber >::type > | operator+ (const Tensor< 0, dim, Number > &p, const Tensor< 0, dim, OtherNumber > &q) |
template<int dim, typename Number , typename OtherNumber > | |
constexpr Tensor< 0, dim, typename ProductType< Number, OtherNumber >::type > | operator- (const Tensor< 0, dim, Number > &p, const Tensor< 0, dim, OtherNumber > &q) |
template<int rank, int dim, typename Number , typename OtherNumber > | |
constexpr Tensor< rank, dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > | operator* (const Tensor< rank, dim, Number > &t, const OtherNumber &factor) |
template<int rank, int dim, typename Number , typename OtherNumber > | |
constexpr Tensor< rank, dim, typename ProductType< typename EnableIfScalar< Number >::type, OtherNumber >::type > | operator* (const Number &factor, const Tensor< rank, dim, OtherNumber > &t) |
template<int rank, int dim, typename Number , typename OtherNumber > | |
constexpr Tensor< rank, dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > | operator/ (const Tensor< rank, dim, Number > &t, const OtherNumber &factor) |
template<int rank, int dim, typename Number , typename OtherNumber > | |
constexpr Tensor< rank, dim, typename ProductType< Number, OtherNumber >::type > | operator+ (const Tensor< rank, dim, Number > &p, const Tensor< rank, dim, OtherNumber > &q) |
template<int rank, int dim, typename Number , typename OtherNumber > | |
constexpr Tensor< rank, dim, typename ProductType< Number, OtherNumber >::type > | operator- (const Tensor< rank, dim, Number > &p, const Tensor< rank, dim, OtherNumber > &q) |
template<int dim, typename Number , typename OtherNumber > | |
constexpr Tensor< 0, dim, typename ProductType< Number, OtherNumber >::type > | schur_product (const Tensor< 0, dim, Number > &src1, const Tensor< 0, dim, OtherNumber > &src2) |
template<int rank, int dim, typename Number , typename OtherNumber > | |
constexpr Tensor< rank, dim, typename ProductType< Number, OtherNumber >::type > | schur_product (const Tensor< rank, dim, Number > &src1, const Tensor< rank, dim, OtherNumber > &src2) |
Contraction operations and the outer product for tensor objects | |
template<int rank_1, int rank_2, int dim, typename Number , typename OtherNumber , typename = typename std::enable_if<rank_1 >= 1 && rank_2> | |
OtherNumber ::type ::tensor_type | operator* (const Tensor< rank_1, dim, Number > &src1, const Tensor< rank_2, dim, OtherNumber > &src2) |
template<int index_1, int index_2, int rank_1, int rank_2, int dim, typename Number , typename OtherNumber > | |
constexpr Tensor< rank_1+rank_2-2, dim, typenameProductType< Number, OtherNumber >::type >::tensor_type | contract (const Tensor< rank_1, dim, Number > &src1, const Tensor< rank_2, dim, OtherNumber > &src2) |
template<int index_1, int index_2, int index_3, int index_4, int rank_1, int rank_2, int dim, typename Number , typename OtherNumber > | |
constexpr Tensor< rank_1+rank_2-4, dim, typenameProductType< Number, OtherNumber >::type >::tensor_type | double_contract (const Tensor< rank_1, dim, Number > &src1, const Tensor< rank_2, dim, OtherNumber > &src2) |
template<int rank, int dim, typename Number , typename OtherNumber > | |
constexpr ProductType< Number, OtherNumber >::type | scalar_product (const Tensor< rank, dim, Number > &left, const Tensor< rank, dim, OtherNumber > &right) |
template<template< int, int, typename > class TensorT1, template< int, int, typename > class TensorT2, template< int, int, typename > class TensorT3, int rank_1, int rank_2, int dim, typename T1 , typename T2 , typename T3 > | |
constexpr ProductType< T1, typenameProductType< T2, T3 >::type >::type | contract3 (const TensorT1< rank_1, dim, T1 > &left, const TensorT2< rank_1+rank_2, dim, T2 > &middle, const TensorT3< rank_2, dim, T3 > &right) |
template<int rank_1, int rank_2, int dim, typename Number , typename OtherNumber > | |
constexpr Tensor< rank_1+rank_2, dim, typename ProductType< Number, OtherNumber >::type > | outer_product (const Tensor< rank_1, dim, Number > &src1, const Tensor< rank_2, dim, OtherNumber > &src2) |
Special operations on tensors of rank 1 | |
template<int dim, typename Number > | |
constexpr Tensor< 1, dim, Number > | cross_product_2d (const Tensor< 1, dim, Number > &src) |
template<int dim, typename Number1 , typename Number2 > | |
constexpr Tensor< 1, dim, typename ProductType< Number1, Number2 >::type > | cross_product_3d (const Tensor< 1, dim, Number1 > &src1, const Tensor< 1, dim, Number2 > &src2) |
Special operations on tensors of rank 2 | |
template<int dim, typename Number > | |
constexpr Number | determinant (const Tensor< 2, dim, Number > &t) |
template<typename Number > | |
constexpr Number | determinant (const Tensor< 2, 1, Number > &t) |
template<typename Number > | |
constexpr Number | determinant (const Tensor< 2, 2, Number > &t) |
template<typename Number > | |
constexpr Number | determinant (const Tensor< 2, 3, Number > &t) |
template<int dim, typename Number > | |
constexpr Number | trace (const Tensor< 2, dim, Number > &d) |
template<int dim, typename Number > | |
constexpr Tensor< 2, dim, Number > | invert (const Tensor< 2, dim, Number > &) |
template<int dim, typename Number > | |
constexpr Tensor< 2, dim, Number > | transpose (const Tensor< 2, dim, Number > &t) |
template<int dim, typename Number > | |
constexpr Tensor< 2, dim, Number > | adjugate (const Tensor< 2, dim, Number > &t) |
template<int dim, typename Number > | |
constexpr Tensor< 2, dim, Number > | cofactor (const Tensor< 2, dim, Number > &t) |
template<int dim, typename Number > | |
Tensor< 2, dim, Number > | project_onto_orthogonal_tensors (const Tensor< 2, dim, Number > &A) |
template<int dim, typename Number > | |
Number | l1_norm (const Tensor< 2, dim, Number > &t) |
template<int dim, typename Number > | |
Number | linfty_norm (const Tensor< 2, dim, Number > &t) |
|
inline |
|
inlineconstexpr |
|
inlineconstexpr |
|
constexpr |
Scalar multiplication of two tensors of rank 0.
This function unwraps the underlying objects of type Number
and OtherNumber
that are stored within the Tensor and multiplies them. It returns an unwrapped number of product type.
|
constexpr |
|
constexpr |
|
constexpr |
|
inlineconstexpr |
Multiplication of a tensor of general rank with a scalar number from the right.
Only multiplication with a scalar number type (i.e., a floating point number, a complex floating point number, etc.) is allowed, see the documentation of EnableIfScalar for details.
|
inlineconstexpr |
Multiplication of a tensor of general rank with a scalar number from the left.
Only multiplication with a scalar number type (i.e., a floating point number, a complex floating point number, etc.) is allowed, see the documentation of EnableIfScalar for details.
|
inlineconstexpr |
Division of a tensor of general rank with a scalar number. See the discussion on operator*() above for more information about template arguments and the return type.
|
inlineconstexpr |
|
inlineconstexpr |
|
inlineconstexpr |
|
inlineconstexpr |
Entrywise multiplication of two tensor objects of general rank.
This multiplication is also called "Hadamard-product" (c.f. https://en.wikipedia.org/wiki/Hadamard_product_(matrices)), and generates a new tensor of size <rank, dim>:
\[ \text{result}_{i, j} = \text{left}_{i, j}\circ \text{right}_{i, j} \]
rank | The rank of both tensors. |
|
inlineconstexpr |
Generic contraction of a pair of indices of two tensors of arbitrary rank: Return a tensor of rank \((\text{rank}_1 + \text{rank}_2 - 2)\) that is the contraction of index index_1
of a tensor src1
of rank rank_1
with the index index_2
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,k,\ldots,i_{r1}} \text{right}_{j_1,\ldots,k,\ldots,j_{r2}} \]
If for example the first index (index_1==0
) of a tensor t1
shall be contracted with the third index (index_2==2
) of a tensor t2
, this function should be invoked as
|
inlineconstexpr |
Generic contraction of two pairs of indices of two tensors of arbitrary rank: Return a tensor of rank \((\text{rank}_1 + \text{rank}_2 - 4)\) that is the contraction of index index_1
with index index_2
, and index index_3
with index index_4
of a tensor src1
of rank rank_1
and a tensor src2
of rank rank_2:
\[ \text{result}_{i_1,\ldots,i_{r1},j_1,\ldots,j_{r2}} = \sum_{k, l} \text{left}_{i_1,\ldots,k,\ldots,l,\ldots,i_{r1}} \text{right}_{j_1,\ldots,k,\ldots,l\ldots,j_{r2}} \]
If for example the first index (index_1==0
) shall be contracted with the third index (index_2==2
), and the second index (index_3==1
) with the first index (index_4==0
) of a tensor t2
, this function should be invoked as
|
inlineconstexpr |
The scalar product, or (generalized) Frobenius inner product of two tensors of equal rank: Return a scalar number that is the result of a full contraction of a tensor left
and right:
\[ \sum_{i_1,\ldots,i_r} \text{left}_{i_1,\ldots,i_r} \text{right}_{i_1,\ldots,i_r} \]
|
inlineconstexpr |
Full contraction of three tensors: Return a scalar number that is the result of a full contraction of a tensor left
of rank rank_1
, a tensor middle
of rank \((\text{rank}_1+\text{rank}_2)\) and a tensor right
of rank rank_2:
\[ \sum_{i_1,\ldots,i_{r1},j_1,\ldots,j_{r2}} \text{left}_{i_1,\ldots,i_{r1}} \text{middle}_{i_1,\ldots,i_{r1},j_1,\ldots,j_{r2}} \text{right}_{j_1,\ldots,j_{r2}} \]
|
inlineconstexpr |
|
inlineconstexpr |
Return the cross product in 2d. This is just a rotation by 90 degrees clockwise to compute the outer normal from a tangential vector. This function is defined for all space dimensions to allow for dimension independent programming (e.g. within switches over the space dimension), but may only be called if the actual dimension of the arguments is two (e.g. from the dim==2
case in the switch).
|
inlineconstexpr |
Return the cross product of 2 vectors in 3d. This function is defined for all space dimensions to allow for dimension independent programming (e.g. within switches over the space dimension), but may only be called if the actual dimension of the arguments is three (e.g. from the dim==3
case in the switch).
|
constexpr |
|
constexpr |
|
constexpr |
|
constexpr |
Return the cofactor of the given tensor of rank 2. The cofactor of a tensor \(\mathbf A\) is defined as
\[ \textrm{cof}\mathbf A \dealcoloneq \textrm{det}\mathbf A \; \mathbf{A}^{-T} = \left[ \textrm{adj}\mathbf A \right]^{T} \; . \]
Tensor< 2, dim, Number > project_onto_orthogonal_tensors | ( | const Tensor< 2, dim, Number > & | A | ) |
Return the nearest orthogonal matrix \(\hat {\mathbf A}=\mathbf U \mathbf{V}^T\) by combining the products of the singular value decomposition (SVD) \({\mathbf A}=\mathbf U \mathbf S \mathbf V^T\) for a given input \({\mathbf A}\), effectively replacing \(\mathbf S\) with the identity matrix.
This is a (nonlinear) projection operation since when applied twice, we have \(\hat{\hat{\mathbf A}}=\hat{\mathbf A}\) as is easy to see. (That is because the SVD of \(\hat {\mathbf A}\) is simply \(\mathbf U \mathbf I \mathbf{V}^T\).) Furthermore, \(\hat {\mathbf A}\) is really an orthogonal matrix because orthogonal matrices have to satisfy \({\hat {\mathbf A}}^T \hat {\mathbf A}={\mathbf I}\), which here implies that
\begin{align*} {\hat {\mathbf A}}^T \hat {\mathbf A} &= \left(\mathbf U \mathbf{V}^T\right)^T\left(\mathbf U \mathbf{V}^T\right) \\ &= \mathbf V \mathbf{U}^T \mathbf U \mathbf{V}^T \\ &= \mathbf V \left(\mathbf{U}^T \mathbf U\right) \mathbf{V}^T \\ &= \mathbf V \mathbf I \mathbf{V}^T \\ &= \mathbf V \mathbf{V}^T \\ &= \mathbf I \end{align*}
due to the fact that the \(\mathbf U\) and \(\mathbf V\) factors that come out of the SVD are themselves orthogonal matrices.
A | The tensor for which to find the closest orthogonal tensor. |
Number | The type used to store the entries of the tensor. Must be either float or double . |
A
must not be singular. This is because, conceptually, the problem to be solved here is trying to find a matrix \(\hat{\mathbf A}\) that minimizes some kind of distance from \(\mathbf A\) while satisfying the quadratic constraint \({\hat {\mathbf A}}^T \hat {\mathbf A}={\mathbf I}\). This is not so dissimilar to the kind of problem where one wants to find a vector \(\hat{\mathbf x}\in{\mathbb R}^n\) that minimizes the quadratic objective function \(\|\hat {\mathbf x} - \mathbf x\|^2\) for a given \(\mathbf x\) subject to the constraint \(\|\mathbf x\|^2=1\) – in other words, we are seeking the point \(\hat{\mathbf x}\) on the unit sphere that is closest to \(\mathbf x\). This problem has a solution for all \(\mathbf x\) except if \(\mathbf x=0\). The corresponding condition for the problem we are considering here is that \(\mathbf A\) must not have a zero eigenvalue.