Reference documentation for deal.II version 9.1.1
|
#include <deal.II/base/tensor.h>
Public Types | |
using | value_type = typename Tensor< rank_ - 1, dim, Number >::tensor_type |
using | array_type = typename Tensor< rank_ - 1, dim, Number >::array_type[(dim !=0) ? dim :1] |
using | tensor_type = Tensor< rank_, dim, Number > |
Public Member Functions | |
constexpr | Tensor ()=default |
Tensor (const array_type &initializer) | |
template<typename OtherNumber > | |
Tensor (const Tensor< rank_, dim, OtherNumber > &initializer) | |
template<typename OtherNumber > | |
Tensor (const Tensor< 1, dim, Tensor< rank_ - 1, dim, OtherNumber >> &initializer) | |
template<typename OtherNumber > | |
constexpr | operator Tensor< 1, dim, Tensor< rank_ - 1, dim, OtherNumber >> () const |
value_type & | operator[] (const unsigned int i) |
constexpr const value_type & | operator[] (const unsigned int i) const |
const Number & | operator[] (const TableIndices< rank_ > &indices) const |
Number & | operator[] (const TableIndices< rank_ > &indices) |
Number * | begin_raw () |
constexpr const Number * | begin_raw () const |
Number * | end_raw () |
constexpr const Number * | end_raw () const |
template<typename OtherNumber > | |
Tensor & | operator= (const Tensor< rank_, dim, OtherNumber > &rhs) |
Tensor & | operator= (const Number &d) |
template<typename OtherNumber > | |
bool | operator== (const Tensor< rank_, dim, OtherNumber > &) const |
template<typename OtherNumber > | |
constexpr bool | operator!= (const Tensor< rank_, dim, OtherNumber > &) const |
template<typename OtherNumber > | |
Tensor & | operator+= (const Tensor< rank_, dim, OtherNumber > &) |
template<typename OtherNumber > | |
Tensor & | operator-= (const Tensor< rank_, dim, OtherNumber > &) |
template<typename OtherNumber > | |
Tensor & | operator*= (const OtherNumber &factor) |
template<typename OtherNumber > | |
Tensor & | operator/= (const OtherNumber &factor) |
Tensor | operator- () const |
void | clear () |
numbers::NumberTraits< Number >::real_type | norm () const |
numbers::NumberTraits< Number >::real_type | norm_square () const |
template<typename OtherNumber > | |
void | unroll (Vector< OtherNumber > &result) const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Static Public Member Functions | |
static unsigned int | component_to_unrolled_index (const TableIndices< rank_ > &indices) |
static TableIndices< rank_ > | unrolled_to_component_indices (const unsigned int i) |
static constexpr std::size_t | memory_consumption () |
Static Public Attributes | |
static constexpr unsigned int | dimension = dim |
static constexpr unsigned int | rank = rank_ |
static constexpr unsigned int | n_independent_components |
Private Member Functions | |
template<typename OtherNumber > | |
void | unroll_recursion (Vector< OtherNumber > &result, unsigned int &start_index) const |
Private Attributes | |
Tensor< rank_ - 1, dim, Number > | values [(dim !=0) ? dim :1] |
Friends | |
template<int , int , typename > | |
class | Tensor |
class | Point< dim, Number > |
Related Functions | |
(Note that these are not member functions.) | |
template<int rank, int dim, typename Number > | |
Tensor< rank, dim, Number > | sum (const Tensor< rank, dim, Number > &local, const MPI_Comm &mpi_communicator) |
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 > | |
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 > | |
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 > | |
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 > | |
Tensor< rank, dim, typename ProductType< Number, OtherNumber >::type > | operator- (const Tensor< rank, dim, Number > &p, const Tensor< rank, dim, OtherNumber > &q) |
Contraction operations and the outer product for tensor objects | |
template<int rank_1, int rank_2, int dim, typename Number , typename OtherNumber > | |
Tensor< rank_1+rank_2 - 2, dim, typename ProductType< Number, 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 > | |
Tensor< rank_1+rank_2 - 2, dim, typename ProductType< 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 > | |
Tensor< rank_1+rank_2 - 4, dim, typename ProductType< 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 > | |
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 > | |
ProductType< T1, typename ProductType< 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 > | |
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 > | |
Tensor< 1, dim, Number > | cross_product_2d (const Tensor< 1, dim, Number > &src) |
template<int dim, typename Number > | |
Tensor< 1, dim, Number > | cross_product_3d (const Tensor< 1, dim, Number > &src1, const Tensor< 1, dim, Number > &src2) |
Special operations on tensors of rank 2 | |
template<int dim, typename Number > | |
Number | determinant (const Tensor< 2, dim, Number > &t) |
template<typename Number > | |
constexpr Number | determinant (const Tensor< 2, 1, Number > &t) |
template<int dim, typename Number > | |
Number | trace (const Tensor< 2, dim, Number > &d) |
template<int dim, typename Number > | |
Tensor< 2, dim, Number > | invert (const Tensor< 2, dim, Number > &) |
template<int dim, typename Number > | |
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 > | |
Number | l1_norm (const Tensor< 2, dim, Number > &t) |
template<int dim, typename Number > | |
Number | linfty_norm (const Tensor< 2, dim, Number > &t) |
Deprecated Tensor operations | |
template<int dim, typename Number > | |
Number | double_contract (const Tensor< 2, dim, Number > &src1, const Tensor< 2, dim, Number > &src2) |
template<int dim, typename Number > | |
void | double_contract (Tensor< 2, dim, Number > &dest, const Tensor< 4, dim, Number > &src1, const Tensor< 2, dim, Number > &src2) |
template<int dim, typename Number > | |
void | contract (Tensor< 2, dim, Number > &dest, const Tensor< 2, dim, Number > &src1, const unsigned int index1, const Tensor< 2, dim, Number > &src2, const unsigned int index3) |
template<int dim, typename Number > | |
void | contract (Tensor< 2, dim, Number > &dest, const Tensor< 3, dim, Number > &src1, const unsigned int index1, const Tensor< 1, dim, Number > &src2) |
template<int dim, typename Number > | |
void | contract (Tensor< 3, dim, Number > &dest, const Tensor< 3, dim, Number > &src1, const unsigned int index1, const Tensor< 2, dim, Number > &src2, const unsigned int index2) |
template<int rank_1, int rank_2, int dim, typename Number > | |
void | contract (Tensor< rank_1+rank_2 - 2, dim, Number > &dest, const Tensor< rank_1, dim, Number > &src1, const Tensor< rank_2, dim, Number > &src2) |
template<int dim, typename Number , typename OtherNumber > | |
ProductType< Number, OtherNumber >::type | contract (const Tensor< 1, dim, Number > &src1, const Tensor< 1, dim, OtherNumber > &src2) |
template<int dim, typename Number > | |
void | cross_product (Tensor< 1, dim, Number > &dst, const Tensor< 1, dim, Number > &src) |
template<int dim, typename Number > | |
void | cross_product (Tensor< 1, dim, Number > &dst, const Tensor< 1, dim, Number > &src1, const Tensor< 1, dim, Number > &src2) |
template<int rank_1, int rank_2, int dim, typename Number > | |
void | outer_product (Tensor< rank_1+rank_2, dim, Number > &dst, const Tensor< rank_1, dim, Number > &src1, const Tensor< rank_2, dim, Number > &src2) |
template<int dim, typename Number > | |
void | outer_product (Tensor< 1, dim, Number > &dst, const Number src1, const Tensor< 1, dim, Number > &src2) |
template<int dim, typename Number > | |
void | outer_product (Tensor< 1, dim, Number > &dst, const Tensor< 1, dim, Number > src1, const Number src2) |
template<int rank, typename Number > | |
Number | determinant (const Tensor< rank, 1, Number > &t) |
template<typename Number > | |
Number | determinant (const Tensor< 1, 1, Number > &t) |
A general tensor class with an arbitrary rank, i.e. with an arbitrary number of indices. The Tensor class provides an indexing operator and a bit of infrastructure, but most functionality is recursively handed down to tensors of rank 1 or put into external templated functions, e.g. the contract
family.
Using this tensor class for objects of rank 2 has advantages over matrices in many cases since the dimension is known to the compiler as well as the location of the data. It is therefore possible to produce far more efficient code than for matrices with runtime-dependent dimension. It also makes the code easier to read because of the semantic difference between a tensor (an object that relates to a coordinate system and has transformation properties with regard to coordinate rotations and transforms) and matrices (which we consider as operators on arbitrary vector spaces related to linear algebra things).
rank_ | An integer that denotes the rank of this tensor. A rank-0 tensor is a scalar, a rank-1 tensor is a vector with dim components, a rank-2 tensor is a matrix with dim-by-dim components, etc. There are specializations of this class for rank-0 and rank-1 tensors. There is also a related class SymmetricTensor for tensors of even rank whose elements are symmetric. |
dim | An integer that denotes the dimension of the space in which this tensor operates. This of course equals the number of coordinates that identify a point and rank-1 tensor. |
Number | The data type in which the tensor elements are to be stored. This will, in almost all cases, simply be the default double , but there are cases where one may want to store elements in a different (and always scalar) type. It can be used to base tensors on float or complex numbers or any other data type that implements basic arithmetic operations. Another example would be a type that allows for Automatic Differentiation (see, for example, the Sacado type used in step-33) and thereby can generate analytic (spatial) derivatives of a function that takes a tensor as argument. |
using Tensor< rank_, dim, Number >::value_type = typename Tensor<rank_ - 1, dim, Number>::tensor_type |
Type of objects encapsulated by this container and returned by operator[](). This is a tensor of lower rank for a general tensor, and a scalar number type for Tensor<1,dim,Number>.
using Tensor< rank_, dim, Number >::array_type = typename Tensor<rank_ - 1, dim, Number>::array_type[(dim != 0) ? dim : 1] |
using Tensor< rank_, dim, Number >::tensor_type = Tensor<rank_, dim, Number> |
Internal type declaration that is used to specialize the return type of operator[]() for Tensor<1,dim,Number>
|
default |
Constructor. Initialize all entries to zero.
|
inlineexplicit |
|
inline |
constexpr const Tensor< rank_, dim, Number >::value_type & Tensor< rank_, dim, Number >::operator[] | ( | const unsigned int | i | ) | const |
|
inline |
Read access using TableIndices indices
|
inline |
Read and write access using TableIndices indices
|
inline |
constexpr const Number * Tensor< rank_, dim, Number >::begin_raw | ( | ) | const |
|
inline |
constexpr const Number * Tensor< rank_, dim, Number >::end_raw | ( | ) | const |
Tensor& Tensor< rank_, dim, Number >::operator= | ( | const Tensor< rank_, dim, OtherNumber > & | rhs | ) |
Assignment operator from tensors with different underlying scalar type. This obviously requires that the OtherNumber
type is convertible to Number
.
Tensor& Tensor< rank_, dim, Number >::operator+= | ( | const Tensor< rank_, dim, OtherNumber > & | ) |
Add another tensor.
Tensor& Tensor< rank_, dim, Number >::operator-= | ( | const Tensor< rank_, dim, OtherNumber > & | ) |
Subtract another tensor.
Tensor& Tensor< rank_, dim, Number >::operator*= | ( | const OtherNumber & | factor | ) |
Scale the tensor by factor
, i.e. multiply all components by factor
.
Tensor& Tensor< rank_, dim, Number >::operator/= | ( | const OtherNumber & | factor | ) |
Scale the vector by 1/factor
.
|
inline |
Reset all values to zero.
Note that this is partly inconsistent with the semantics of the clear()
member functions of the standard library containers and of several other classes within deal.II, which not only reset the values of stored elements to zero, but release all memory and return the object into a virginial state. However, since the size of objects of the present type is determined by its template parameters, resizing is not an option, and indeed the state where all elements have a zero value is the state right after construction of such an object.
|
inline |
|
inline |
|
inlinestatic |
|
inlinestatic |
|
static |
|
inline |
|
friend |
|
friend |
|
related |
Perform an MPI sum of the entries of a tensor.
|
related |
|
related |
|
related |
|
related |
|
related |
|
related |
|
related |
|
related |
|
related |
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.
|
related |
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.
|
related |
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.
|
related |
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}} \]
|
related |
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
, the invocation of this function is
|
related |
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
) the invocation of this function is this function is
|
related |
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} \]
|
related |
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}} \]
|
related |
|
related |
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).
|
related |
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).
|
related |
|
related |
|
related |
|
related |
Return the adjugate of the given tensor of rank 2. The adjugate of a tensor \(\left(\bullet\right)\) is defined as
\[ \textrm{adj}\left(\bullet\right) \dealcoloneq \textrm{det}\left(\bullet\right) \; \left(\bullet\right)^{-1} \; . \]
|
related |
Return the cofactor of the given tensor of rank 2. The cofactor of a tensor \(\left(\bullet\right)\) is defined as
\[ \textrm{cof}\left(\bullet\right) \dealcoloneq \textrm{det}\left(\bullet\right) \; \left(\bullet\right)^{-T} = \left[ \textrm{adj}\left(\bullet\right) \right]^{T} \; . \]
|
related |
|
related |
|
related |
Double contract two tensors of rank 2, thus computing the Frobenius inner product sumi,j src1[i][j]*src2[i][j]
.
Definition at line 228 of file tensor_deprecated.h.
|
related |
Contract the last two indices of src1
with the two indices src2
, creating a rank-2 tensor. This is the matrix-vector product analog operation between tensors of rank 4 and rank 2.
Definition at line 239 of file tensor_deprecated.h.
|
related |
Contract a tensor of rank 2 with a tensor of rank 2. The contraction is performed over index index1
of the first tensor, and index2
of the second tensor. Note that the number of the index is counted from 1 on, not from zero as usual.
Definition at line 252 of file tensor_deprecated.h.
|
related |
Contract a tensor of rank 3 with a tensor of rank 1. The contraction is performed over index index1
of the first tensor. Note that the number of the index is counted from 1 on, not from zero as usual.
Definition at line 309 of file tensor_deprecated.h.
|
related |
Contract a tensor of rank 3 with a tensor of rank 2. The contraction is performed over index index1
of the first tensor, and index2
of the second tensor. Note that the number of the index is counted from 1 on, not from zero as usual.
Definition at line 345 of file tensor_deprecated.h.
|
related |
Single contraction for tensors: contract the last index of a tensor src1
of rank rank_1
with the first index of a tensor src2
of rank rank_2
.
Definition at line 427 of file tensor_deprecated.h.
|
related |
Contract a tensor of rank 1 with a tensor of rank 1 and return the result.
Definition at line 439 of file tensor_deprecated.h.
|
related |
The cross product of one vector in 2d. This is just a rotation by 90 degrees.
Definition at line 451 of file tensor_deprecated.h.
|
related |
The cross product of 2 vectors in 3d.
Definition at line 458 of file tensor_deprecated.h.
|
related |
Form the outer product of two tensors.
Definition at line 467 of file tensor_deprecated.h.
|
related |
Multiply a Tensor<1,dim,Number> with a Number.
Definition at line 475 of file tensor_deprecated.h.
|
related |
Multiply a Tensor<1,dim,Number> with a Number.
Definition at line 484 of file tensor_deprecated.h.
|
related |
Definition at line 494 of file tensor_deprecated.h.
|
related |
Definition at line 501 of file tensor_deprecated.h.
|
static |
Provide a way to get the dimension of an object without explicit knowledge of it's data type. Implementation is this way instead of providing a function dimension()
because now it is possible to get the dimension at compile time without the expansion and preevaluation of an inlined function; the compiler may therefore produce more efficient code and you may use this value to declare other data types.
|
static |
|
static |
Number of independent components of a tensor of current rank. This is dim times the number of independent components of each sub-tensor.