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 Attributes | Private Member Functions | Private Attributes | Friends | List of all members
Tensor< 0, dim, Number > Class Template Reference

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

Public Types

using real_type = typename numbers::NumberTraits< Number >::real_type
 
using value_type = Number
 
using array_type = Number
 
using tensor_type = Number
 

Public Member Functions

constexpr Tensor ()
 
template<typename OtherNumber >
constexpr Tensor (const Tensor< 0, dim, OtherNumber > &initializer)
 
template<typename OtherNumber >
constexpr Tensor (const OtherNumber &initializer)
 
Number * begin_raw ()
 
const Number * begin_raw () const
 
Number * end_raw ()
 
const Number * end_raw () const
 
constexpr operator Number & ()
 
constexpr operator const Number & () const
 
template<typename OtherNumber >
constexpr Tensoroperator= (const Tensor< 0, dim, OtherNumber > &rhs)
 
template<typename OtherNumber >
constexpr Tensoroperator= (const OtherNumber &d) &
 
template<typename OtherNumber >
constexpr Tensoroperator= (const OtherNumber &d) &&=delete
 
template<typename OtherNumber >
constexpr bool operator== (const Tensor< 0, dim, OtherNumber > &rhs) const
 
template<typename OtherNumber >
constexpr bool operator!= (const Tensor< 0, dim, OtherNumber > &rhs) const
 
template<typename OtherNumber >
constexpr Tensoroperator+= (const Tensor< 0, dim, OtherNumber > &rhs)
 
template<typename OtherNumber >
constexpr Tensoroperator-= (const Tensor< 0, dim, OtherNumber > &rhs)
 
template<typename OtherNumber >
constexpr Tensoroperator*= (const OtherNumber &factor)
 
template<typename OtherNumber >
constexpr Tensoroperator/= (const OtherNumber &factor)
 
constexpr Tensor operator- () const
 
constexpr void clear ()
 
real_type norm () const
 
constexpr real_type norm_square () 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 Attributes

static constexpr unsigned int dimension = dim
 
static constexpr unsigned int rank = 0
 
static constexpr unsigned int n_independent_components = 1
 

Private Member Functions

template<typename Iterator >
Iterator unroll_recursion (const Iterator current, const Iterator end) const
 

Private Attributes

Number value
 

Friends

template<int , int , typename >
class Tensor
 

Detailed Description

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

This class is a specialized version of the Tensor<rank,dim,Number> class. It handles tensors of rank zero, i.e. scalars. The second template argument dim is ignored.

This class exists because in some cases we want to construct objects of type Tensor<spacedim-dim,dim,Number>, which should expand to scalars, vectors, matrices, etc, depending on the values of the template arguments dim and spacedim. We therefore need a class that acts as a scalar (i.e. Number) for all purposes but is part of the Tensor template family.

Template Parameters
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. Since the current object is a rank-0 tensor (a scalar), this template argument has no meaning for this class.
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 92 of file tensor.h.

Member Typedef Documentation

◆ real_type

template<int dim, typename Number >
using Tensor< 0, dim, Number >::real_type = typename numbers::NumberTraits<Number>::real_type

Declare a type that has holds real-valued numbers with the same precision as the template argument to this class. For std::complex<number>, this corresponds to type number, and it is equal to Number for all other cases. See also the respective field in Vector<Number>.

This alias is used to represent the return type of norms.

Definition at line 126 of file tensor.h.

◆ value_type

template<int dim, typename Number >
using Tensor< 0, dim, Number >::value_type = Number

Type of objects encapsulated by this container and returned by operator[](). This is a scalar number type for a rank 0 tensor.

Definition at line 132 of file tensor.h.

◆ array_type

template<int dim, typename Number >
using Tensor< 0, dim, Number >::array_type = Number

Declare an array type which can be used to initialize an object of this type statically. In case of a tensor of rank 0 this is just the scalar number type Number.

Definition at line 139 of file tensor.h.

◆ tensor_type

template<int dim, typename Number >
using Tensor< 0, dim, Number >::tensor_type = Number

Internal type declaration that is used to specialize the return type of operator[]() for Tensor<1,dim,Number>

Definition at line 419 of file tensor.h.

Constructor & Destructor Documentation

◆ Tensor() [1/3]

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

Constructor. Set to zero.

Note
This function can also be used in device code.

◆ Tensor() [2/3]

template<int dim, typename Number >
template<typename OtherNumber >
constexpr Tensor< 0, dim, Number >::Tensor ( const Tensor< 0, 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() [3/3]

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

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

Note
This function can also be used in device code.

Member Function Documentation

◆ begin_raw() [1/2]

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

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

Deprecated:
This function suggests that the elements of a Tensor object are stored as a contiguous array, but this is not in fact true and one should not pretend that this so. As a consequence, this function is deprecated.

◆ begin_raw() [2/2]

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

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

Deprecated:
This function suggests that the elements of a Tensor object are stored as a contiguous array, but this is not in fact true and one should not pretend that this so. As a consequence, this function is deprecated.

◆ end_raw() [1/2]

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

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

Deprecated:
This function suggests that the elements of a Tensor object are stored as a contiguous array, but this is not in fact true and one should not pretend that this so. As a consequence, this function is deprecated.

◆ end_raw() [2/2]

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

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

Deprecated:
This function suggests that the elements of a Tensor object are stored as a contiguous array, but this is not in fact true and one should not pretend that this so. As a consequence, this function is deprecated.

◆ operator Number &()

template<int dim, typename Number >
constexpr Tensor< 0, dim, Number >::operator Number & ( )
constexpr

Return a reference to the encapsulated Number object. Since rank-0 tensors are scalars, this is a natural operation.

This is the non-const conversion operator that returns a writable reference.

Note
This function can also be used in device code.

◆ operator const Number &()

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

Return a reference to the encapsulated Number object. Since rank-0 tensors are scalars, this is a natural operation.

This is the const conversion operator that returns a read-only reference.

Note
This function can also be used in device code.

◆ operator=() [1/3]

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

Assignment 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 dim, typename Number >
template<typename OtherNumber >
constexpr Tensor & Tensor< 0, dim, Number >::operator= ( const OtherNumber &  d) &
constexpr

This operator assigns a scalar to a tensor. This obviously requires that the OtherNumber type is convertible to Number.

Note
This function can also be used in device code.

◆ operator=() [3/3]

template<int dim, typename Number >
template<typename OtherNumber >
constexpr Tensor & Tensor< 0, dim, Number >::operator= ( const OtherNumber &  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 dim, typename Number >
template<typename OtherNumber >
constexpr bool Tensor< 0, dim, Number >::operator== ( const Tensor< 0, dim, OtherNumber > &  rhs) const
constexpr

Test for equality of two tensors.

◆ operator!=()

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

Test for inequality of two tensors.

◆ operator+=()

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

Add another scalar.

Note
This function can also be used in device code.

◆ operator-=()

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

Subtract another scalar.

Note
This function can also be used in device code.

◆ operator*=()

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

Multiply the scalar with a factor.

Note
This function can also be used in device code.

◆ operator/=()

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

Divide the scalar by factor.

Note
This function can also be used in device code.

◆ operator-()

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

Tensor with inverted entries.

Note
This function can also be used in device code.

◆ clear()

template<int dim, typename Number >
constexpr void Tensor< 0, 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 dim, typename Number >
real_type Tensor< 0, 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.

◆ norm_square()

template<int dim, typename Number >
constexpr real_type Tensor< 0, 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()

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

Fill a range with all tensor elements. Since this type of Tensor only has one entry this just copies the value of this tensor into *begin.

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

◆ serialize()

template<int dim, typename Number >
template<class Archive >
void Tensor< 0, 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 dim, typename Number >
template<typename Iterator >
Iterator Tensor< 0, 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 dim, typename Number >
template<int , int , typename >
friend class Tensor
friend

Definition at line 436 of file tensor.h.

Member Data Documentation

◆ dimension

template<int dim, typename Number >
constexpr unsigned int Tensor< 0, 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 106 of file tensor.h.

◆ rank

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

Publish the rank of this tensor to the outside world.

Definition at line 111 of file tensor.h.

◆ n_independent_components

template<int dim, typename Number >
constexpr unsigned int Tensor< 0, dim, Number >::n_independent_components = 1
staticconstexpr

Number of independent components of a tensor of rank 0.

Definition at line 116 of file tensor.h.

◆ value

template<int dim, typename Number >
Number Tensor< 0, dim, Number >::value
private

The value of this scalar object.

Definition at line 425 of file tensor.h.


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