Reference documentation for deal.II version 9.2.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\}}\)
Classes | Namespaces | Functions
tensor.h File Reference
#include <deal.II/base/config.h>
#include <deal.II/base/exceptions.h>
#include <deal.II/base/numbers.h>
#include <deal.II/base/std_cxx14/utility.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 <deal.II/lac/lapack_full_matrix.h>
#include <adolc/adouble.h>
#include <cmath>
#include <ostream>
#include <vector>

Go to the source code of this file.

Classes

class  Tensor< 0, dim, Number >
 
class  Tensor< rank_, dim, Number >
 

Namespaces

 internal
 
 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, 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 >
constexpr 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 >
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, 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 >
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<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 > &tensor)
 
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 1661 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 1682 of file tensor.h.

◆ operator*() [1/6]

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 CUDA device code.

Definition at line 1709 of file tensor.h.

◆ operator*() [2/6]

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 CUDA device code.

Definition at line 1729 of file tensor.h.

◆ operator*() [3/6]

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 CUDA device code.

Definition at line 1749 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 CUDA device code.

Definition at line 1770 of file tensor.h.

◆ operator+() [1/2]

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

Add two tensors of rank 0.

Note
This function can also be used in CUDA device code.

Definition at line 1786 of file tensor.h.

◆ operator-() [1/2]

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

Subtract two tensors of rank 0.

Note
This function can also be used in CUDA device code.

Definition at line 1803 of file tensor.h.

◆ operator*() [4/6]

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 CUDA device code.

Definition at line 1828 of file tensor.h.

◆ operator*() [5/6]

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 CUDA device code.

Definition at line 1856 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 CUDA device code.

Definition at line 1927 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 CUDA device code.

Definition at line 1945 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 CUDA device code.

Definition at line 1969 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 1989 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 2018 of file tensor.h.

◆ operator*() [6/6]

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 
)

Definition at line 2070 of file tensor.h.

◆ contract()

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, typename ProductType< Number, OtherNumber >::type >::tensor_type contract ( const Tensor< rank_1, dim, Number > &  src1,
const Tensor< rank_2, dim, OtherNumber > &  src2 
)
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

contract<0, 2>(t1, t2);
Note
The position of the index is counted from 0, i.e., \(0\le\text{index}_i<\text{range}_i\).
In case the contraction yields a tensor of rank 0 the scalar number is returned as an unwrapped number type.
Author
Matthias Maier, 2015

Definition at line 2127 of file tensor.h.

◆ double_contract()

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, typename ProductType< Number, OtherNumber >::type >::tensor_type double_contract ( const Tensor< rank_1, dim, Number > &  src1,
const Tensor< rank_2, dim, OtherNumber > &  src2 
)
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

double_contract<0, 2, 1, 0>(t1, t2);
Note
The position of the index is counted from 0, i.e., \(0\le\text{index}_i<\text{range}_i\).
In case the contraction yields a tensor of rank 0 the scalar number is returned as an unwrapped number type.
Author
Matthias Maier, 2015

Definition at line 2201 of file tensor.h.

◆ scalar_product()

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

Author
Matthias Maier, 2015

Definition at line 2281 of file tensor.h.

◆ contract3()

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

Note
Each of the three input tensors can be either a Tensor or SymmetricTensor.
Author
Matthias Maier, 2015, Jean-Paul Pelteret 2017

Definition at line 2319 of file tensor.h.

◆ outer_product()

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

The outer product of two tensors of rank_1 and rank_2: Returns a tensor of rank \((\text{rank}_1 + \text{rank}_2)\):

\[ \text{result}_{i_1,\ldots,i_{r1},j_1,\ldots,j_{r2}} = \text{left}_{i_1,\ldots,i_{r1}}\,\text{right}_{j_1,\ldots,j_{r2}.} \]

Author
Matthias Maier, 2015

Definition at line 2349 of file tensor.h.

◆ cross_product_2d()

template<int dim, typename Number >
constexpr Tensor< 1, dim, Number > cross_product_2d ( const Tensor< 1, dim, Number > &  src)
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).

Author
Guido Kanschat, 2001

Definition at line 2381 of file tensor.h.

◆ cross_product_3d()

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

Author
Guido Kanschat, 2001

Definition at line 2407 of file tensor.h.

◆ determinant() [1/2]

template<int dim, typename Number >
constexpr Number determinant ( const Tensor< 2, dim, Number > &  t)
inlineconstexpr

Compute the determinant of a tensor or rank 2.

Author
Wolfgang Bangerth, 2009

Definition at line 2442 of file tensor.h.

◆ determinant() [2/2]

template<typename Number >
constexpr Number determinant ( const Tensor< 2, 1, Number > &  t)
constexpr

Specialization for dim==1.

Definition at line 2470 of file tensor.h.

◆ trace()

template<int dim, typename Number >
constexpr Number trace ( const Tensor< 2, dim, Number > &  d)
inlineconstexpr

Compute and return the trace of a tensor of rank 2, i.e. the sum of its diagonal entries.

Author
Wolfgang Bangerth, 2001

Definition at line 2485 of file tensor.h.

◆ invert()

template<int dim, typename Number >
constexpr Tensor< 2, dim, Number > invert ( const Tensor< 2, dim, Number > &  )
inlineconstexpr

Compute and return the inverse of the given tensor. Since the compiler can perform the return value optimization, and since the size of the return object is known, it is acceptable to return the result by value, rather than by reference as a parameter.

Author
Wolfgang Bangerth, 2000

Definition at line 2505 of file tensor.h.

◆ transpose()

template<int dim, typename Number >
constexpr Tensor< 2, dim, Number > transpose ( const Tensor< 2, dim, Number > &  t)
inlineconstexpr

Return the transpose of the given tensor.

Author
Wolfgang Bangerth, 2002

Definition at line 2599 of file tensor.h.

◆ adjugate()

template<int dim, typename Number >
constexpr Tensor< 2, dim, Number > adjugate ( const Tensor< 2, dim, Number > &  t)
constexpr

Return the adjugate of the given tensor of rank 2. The adjugate of a tensor \(\mathbf A\) is defined as

\[ \textrm{adj}\mathbf A \dealcoloneq \textrm{det}\mathbf A \; \mathbf{A}^{-1} \; . \]

Note
This requires that the tensor is invertible.
Author
Jean-Paul Pelteret, 2016

Definition at line 2631 of file tensor.h.

◆ cofactor()

template<int dim, typename Number >
constexpr Tensor< 2, dim, Number > cofactor ( const Tensor< 2, dim, Number > &  t)
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} \; . \]

Note
This requires that the tensor is invertible.
Author
Jean-Paul Pelteret, 2016

Definition at line 2653 of file tensor.h.

◆ project_onto_orthogonal_tensors()

template<int dim, typename Number >
Tensor< 2, dim, Number > project_onto_orthogonal_tensors ( const Tensor< 2, dim, Number > &  tensor)

Return the nearest orthogonal matrix by combining the products of the SVD decomposition: \(\mathbf U \mathbf{V}^T\), where \(\mathbf U\) and \(\mathbf V\) are computed from the SVD decomposition: \(\mathbf U \mathbf S \mathbf V^T\), effectively replacing \(\mathbf S\) with the identity matrix.

Parameters
tensorThe tensor which to find the closest orthogonal tensor to.

Definition at line 2671 of file tensor.h.

◆ l1_norm()

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

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

Author
Wolfgang Bangerth, 2012

Definition at line 2705 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).

Author
Wolfgang Bangerth, 2012

Definition at line 2732 of file tensor.h.