Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Member Functions | Private Attributes | Friends | Related Functions | List of all members
Tensor< rank_, dim, Number > Class Template Reference

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

Inheritance diagram for Tensor< rank_, dim, Number >:
[legend]

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_typeoperator[] (const unsigned int i)
 
constexpr const value_typeoperator[] (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 >
Tensoroperator= (const Tensor< rank_, dim, OtherNumber > &rhs)
 
Tensoroperator= (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 >
Tensoroperator+= (const Tensor< rank_, dim, OtherNumber > &)
 
template<typename OtherNumber >
Tensoroperator-= (const Tensor< rank_, dim, OtherNumber > &)
 
template<typename OtherNumber >
Tensoroperator*= (const OtherNumber &factor)
 
template<typename OtherNumber >
Tensoroperator/= (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)
 

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.

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 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.
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.
Author
Wolfgang Bangerth, 1998-2005, Matthias Maier, 2015

Definition at line 90 of file mpi.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 427 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 434 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 666 of file tensor.h.

Constructor & Destructor Documentation

◆ Tensor() [1/4]

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

Constructor. Initialize all entries to zero.

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

◆ Tensor() [2/4]

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

Constructor, where the data is copied from a C-style array.

Definition at line 1041 of file tensor.h.

◆ Tensor() [3/4]

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

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

Definition at line 1051 of file tensor.h.

◆ Tensor() [4/4]

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

Constructor that converts from a "tensor of tensors".

Definition at line 1062 of file tensor.h.

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

Conversion operator to tensor of tensors.

Definition at line 1073 of file tensor.h.

◆ operator[]() [1/4]

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

Read-Write access operator.

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

Definition at line 1118 of file tensor.h.

◆ operator[]() [2/4]

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

Read-only access operator.

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

Definition at line 1128 of file tensor.h.

◆ operator[]() [3/4]

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

Read access using TableIndices indices

Definition at line 1137 of file tensor.h.

◆ operator[]() [4/4]

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

Read and write access using TableIndices indices

Definition at line 1149 of file tensor.h.

◆ begin_raw() [1/2]

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

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

Definition at line 1161 of file tensor.h.

◆ begin_raw() [2/2]

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

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

Definition at line 1171 of file tensor.h.

◆ end_raw() [1/2]

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

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

Definition at line 1181 of file tensor.h.

◆ end_raw() [2/2]

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

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

Definition at line 1190 of file tensor.h.

◆ operator=() [1/2]

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

◆ operator=() [2/2]

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

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.

Definition at line 1210 of file tensor.h.

◆ operator==()

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

Test for equality of two tensors.

Definition at line 1226 of file tensor.h.

◆ operator!=()

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

Test for inequality of two tensors.

Definition at line 1253 of file tensor.h.

◆ operator+=()

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

Add another tensor.

◆ operator-=()

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

Subtract another tensor.

◆ operator*=()

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

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

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

◆ operator/=()

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

Scale the vector by 1/factor.

◆ operator-()

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

Unary minus operator. Negate all entries of a tensor.

Definition at line 1305 of file tensor.h.

◆ clear()

template<int rank_, int dim, typename Number >
void Tensor< rank_, dim, Number >::clear ( )
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.

Definition at line 1437 of file tensor.h.

◆ norm()

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

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.

Definition at line 1318 of file tensor.h.

◆ norm_square()

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

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

Definition at line 1326 of file tensor.h.

◆ unroll()

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

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.

Definition at line 1340 of file tensor.h.

◆ component_to_unrolled_index()

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

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

Definition at line 1363 of file tensor.h.

◆ unrolled_to_component_indices()

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

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

Definition at line 1416 of file tensor.h.

◆ memory_consumption()

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

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

Definition at line 1446 of file tensor.h.

◆ serialize()

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

Read or write the data of this object to or from a stream for the purpose of serialization

Definition at line 1455 of file tensor.h.

◆ unroll_recursion()

template<int rank_, int dim, typename Number >
template<typename OtherNumber >
void Tensor< rank_, dim, Number >::unroll_recursion ( Vector< OtherNumber > &  result,
unsigned int &  start_index 
) const
inlineprivate

Internal helper function for unroll.

Definition at line 1353 of file tensor.h.

Friends And Related Function Documentation

◆ Tensor

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

Allow an arbitrary Tensor to access the underlying values.

Definition at line 688 of file tensor.h.

◆ Point< dim, Number >

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

Point is allowed access to the coordinates. This is supposed to improve speed.

Definition at line 694 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 1477 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.

<0,dim,Number>

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

<0,dim,Number>

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

<0,dim,Number>

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

<0,dim,Number>

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

<0,dim,Number>

Definition at line 1575 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 
)
related

Add two tensors of rank 0.

<0,dim,Number>

Definition at line 1589 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 
)
related

Subtract two tensors of rank 0.

<0,dim,Number>

Definition at line 1604 of file tensor.h.

◆ operator*() [4/6]

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

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

Definition at line 1653 of file tensor.h.

◆ operator/() [2/2]

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

Definition at line 1673 of file tensor.h.

◆ operator+() [2/2]

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

Addition of two tensors of general rank.

Template Parameters
rankThe rank of both tensors.

Definition at line 1693 of file tensor.h.

◆ operator-() [2/2]

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

Subtraction of two tensors of general rank.

Template Parameters
rankThe rank of both tensors.

Definition at line 1715 of file tensor.h.

◆ operator*() [6/6]

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

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, which does the double contraction.
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 1766 of file tensor.h.

◆ contract() [1/6]

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

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 1823 of file tensor.h.

◆ double_contract() [1/3]

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

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 1897 of file tensor.h.

◆ scalar_product()

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

Author
Matthias Maier, 2015

Definition at line 1976 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 >
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 
)
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}} \]

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 2013 of file tensor.h.

◆ outer_product() [1/4]

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

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 2043 of file tensor.h.

◆ cross_product_2d()

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

Author
Guido Kanschat, 2001

Definition at line 2075 of file tensor.h.

◆ cross_product_3d()

template<int dim, typename Number >
Tensor< 1, dim, Number > cross_product_3d ( const Tensor< 1, dim, Number > &  src1,
const Tensor< 1, dim, Number > &  src2 
)
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).

Author
Guido Kanschat, 2001

Definition at line 2100 of file tensor.h.

◆ determinant() [1/4]

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

Compute the determinant of a tensor or rank 2.

Author
Wolfgang Bangerth, 2009

Definition at line 2130 of file tensor.h.

◆ determinant() [2/4]

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

Specialization for dim==1.

Definition at line 2158 of file tensor.h.

◆ trace()

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

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 2173 of file tensor.h.

◆ invert()

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

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 2193 of file tensor.h.

◆ transpose()

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

Return the transpose of the given tensor.

Author
Wolfgang Bangerth, 2002

Definition at line 2287 of file tensor.h.

◆ adjugate()

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

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

Definition at line 2319 of file tensor.h.

◆ cofactor()

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

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

Definition at line 2341 of file tensor.h.

◆ l1_norm()

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

Return the \(l_1\) norm of the given rank-2 tensor, where \(||t||_1 = \max_j \sum_i |t_{ij}|\) (maximum of the sums over columns).

Author
Wolfgang Bangerth, 2012

Definition at line 2356 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 \(||t||_\infty = \max_i \sum_j |t_{ij}|\) (maximum of the sums over rows).

Author
Wolfgang Bangerth, 2012

Definition at line 2382 of file tensor.h.

◆ double_contract() [2/3]

template<int dim, typename Number >
Number double_contract ( const Tensor< 2, dim, Number > &  src1,
const Tensor< 2, dim, Number > &  src2 
)
related

Double contract two tensors of rank 2, thus computing the Frobenius inner product sumi,j src1[i][j]*src2[i][j].

Deprecated:
Use the double_contract() function that takes indices as template arguments and returns its result instead.

Definition at line 228 of file tensor_deprecated.h.

◆ double_contract() [3/3]

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

Deprecated:
Use the double_contract() function that takes indices as template arguments and returns its result instead.

Definition at line 239 of file tensor_deprecated.h.

◆ contract() [2/6]

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

Deprecated:
Use the contract() function that takes indices as template arguments and returns its result instead.

Definition at line 252 of file tensor_deprecated.h.

◆ contract() [3/6]

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

Deprecated:
Use the contract() function that takes indices as template arguments and returns its result instead.

Definition at line 309 of file tensor_deprecated.h.

◆ contract() [4/6]

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

Deprecated:
Use the contract() function that takes indices as template arguments and returns its result instead.

Definition at line 345 of file tensor_deprecated.h.

◆ contract() [5/6]

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

Deprecated:
Use operator* instead. It denotes a single contraction.

Definition at line 427 of file tensor_deprecated.h.

◆ contract() [6/6]

template<int dim, typename Number , typename OtherNumber >
ProductType< Number, OtherNumber >::type contract ( const Tensor< 1, dim, Number > &  src1,
const Tensor< 1, dim, OtherNumber > &  src2 
)
related

Contract a tensor of rank 1 with a tensor of rank 1 and return the result.

Deprecated:
Use operator* instead. It denotes a single contraction.

Definition at line 439 of file tensor_deprecated.h.

◆ cross_product() [1/2]

template<int dim, typename Number >
void cross_product ( Tensor< 1, dim, Number > &  dst,
const Tensor< 1, dim, Number > &  src 
)
related

The cross product of one vector in 2d. This is just a rotation by 90 degrees.

Deprecated:
Use the function cross_product_2d() that returns the value.

Definition at line 451 of file tensor_deprecated.h.

◆ cross_product() [2/2]

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

The cross product of 2 vectors in 3d.

Deprecated:
Use the function cross_product_3d() that returns the value.

Definition at line 458 of file tensor_deprecated.h.

◆ outer_product() [2/4]

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

Form the outer product of two tensors.

Deprecated:
Use the generic version that returns its result instead.

Definition at line 467 of file tensor_deprecated.h.

◆ outer_product() [3/4]

template<int dim, typename Number >
void outer_product ( Tensor< 1, dim, Number > &  dst,
const Number  src1,
const Tensor< 1, dim, Number > &  src2 
)
related

Multiply a Tensor<1,dim,Number> with a Number.

Deprecated:
Use operator* instead.

Definition at line 475 of file tensor_deprecated.h.

◆ outer_product() [4/4]

template<int dim, typename Number >
void outer_product ( Tensor< 1, dim, Number > &  dst,
const Tensor< 1, dim, Number >  src1,
const Number  src2 
)
related

Multiply a Tensor<1,dim,Number> with a Number.

Deprecated:
Use operator* instead.

Definition at line 484 of file tensor_deprecated.h.

◆ determinant() [3/4]

template<int rank, typename Number >
Number determinant ( const Tensor< rank, 1, Number > &  t)
related
Deprecated:
Do not use this function, evaluate the value manually.

Definition at line 494 of file tensor_deprecated.h.

◆ determinant() [4/4]

template<typename Number >
Number determinant ( const Tensor< 1, 1, Number > &  t)
related
Deprecated:
Do not use this function, evaluate the value manually.

Definition at line 501 of file tensor_deprecated.h.

Member Data Documentation

◆ dimension

template<int rank_, int dim, typename Number>
constexpr unsigned int Tensor< rank_, dim, Number >::dimension = dim
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.

Definition at line 408 of file tensor.h.

◆ rank

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

Publish the rank of this tensor to the outside world.

Definition at line 413 of file tensor.h.

◆ n_independent_components

template<int rank_, int dim, typename Number>
constexpr unsigned int Tensor< rank_, dim, Number >::n_independent_components
static
Initial value:
=
Tensor<rank_ - 1, dim>::n_independent_components * dim

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 419 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 672 of file tensor.h.


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