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
Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Member Functions | Private Attributes | Friends | Related Symbols | List of all members
Tensor< rank_, dim, Number > Class Template Reference

#include <deal.II/base/tensor.h>

Inherited by Point< spacedim >, Point< dim, VectorizedArray< Number > >, Point< spacedim, VectorizedArray< Number > >, Point< dim, VectorizedArrayType >, Point< spacedim, VectorizedArrayType >, Point< dim, VectorizedArray< double > >, Point< spacedim, VectorizedArray< double > >, Point< 1 >, Point< spacedim, VectorizedArray< typename VectorType::value_type > >, Point< 3, float >, Point< 3 >, Point< dim - 1, VectorizedArrayType >, Point< dim, double >, Point< spacedim, double >, Point< 2, double >, Point< 2 >, Point< 0 >, and Point< dim, Number >.

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 DEAL_II_HOST_DEVICE_ALWAYS_INLINE Tensor ()
 
constexpr Tensor (const array_type &initializer)
 
template<typename ElementType , typename MemorySpace >
constexpr Tensor (const ArrayView< ElementType, MemorySpace > &initializer)
 
template<typename OtherNumber >
constexpr Tensor (const Tensor< rank_, dim, OtherNumber > &initializer)
 
template<typename OtherNumber >
constexpr 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
 
constexpr value_typeoperator[] (const unsigned int i)
 
constexpr const value_typeoperator[] (const unsigned int i) const
 
constexpr const Number & operator[] (const TableIndices< rank_ > &indices) const
 
constexpr Number & operator[] (const TableIndices< rank_ > &indices)
 
Number * begin_raw ()
 
const Number * begin_raw () const
 
Number * end_raw ()
 
const Number * end_raw () const
 
template<typename OtherNumber >
constexpr Tensoroperator= (const Tensor< rank_, dim, OtherNumber > &rhs)
 
constexpr Tensoroperator= (const Number &d) &
 
constexpr Tensoroperator= (const Number &d) &&=delete
 
template<typename OtherNumber >
constexpr bool operator== (const Tensor< rank_, dim, OtherNumber > &) const
 
template<typename OtherNumber >
constexpr bool operator!= (const Tensor< rank_, dim, OtherNumber > &) const
 
template<typename OtherNumber >
constexpr Tensoroperator+= (const Tensor< rank_, dim, OtherNumber > &)
 
template<typename OtherNumber >
constexpr Tensoroperator-= (const Tensor< rank_, dim, OtherNumber > &)
 
template<typename OtherNumber >
constexpr Tensoroperator*= (const OtherNumber &factor)
 
template<typename OtherNumber >
constexpr Tensoroperator/= (const OtherNumber &factor)
 
constexpr Tensor operator- () const
 
constexpr void clear ()
 
numbers::NumberTraits< Number >::real_type norm () const
 
constexpr numbers::NumberTraits< Number >::real_type norm_square () const
 
template<typename OtherNumber >
void unroll (Vector< OtherNumber > &result) const
 
template<class Iterator >
void unroll (const Iterator begin, const Iterator end) const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Static Public Member Functions

static constexpr unsigned int component_to_unrolled_index (const TableIndices< rank_ > &indices)
 
static constexpr 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 Iterator >
Iterator unroll_recursion (const Iterator current, const Iterator end) const
 
template<typename ArrayLike , std::size_t... Indices>
constexpr Tensor (const ArrayLike &initializer, std::index_sequence< Indices... >)
 

Private Attributes

Tensor< rank_ - 1, dim, Number > values [(dim !=0) ? dim :1]
 

Friends

template<int , int , typename >
class Tensor
 
class Point< dim, Number >
 

Related Symbols

(Note that these are not member symbols.)

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

Detailed Description

template<int rank_, int dim, typename Number>
class Tensor< rank_, dim, Number >

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.

The rank of a tensor specifies which types of physical quantities it can represent:

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

Template Parameters
rank_An integer that denotes the rank of this tensor. A specialization of this class exists for rank-0 tensors.
dimAn 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.
NumberThe 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.

Definition at line 515 of file tensor.h.

Member Typedef Documentation

◆ value_type

template<int rank_, int dim, typename Number >
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>.

Definition at line 549 of file tensor.h.

◆ array_type

template<int rank_, int dim, typename Number >
using Tensor< rank_, dim, Number >::array_type = typename Tensor<rank_ - 1, dim, Number>::array_type[(dim != 0) ? dim : 1]

Declare an array type which can be used to initialize an object of this type statically. For dim == 0, its size is 1. Otherwise, it is dim.

Definition at line 555 of file tensor.h.

◆ tensor_type

template<int rank_, int dim, typename Number >
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>

Definition at line 878 of file tensor.h.

Constructor & Destructor Documentation

◆ Tensor() [1/6]

template<int rank_, int dim, typename Number >
constexpr DEAL_II_HOST_DEVICE_ALWAYS_INLINE Tensor< rank_, dim, Number >::Tensor ( )
constexpr

Constructor. Initialize all entries to zero.

Note
This function can also be used in device code.

◆ Tensor() [2/6]

template<int rank_, int dim, typename Number >
constexpr Tensor< rank_, dim, Number >::Tensor ( const array_type initializer)
explicitconstexpr

A constructor where the data is copied from a C-style array.

Note
This function can also be used in device code.

◆ Tensor() [3/6]

template<int rank_, int dim, typename Number >
template<typename ElementType , typename MemorySpace >
constexpr Tensor< rank_, dim, Number >::Tensor ( const ArrayView< ElementType, MemorySpace > &  initializer)
explicitconstexpr

A constructor where the data is copied from an ArrayView object. Obviously, the ArrayView object must represent a stretch of data of size dimrank. The sequentially ordered elements of the argument initializer are interpreted as described by unrolled_to_component_indices().

This constructor obviously requires that the ElementType type is either equal to Number, or is convertible to Number.

Note
This function can also be used in device code.

◆ Tensor() [4/6]

template<int rank_, int dim, typename Number >
template<typename OtherNumber >
constexpr Tensor< rank_, dim, Number >::Tensor ( const Tensor< rank_, dim, OtherNumber > &  initializer)
constexpr

Constructor from tensors with different underlying scalar type. This obviously requires that the OtherNumber type is convertible to Number.

Note
This function can also be used in device code.

◆ Tensor() [5/6]

template<int rank_, int dim, typename Number >
template<typename OtherNumber >
constexpr Tensor< rank_, dim, Number >::Tensor ( const Tensor< 1, dim, Tensor< rank_ - 1, dim, OtherNumber > > &  initializer)
constexpr

Constructor that converts from a "tensor of tensors".

◆ Tensor() [6/6]

template<int rank_, int dim, typename Number >
template<typename ArrayLike , std::size_t... Indices>
constexpr Tensor< rank_, dim, Number >::Tensor ( const ArrayLike &  initializer,
std::index_sequence< Indices... >   
)
constexprprivate

This constructor is for internal use. It provides a way to create constexpr constructors for Tensor<rank, dim, Number>

Note
This function can also be used in device code.

Member Function Documentation

◆ operator Tensor< 1, dim, Tensor< rank_ - 1, dim, OtherNumber > >()

template<int rank_, int dim, typename Number >
template<typename OtherNumber >
constexpr Tensor< rank_, dim, Number >::operator Tensor< 1, dim, Tensor< rank_ - 1, dim, OtherNumber > > ( ) const
constexpr

Conversion operator to tensor of tensors.

◆ operator[]() [1/4]

template<int rank_, int dim, typename Number >
constexpr value_type & Tensor< rank_, dim, Number >::operator[] ( const unsigned int  i)
constexpr

Read-Write access operator.

Note
This function can also be used in device code.

◆ operator[]() [2/4]

template<int rank_, int dim, typename Number >
constexpr const value_type & Tensor< rank_, dim, Number >::operator[] ( const unsigned int  i) const
constexpr

Read-only access operator.

Note
This function can also be used in device code.

◆ operator[]() [3/4]

template<int rank_, int dim, typename Number >
constexpr const Number & Tensor< rank_, dim, Number >::operator[] ( const TableIndices< rank_ > &  indices) const
constexpr

Read access using TableIndices indices

◆ operator[]() [4/4]

template<int rank_, int dim, typename Number >
constexpr Number & Tensor< rank_, dim, Number >::operator[] ( const TableIndices< rank_ > &  indices)
constexpr

Read and write access using TableIndices indices

◆ begin_raw() [1/2]

template<int rank_, int dim, typename Number >
Number * Tensor< rank_, dim, Number >::begin_raw ( )

Return a pointer to the first element of the underlying storage.

◆ begin_raw() [2/2]

template<int rank_, int dim, typename Number >
const Number * Tensor< rank_, dim, Number >::begin_raw ( ) const

Return a const pointer to the first element of the underlying storage.

◆ end_raw() [1/2]

template<int rank_, int dim, typename Number >
Number * Tensor< rank_, dim, Number >::end_raw ( )

Return a pointer to the element past the end of the underlying storage.

◆ end_raw() [2/2]

template<int rank_, int dim, typename Number >
const Number * Tensor< rank_, dim, Number >::end_raw ( ) const

Return a pointer to the element past the end of the underlying storage.

◆ operator=() [1/3]

template<int rank_, int dim, typename Number >
template<typename OtherNumber >
constexpr Tensor & Tensor< rank_, dim, Number >::operator= ( const Tensor< rank_, dim, OtherNumber > &  rhs)
constexpr

Assignment operator from tensors with different underlying scalar type. This obviously requires that the OtherNumber type is convertible to Number.

Note
This function can also be used in device code.

◆ operator=() [2/3]

template<int rank_, int dim, typename Number >
constexpr Tensor & Tensor< rank_, dim, Number >::operator= ( const Number &  d) &
constexpr

This operator assigns a scalar to a tensor. To avoid confusion with what exactly it means to assign a scalar value to a tensor, zero is the only value allowed for d, allowing the intuitive notation t=0 to reset all elements of the tensor to zero.

◆ operator=() [3/3]

template<int rank_, int dim, typename Number >
constexpr Tensor & Tensor< rank_, dim, Number >::operator= ( const Number &  d) &&
constexprdelete

Assign a scalar to the current object. This overload is used for rvalue references; because it does not make sense to assign something to a temporary, the function is deleted.

◆ operator==()

template<int rank_, int dim, typename Number >
template<typename OtherNumber >
constexpr bool Tensor< rank_, dim, Number >::operator== ( const Tensor< rank_, dim, OtherNumber > &  ) const
constexpr

Test for equality of two tensors.

◆ operator!=()

template<int rank_, int dim, typename Number >
template<typename OtherNumber >
constexpr bool Tensor< rank_, dim, Number >::operator!= ( const Tensor< rank_, dim, OtherNumber > &  ) const
constexpr

Test for inequality of two tensors.

◆ operator+=()

template<int rank_, int dim, typename Number >
template<typename OtherNumber >
constexpr Tensor & Tensor< rank_, dim, Number >::operator+= ( const Tensor< rank_, dim, OtherNumber > &  )
constexpr

Add another tensor.

Note
This function can also be used in device code.

◆ operator-=()

template<int rank_, int dim, typename Number >
template<typename OtherNumber >
constexpr Tensor & Tensor< rank_, dim, Number >::operator-= ( const Tensor< rank_, dim, OtherNumber > &  )
constexpr

Subtract another tensor.

Note
This function can also be used in device code.

◆ operator*=()

template<int rank_, int dim, typename Number >
template<typename OtherNumber >
constexpr Tensor & Tensor< rank_, dim, Number >::operator*= ( const OtherNumber &  factor)
constexpr

Scale the tensor by factor, i.e. multiply all components by factor.

Note
This function can also be used in device code.

◆ operator/=()

template<int rank_, int dim, typename Number >
template<typename OtherNumber >
constexpr Tensor & Tensor< rank_, dim, Number >::operator/= ( const OtherNumber &  factor)
constexpr

Scale the vector by 1/factor.

Note
This function can also be used in device code.

◆ operator-()

template<int rank_, int dim, typename Number >
constexpr Tensor Tensor< rank_, dim, Number >::operator- ( ) const
constexpr

Unary minus operator. Negate all entries of a tensor.

Note
This function can also be used in device code.

◆ clear()

template<int rank_, int dim, typename Number >
constexpr void Tensor< rank_, dim, Number >::clear ( )
constexpr

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.

◆ norm()

template<int rank_, int dim, typename Number >
numbers::NumberTraits< Number >::real_type Tensor< rank_, dim, Number >::norm ( ) const

Return the Frobenius-norm of a tensor, i.e. the square root of the sum of the absolute squares of all entries. For the present case of rank-1 tensors, this equals the usual l2 norm of the vector.

Note
This function can also be used in device code.

◆ norm_square()

template<int rank_, int dim, typename Number >
constexpr numbers::NumberTraits< Number >::real_type Tensor< rank_, dim, Number >::norm_square ( ) const
constexpr

Return the square of the Frobenius-norm of a tensor, i.e. the sum of the absolute squares of all entries.

Note
This function can also be used in device code.

◆ unroll() [1/2]

template<int rank_, int dim, typename Number >
template<typename OtherNumber >
void Tensor< rank_, dim, Number >::unroll ( Vector< OtherNumber > &  result) const

Fill a vector with all tensor elements.

This function unrolls all tensor entries into a single, linearly numbered vector. As usual in C++, the rightmost index of the tensor marches fastest.

Deprecated:
Use the more general function that takes a pair of iterators instead.

◆ unroll() [2/2]

template<int rank_, int dim, typename Number >
template<class Iterator >
void Tensor< rank_, dim, Number >::unroll ( const Iterator  begin,
const Iterator  end 
) const

Fill a range with all tensor elements.

This function unrolls all tensor entries into a single, linearly numbered sequence. The order of the elements is the one given by component_to_unrolled_index().

The template type Number must be convertible to the type of *begin.

◆ component_to_unrolled_index()

template<int rank_, int dim, typename Number >
static constexpr unsigned int Tensor< rank_, dim, Number >::component_to_unrolled_index ( const TableIndices< rank_ > &  indices)
staticconstexpr

Return an unrolled index in the range \([0,\text{dim}^{\text{rank}}-1]\) for the element of the tensor indexed by the argument to the function.

◆ unrolled_to_component_indices()

template<int rank_, int dim, typename Number >
static constexpr TableIndices< rank_ > Tensor< rank_, dim, Number >::unrolled_to_component_indices ( const unsigned int  i)
staticconstexpr

Opposite of component_to_unrolled_index: For an index in the range \([0, \text{dim}^{\text{rank}}-1]\), return which set of indices it would correspond to.

◆ memory_consumption()

template<int rank_, int dim, typename Number >
static constexpr std::size_t Tensor< rank_, dim, Number >::memory_consumption ( )
staticconstexpr

Determine an estimate for the memory consumption (in bytes) of this object.

◆ serialize()

template<int rank_, int dim, typename Number >
template<class Archive >
void Tensor< rank_, dim, Number >::serialize ( Archive &  ar,
const unsigned int  version 
)

Read or write the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.

◆ unroll_recursion()

template<int rank_, int dim, typename Number >
template<typename Iterator >
Iterator Tensor< rank_, dim, Number >::unroll_recursion ( const Iterator  current,
const Iterator  end 
) const
private

Internal helper function for unroll.

Friends And Related Symbol Documentation

◆ Tensor

template<int rank_, int dim, typename Number >
template<int , int , typename >
friend class Tensor
friend

Definition at line 907 of file tensor.h.

◆ Point< dim, Number >

template<int rank_, int dim, typename Number >
friend class Point< dim, Number >
friend

Definition at line 907 of file tensor.h.

◆ sum()

template<int rank, int dim, typename Number >
Tensor< rank, dim, Number > sum ( const Tensor< rank, dim, Number > &  local,
const MPI_Comm  mpi_communicator 
)
related

Perform an MPI sum of the entries of a tensor.

◆ operator<<() [1/2]

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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 rank_1, int rank_2, int dim, typename Number , typename OtherNumber , typename = std::enable_if_t<rank_1 >= 1 && rank_2>
template<int dim, typename Number >
Number l1_norm ( const Tensor< 2, dim, Number > &  t)
related

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

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.

Member Data Documentation

◆ dimension

template<int rank_, int dim, typename Number >
constexpr unsigned int Tensor< rank_, dim, Number >::dimension = dim
staticconstexpr

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.

Definition at line 530 of file tensor.h.

◆ rank

template<int rank_, int dim, typename Number >
constexpr unsigned int Tensor< rank_, dim, Number >::rank = rank_
staticconstexpr

Publish the rank of this tensor to the outside world.

Definition at line 535 of file tensor.h.

◆ n_independent_components

template<int rank_, int dim, typename Number >
constexpr unsigned int Tensor< rank_, dim, Number >::n_independent_components
staticconstexpr
Initial value:
=
Tensor<rank_ - 1, dim>::n_independent_components * dim
static constexpr unsigned int n_independent_components
Definition tensor.h:541

Number of independent components of a tensor of current rank. This is dim times the number of independent components of each sub-tensor.

Definition at line 541 of file tensor.h.

◆ values

template<int rank_, int dim, typename Number >
Tensor<rank_ - 1, dim, Number> Tensor< rank_, dim, Number >::values[(dim !=0) ? dim :1]
private

Array of tensors holding the subelements.

Definition at line 884 of file tensor.h.


The documentation for this class was generated from the following files: