Reference documentation for deal.II version 9.5.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Classes | Namespaces | Functions
tensor.h File Reference
#include <deal.II/base/config.h>
#include <deal.II/base/exceptions.h>
#include <deal.II/base/kokkos.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 <adolc/adouble.h>
#include <cmath>
#include <ostream>

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 , std::enable_if_t< !std::is_integral< typename ProductType< Number, OtherNumber >::type >::value, int > = 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 DEAL_II_HOST_DEVICE_ALWAYS_INLINE 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 DEAL_II_HOST_DEVICE_ALWAYS_INLINE 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 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)
 

Function Documentation

◆ operator<<() [1/2]

template<int rank_, int dim, typename Number >
std::ostream & operator<< ( std::ostream &  out,
const Tensor< rank_, dim, Number > &  p 
)
inline

Output operator for tensors. Print the elements consecutively, with a space in between, two spaces between rank 1 subtensors, three between rank 2 and so on.

Definition at line 1925 of file tensor.h.

◆ operator<<() [2/2]

template<int dim, typename Number >
std::ostream & operator<< ( std::ostream &  out,
const Tensor< 0, dim, Number > &  p 
)
inline

Output operator for tensors of rank 0. Since such tensors are scalars, we simply print this one value.

Definition at line 1946 of file tensor.h.

◆ operator*() [1/5]

template<int dim, typename Number , typename Other >
constexpr ProductType< Other, Number >::type operator* ( const Other &  object,
const Tensor< 0, dim, Number > &  t 
)
inlineconstexpr

Scalar multiplication of a tensor of rank 0 with an object from the left.

This function unwraps the underlying Number stored in the Tensor and multiplies object with it.

Note
This function can also be used in device code.

Definition at line 1975 of file tensor.h.

◆ operator*() [2/5]

template<int dim, typename Number , typename Other >
constexpr ProductType< Number, Other >::type operator* ( const Tensor< 0, dim, Number > &  t,
const Other &  object 
)
inlineconstexpr

Scalar multiplication of a tensor of rank 0 with an object from the right.

This function unwraps the underlying Number stored in the Tensor and multiplies object with it.

Note
This function can also be used in device code.

Definition at line 1995 of file tensor.h.

◆ operator*() [3/5]

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 
)
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.

Note
This function can also be used in device code.

Definition at line 2015 of file tensor.h.

◆ operator/() [1/2]

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 
)
constexpr

Division of a tensor of rank 0 by a scalar number.

Note
This function can also be used in device code.

Definition at line 2036 of file tensor.h.

◆ operator+() [1/2]

template<int dim, typename Number , typename OtherNumber >
constexpr DEAL_II_HOST_DEVICE_ALWAYS_INLINE Tensor< 0, dim, typename ProductType< Number, OtherNumber >::type > operator+ ( const Tensor< 0, dim, Number > &  p,
const Tensor< 0, dim, OtherNumber > &  q 
)
constexpr

Add two tensors of rank 0.

Note
This function can also be used in device code.

Definition at line 2052 of file tensor.h.

◆ operator-() [1/2]

template<int dim, typename Number , typename OtherNumber >
constexpr DEAL_II_HOST_DEVICE_ALWAYS_INLINE Tensor< 0, dim, typename ProductType< Number, OtherNumber >::type > operator- ( const Tensor< 0, dim, Number > &  p,
const Tensor< 0, dim, OtherNumber > &  q 
)
constexpr

Subtract two tensors of rank 0.

Note
This function can also be used in device code.

Definition at line 2069 of file tensor.h.

◆ operator*() [4/5]

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 
)
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.

Note
This function can also be used in device code.

Definition at line 2094 of file tensor.h.

◆ operator*() [5/5]

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 
)
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.

Note
This function can also be used in device code.

Definition at line 2122 of file tensor.h.

◆ operator/() [2/2]

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 
)
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.

Note
This function can also be used in device code.

Definition at line 2193 of file tensor.h.

◆ operator+() [2/2]

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 
)
inlineconstexpr

Addition of two tensors of general rank.

Template Parameters
rankThe rank of both tensors.
Note
This function can also be used in device code.

Definition at line 2211 of file tensor.h.

◆ operator-() [2/2]

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 
)
inlineconstexpr

Subtraction of two tensors of general rank.

Template Parameters
rankThe rank of both tensors.
Note
This function can also be used in device code.

Definition at line 2235 of file tensor.h.

◆ schur_product() [1/2]

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 
)
inlineconstexpr

Entrywise multiplication of two tensor objects of rank 0 (i.e. the multiplication of two scalar values).

Definition at line 2255 of file tensor.h.

◆ schur_product() [2/2]

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 
)
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} \]

Template Parameters
rankThe rank of both tensors.

Definition at line 2284 of file tensor.h.

◆ l1_norm()

template<int dim, typename Number >
Number l1_norm ( const Tensor< 2, dim, Number > &  t)
inline

The dot product (single contraction) for tensors. This function 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}} \]

Note
For the Tensor class, the multiplication operator only performs a contraction over a single pair of indices. This is in contrast to the multiplication operator for SymmetricTensor, for which the corresponding operator*() performs a double contraction. The origin of the difference in how operator*() is implemented between Tensor and SymmetricTensor is that for the former, the product between two Tensor objects of same rank and dimension results in another Tensor object – that it, operator*() corresponds to the multiplicative group action within the group of tensors. On the other hand, there is no corresponding multiplicative group action with the set of symmetric tensors because, in general, the product of two symmetric tensors is a nonsymmetric tensor. As a consequence, for a mathematician, it is clear that operator*() for symmetric tensors must have a different meaning: namely the dot or scalar product that maps two symmetric tensors of rank 2 to a scalar. This corresponds to the double-dot (colon) operator whose meaning is then extended to the product of any two even-ranked symmetric tensors.
In case the contraction yields a tensor of rank 0, that is, if rank_1==rank_2==1, then a scalar number is returned as an unwrapped number type. Return the \(l_1\) norm of the given rank-2 tensor, where \(\|\mathbf T\|_1 = \max_j \sum_i |T_{ij}|\) (maximum of the sums over columns).

Definition at line 3035 of file tensor.h.

◆ linfty_norm()

template<int dim, typename Number >
Number linfty_norm ( const Tensor< 2, dim, Number > &  t)
inline

Return the \(l_\infty\) norm of the given rank-2 tensor, where \(\|\mathbf T\|_\infty = \max_i \sum_j |T_{ij}|\) (maximum of the sums over rows).

Definition at line 3061 of file tensor.h.