Reference documentation for deal.II version 9.0.0
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

typedef numbers::NumberTraits< Number >::real_type real_type
 
typedef Number value_type
 
typedef Number array_type
 
typedef Number tensor_type
 

Public Member Functions

 Tensor ()
 
template<typename OtherNumber >
 Tensor (const Tensor< 0, dim, OtherNumber > &initializer)
 
template<typename OtherNumber >
 Tensor (const OtherNumber &initializer)
 
Number * begin_raw ()
 
const Number * begin_raw () const
 
Number * end_raw ()
 
const Number * end_raw () const
 
 operator Number & ()
 
 operator const Number & () const
 
template<typename OtherNumber >
Tensoroperator= (const Tensor< 0, dim, OtherNumber > &rhs)
 
template<typename OtherNumber >
Tensoroperator= (const OtherNumber &d)
 
template<typename OtherNumber >
bool operator== (const Tensor< 0, dim, OtherNumber > &rhs) const
 
template<typename OtherNumber >
bool operator!= (const Tensor< 0, dim, OtherNumber > &rhs) const
 
template<typename OtherNumber >
Tensoroperator+= (const Tensor< 0, dim, OtherNumber > &rhs)
 
template<typename OtherNumber >
Tensoroperator-= (const Tensor< 0, dim, OtherNumber > &rhs)
 
template<typename OtherNumber >
Tensoroperator*= (const OtherNumber &factor)
 
template<typename OtherNumber >
Tensoroperator/= (const OtherNumber &factor)
 
Tensor operator- () const
 
void clear ()
 
real_type norm () const
 
real_type norm_square () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Static Public Attributes

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

Private Member Functions

template<typename OtherNumber >
void unroll_recursion (Vector< OtherNumber > &result, unsigned int &start_index) 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.
Author
Wolfgang Bangerth, 2009, Matthias Maier, 2015

Definition at line 91 of file tensor.h.

Member Typedef Documentation

◆ real_type

template<int dim, typename Number >
typedef numbers::NumberTraits<Number>::real_type Tensor< 0, dim, 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 typedef is used to represent the return type of norms.

Definition at line 122 of file tensor.h.

◆ value_type

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

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

◆ array_type

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

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

◆ tensor_type

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

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

Definition at line 312 of file tensor.h.

Constructor & Destructor Documentation

◆ Tensor() [1/3]

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

Constructor. Set to zero.

See also
CUDAWrappers

Definition at line 718 of file tensor.h.

◆ Tensor() [2/3]

template<int dim, typename Number >
template<typename OtherNumber >
Tensor< 0, dim, Number >::Tensor ( const Tensor< 0, 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 740 of file tensor.h.

◆ Tensor() [3/3]

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

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

Definition at line 730 of file tensor.h.

Member Function Documentation

◆ begin_raw() [1/2]

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

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

Definition at line 750 of file tensor.h.

◆ begin_raw() [2/2]

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

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

Definition at line 760 of file tensor.h.

◆ end_raw() [1/2]

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

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

Definition at line 770 of file tensor.h.

◆ end_raw() [2/2]

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

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

Definition at line 780 of file tensor.h.

◆ operator Number &()

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

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.

See also
CUDAWrappers

Definition at line 789 of file tensor.h.

◆ operator const Number &()

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

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.

See also
CUDAWrappers

Definition at line 801 of file tensor.h.

◆ operator=() [1/2]

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

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

◆ operator=() [2/2]

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

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

◆ operator==()

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

Test for equality of two tensors.

Definition at line 845 of file tensor.h.

◆ operator!=()

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

Test for inequality of two tensors.

Definition at line 860 of file tensor.h.

◆ operator+=()

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

Add another scalar

◆ operator-=()

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

Subtract another scalar.

◆ operator*=()

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

Multiply the scalar with a factor.

See also
CUDAWrappers

◆ operator/=()

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

Divide the scalar by factor.

◆ operator-()

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

Tensor with inverted entries.

Definition at line 908 of file tensor.h.

◆ clear()

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

◆ norm()

template<int dim, typename Number >
Tensor< 0, dim, Number >::real_type Tensor< 0, 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 917 of file tensor.h.

◆ norm_square()

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

See also
CUDAWrappers

Definition at line 927 of file tensor.h.

◆ serialize()

template<int dim, typename Number >
template<class Archive >
void Tensor< 0, 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 963 of file tensor.h.

◆ unroll_recursion()

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

Internal helper function for unroll.

Definition at line 941 of file tensor.h.

Friends And Related Function Documentation

◆ Tensor

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

Allow an arbitrary Tensor to access the underlying values.

Definition at line 330 of file tensor.h.

Member Data Documentation

◆ dimension

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

◆ rank

template<int dim, typename Number >
const unsigned int Tensor< 0, dim, Number >::rank = 0
static

Publish the rank of this tensor to the outside world.

Definition at line 107 of file tensor.h.

◆ n_independent_components

template<int dim, typename Number >
const unsigned int Tensor< 0, dim, Number >::n_independent_components = 1
static

Number of independent components of a tensor of rank 0.

Definition at line 112 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 318 of file tensor.h.


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