Reference documentation for deal.II version 9.3.3
\(\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\}}\)
Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Member Functions | Private Attributes | Related Functions | List of all members
Point< dim, Number > Class Template Reference

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

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

Public Types

using value_type = typename Tensor< rank_ - 1, dim, double >::tensor_type
 
using array_type = typename Tensor< rank_ - 1, dim, double >::array_type[(dim !=0) ? dim :1]
 
using tensor_type = Tensor< rank_, dim, double >
 

Public Member Functions

 Point ()
 
 Point (const Tensor< 1, dim, Number > &)
 
 Point (const Number x)
 
 Point (const Number x, const Number y)
 
 Point (const Number x, const Number y, const Number z)
 
template<std::size_t dummy_dim, typename std::enable_if<(dim==dummy_dim) &&(dummy_dim !=0), int >::type = 0>
 Point (const boost::geometry::model::point< Number, dummy_dim, boost::geometry::cs::cartesian > &boost_pt)
 
Number operator() (const unsigned int index) const
 
Number & operator() (const unsigned int index)
 
template<typename OtherNumber >
Point< dim, Number > & operator= (const Tensor< 1, dim, OtherNumber > &p)
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
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 double & operator[] (const TableIndices< rank_ > &indices) const
 
constexpr double & operator[] (const TableIndices< rank_ > &indices)
 
double * begin_raw ()
 
const double * begin_raw () const
 
double * end_raw ()
 
const double * end_raw () const
 
constexpr bool operator== (const Tensor< rank_, dim, OtherNumber > &) const
 
constexpr bool operator!= (const Tensor< rank_, dim, OtherNumber > &) const
 
constexpr Tensoroperator+= (const Tensor< rank_, dim, OtherNumber > &)
 
constexpr Tensoroperator-= (const Tensor< rank_, dim, OtherNumber > &)
 
constexpr Tensoroperator*= (const OtherNumber &factor)
 
constexpr Tensoroperator/= (const OtherNumber &factor)
 
constexpr void clear ()
 
numbers::NumberTraits< double >::real_type norm () const
 
constexpr numbers::NumberTraits< double >::real_type norm_square () const
 
void unroll (Vector< OtherNumber > &result) const
 
Addition and subtraction of points.
Point< dim, Number > operator+ (const Tensor< 1, dim, Number > &) const
 
Tensor< 1, dim, Number > operator- (const Point< dim, Number > &) const
 
Point< dim, Number > operator- (const Tensor< 1, dim, Number > &) const
 
Point< dim, Number > operator- () const
 

Static Public Member Functions

static Point< dim, Number > unit_vector (const unsigned int i)
 
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
 
static constexpr unsigned int rank
 
static constexpr unsigned int n_independent_components
 

Private Member Functions

void unroll_recursion (Vector< OtherNumber > &result, unsigned int &start_index) const
 

Private Attributes

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

Related Functions

(Note that these are not member functions.)

template<int dim, typename Number >
std::ostream & operator<< (std::ostream &out, const Point< dim, Number > &p)
 
template<int dim, typename Number >
std::istream & operator>> (std::istream &in, Point< dim, Number > &p)
 
Tensor< rank, dim, double > sum (const Tensor< rank, dim, double > &local, const MPI_Comm &mpi_communicator)
 
Vector space operations on Tensor objects:
constexpr Tensor< 0, dim, typename ProductType< double, OtherNumber >::type > operator- (const Tensor< 0, dim, double > &p, const Tensor< 0, dim, OtherNumber > &q)
 
constexpr Tensor< rank, dim, typename ProductType< double, OtherNumber >::type > operator- (const Tensor< rank, dim, double > &p, const Tensor< rank, dim, OtherNumber > &q)
 
constexpr ProductType< Other, double >::type operator* (const Other &object, const Tensor< 0, dim, double > &t)
 
constexpr ProductType< double, Other >::type operator* (const Tensor< 0, dim, double > &t, const Other &object)
 
constexpr ProductType< double, OtherNumber >::type operator* (const Tensor< 0, dim, double > &src1, const Tensor< 0, dim, OtherNumber > &src2)
 
constexpr Tensor< rank, dim, typename ProductType< double, typename EnableIfScalar< OtherNumber >::type >::type > operator* (const Tensor< rank, dim, double > &t, const OtherNumber &factor)
 
constexpr Tensor< rank, dim, typename ProductType< typename EnableIfScalar< double >::type, OtherNumber >::type > operator* (const double &factor, const Tensor< rank, dim, OtherNumber > &t)
 
constexpr Tensor< 0, dim, typename ProductType< double, typename EnableIfScalar< OtherNumber >::type >::type > operator/ (const Tensor< 0, dim, double > &t, const OtherNumber &factor)
 
constexpr Tensor< rank, dim, typename ProductType< double, typename EnableIfScalar< OtherNumber >::type >::type > operator/ (const Tensor< rank, dim, double > &t, const OtherNumber &factor)
 
constexpr Tensor< 0, dim, typename ProductType< double, OtherNumber >::type > operator+ (const Tensor< 0, dim, double > &p, const Tensor< 0, dim, OtherNumber > &q)
 
constexpr Tensor< rank, dim, typename ProductType< double, OtherNumber >::type > operator+ (const Tensor< rank, dim, double > &p, const Tensor< rank, dim, OtherNumber > &q)
 
constexpr Tensor< 0, dim, typename ProductType< double, OtherNumber >::type > schur_product (const Tensor< 0, dim, double > &src1, const Tensor< 0, dim, OtherNumber > &src2)
 
constexpr Tensor< rank, dim, typename ProductType< double, OtherNumber >::type > schur_product (const Tensor< rank, dim, double > &src1, const Tensor< rank, dim, OtherNumber > &src2)
 
Output functions for Tensor objects
std::ostream & operator<< (std::ostream &out, const Tensor< rank_, dim, double > &p)
 
std::ostream & operator<< (std::ostream &out, const Tensor< 0, dim, double > &p)
 
Contraction operations and the outer product for tensor objects
OtherNumber::type::tensor_type operator* (const Tensor< rank_1, dim, double > &src1, const Tensor< rank_2, dim, OtherNumber > &src2)
 
constexpr Tensor< rank_1+rank_2-2, dim, typenameProductType< double, OtherNumber >::type >::tensor_type contract (const Tensor< rank_1, dim, double > &src1, const Tensor< rank_2, dim, OtherNumber > &src2)
 
constexpr Tensor< rank_1+rank_2-4, dim, typenameProductType< double, OtherNumber >::type >::tensor_type double_contract (const Tensor< rank_1, dim, double > &src1, const Tensor< rank_2, dim, OtherNumber > &src2)
 
constexpr ProductType< double, OtherNumber >::type scalar_product (const Tensor< rank, dim, double > &left, const Tensor< rank, dim, OtherNumber > &right)
 
constexpr ProductType< T1, typenameProductType< 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)
 
constexpr Tensor< rank_1+rank_2, dim, typename ProductType< double, OtherNumber >::type > outer_product (const Tensor< rank_1, dim, double > &src1, const Tensor< rank_2, dim, OtherNumber > &src2)
 
Special operations on tensors of rank 1
constexpr Tensor< 1, dim, double > cross_product_2d (const Tensor< 1, dim, double > &src)
 
constexpr Tensor< 1, dim, typename ProductType< Number1, Number2 >::type > cross_product_3d (const Tensor< 1, dim, Number1 > &src1, const Tensor< 1, dim, Number2 > &src2)
 
Special operations on tensors of rank 2
constexpr double determinant (const Tensor< 2, dim, double > &t)
 
constexpr double determinant (const Tensor< 2, 1, double > &t)
 
constexpr double determinant (const Tensor< 2, 2, double > &t)
 
constexpr double determinant (const Tensor< 2, 3, double > &t)
 
constexpr double trace (const Tensor< 2, dim, double > &d)
 
constexpr Tensor< 2, dim, double > invert (const Tensor< 2, dim, double > &)
 
constexpr Tensor< 2, dim, double > transpose (const Tensor< 2, dim, double > &t)
 
constexpr Tensor< 2, dim, double > adjugate (const Tensor< 2, dim, double > &t)
 
constexpr Tensor< 2, dim, double > cofactor (const Tensor< 2, dim, double > &t)
 
Tensor< 2, dim, double > project_onto_orthogonal_tensors (const Tensor< 2, dim, double > &A)
 
double l1_norm (const Tensor< 2, dim, double > &t)
 
double linfty_norm (const Tensor< 2, dim, double > &t)
 

Multiplication and scaling of points. Dot products. Norms.

template<typename OtherNumber >
Point< dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > operator/ (const OtherNumber) const
 
Number operator* (const Tensor< 1, dim, Number > &p) const
 
numbers::NumberTraits< Number >::real_type square () const
 
numbers::NumberTraits< Number >::real_type distance (const Point< dim, Number > &p) const
 
numbers::NumberTraits< Number >::real_type distance_square (const Point< dim, Number > &p) const
 
template<typename OtherNumber >
Point< dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > operator* (const OtherNumber) const
 

Detailed Description

template<int dim, typename Number = double>
class Point< dim, Number >

A class that represents a point in a Cartesian space of dimension dim .

Objects of this class are used to represent points (i.e., vectors anchored at the origin) of a vector space equipped with a Cartesian coordinate system. They are, among other uses, passed to functions that operate on points in spaces of a priori fixed dimension: rather than using functions like double f(const double x) and double f(const double x, const double y), you can use double f(const Point<dim> &p) instead as it allows writing dimension independent code.

deal.II specifically uses Point objects as indicating points that are represented by Cartesian coordinates, i.e., where a point in dim space dimensions is characterized by signed distances along the axes of a coordinate system spanned by dim mutually orthogonal unit vectors (called the "coordinate axes"). This choice of representing a vector makes addition and scaling of vectors particularly simple: one only has to add or multiply each coordinate value. On the other hand, adding or scaling vectors is not nearly as simple when a vector is represented in other kinds of coordinate systems (e.g., spherical coordinate systems).

What's a Point<dim> and what is a Tensor<1,dim>?

The Point class is derived from Tensor<1,dim> and consequently shares the latter's member functions and other attributes. In fact, it has relatively few additional functions itself (the most notable exception being the distance() function to compute the Euclidean distance between two points in space), and these two classes can therefore often be used interchangeably.

Nonetheless, there are semantic differences that make us use these classes in different and well-defined contexts. Within deal.II, we use the Point class to denote points in space, i.e., for vectors (rank-1 tensors) that are anchored at the origin. On the other hand, vectors that are anchored elsewhere (and consequently do not represent points in the common usage of the word) are represented by objects of type Tensor<1,dim>. In particular, this is the case for direction vectors, normal vectors, gradients, and the differences between two points (i.e., what you get when you subtract one point from another): all of these are represented by Tensor<1,dim> objects rather than Point<dim>.

Furthermore, the Point class is only used where the coordinates of an object can be thought to possess the dimension of a length. An object that represents the weight, height, and cost of an object is neither a point nor a tensor (because it lacks the transformation properties under rotation of the coordinate system) and should consequently not be represented by either of these classes. Use an array of size 3 in this case, or the std::array class. Alternatively, as in the case of vector-valued functions, you can use objects of type Vector or std::vector.

Template Parameters
dimAn integer that denotes the dimension of the space in which a point lies. This of course equals the number of coordinates that identify a point.
NumberThe data type in which the coordinates values 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 coordinates in a different (and always scalar) type. An example would be an interval type that can store the value of a coordinate as well as its uncertainty. 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 when passed a Point object whose coordinates are stored in such a type.

Definition at line 110 of file point.h.

Member Typedef Documentation

◆ value_type

using Tensor< rank_, dim, double >::value_type = typename Tensor<rank_ - 1, dim, double >::tensor_type
inherited

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

◆ array_type

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

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

◆ tensor_type

using Tensor< rank_, dim, double >::tensor_type = Tensor<rank_, dim, double >
inherited

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

Definition at line 808 of file tensor.h.

Constructor & Destructor Documentation

◆ Point() [1/6]

template<int dim, typename Number = double>
Point< dim, Number >::Point ( )

Standard constructor. Creates an object that corresponds to the origin, i.e., all coordinates are set to zero.

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

◆ Point() [2/6]

template<int dim, typename Number = double>
Point< dim, Number >::Point ( const Tensor< 1, dim, Number > &  )
explicit

Convert a tensor to a point.

◆ Point() [3/6]

template<int dim, typename Number = double>
Point< dim, Number >::Point ( const Number  x)
explicit

Constructor for one dimensional points. This function is only implemented for dim==1 since the usage is considered unsafe for points with dim!=1 as it would leave some components of the point coordinates uninitialized.

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

◆ Point() [4/6]

template<int dim, typename Number = double>
Point< dim, Number >::Point ( const Number  x,
const Number  y 
)

Constructor for two dimensional points. This function is only implemented for dim==2 since the usage is considered unsafe for points with dim!=2 as it would leave some components of the point coordinates uninitialized (if dim>2) or would not use some arguments (if dim<2).

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

◆ Point() [5/6]

template<int dim, typename Number = double>
Point< dim, Number >::Point ( const Number  x,
const Number  y,
const Number  z 
)

Constructor for three dimensional points. This function is only implemented for dim==3 since the usage is considered unsafe for points with dim!=3 as it would leave some components of the point coordinates uninitialized (if dim>3) or would not use some arguments (if dim<3).

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

◆ Point() [6/6]

template<int dim, typename Number = double>
template<std::size_t dummy_dim, typename std::enable_if<(dim==dummy_dim) &&(dummy_dim !=0), int >::type = 0>
Point< dim, Number >::Point ( const boost::geometry::model::point< Number, dummy_dim, boost::geometry::cs::cartesian > &  boost_pt)

Convert a boost::geometry::point to a Point.

Member Function Documentation

◆ unit_vector()

template<int dim, typename Number = double>
static Point< dim, Number > Point< dim, Number >::unit_vector ( const unsigned int  i)
static

Return a unit vector in coordinate direction i, i.e., a vector that is zero in all coordinates except for a single 1 in the ith coordinate.

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

◆ operator()() [1/2]

template<int dim, typename Number = double>
Number Point< dim, Number >::operator() ( const unsigned int  index) const

Read access to the indexth coordinate.

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

◆ operator()() [2/2]

template<int dim, typename Number = double>
Number & Point< dim, Number >::operator() ( const unsigned int  index)

Read and write access to the indexth coordinate.

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

◆ operator=()

template<int dim, typename Number = double>
template<typename OtherNumber >
Point< dim, Number > & Point< dim, Number >::operator= ( const Tensor< 1, dim, OtherNumber > &  p)

Assignment operator from Tensor<1, dim, Number> with different underlying scalar type. This obviously requires that the OtherNumber type is convertible to Number.

◆ operator+()

template<int dim, typename Number = double>
Point< dim, Number > Point< dim, Number >::operator+ ( const Tensor< 1, dim, Number > &  ) const

Add an offset given as Tensor<1,dim,Number> to a point.

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

◆ operator-() [1/3]

template<int dim, typename Number = double>
Tensor< 1, dim, Number > Point< dim, Number >::operator- ( const Point< dim, Number > &  ) const

Subtract two points, i.e., obtain the vector that connects the two. As discussed in the documentation of this class, subtracting two points results in a vector anchored at one of the two points (rather than at the origin) and, consequently, the result is returned as a Tensor<1,dim> rather than as a Point<dim>.

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

◆ operator-() [2/3]

template<int dim, typename Number = double>
Point< dim, Number > Point< dim, Number >::operator- ( const Tensor< 1, dim, Number > &  ) const

Subtract a difference vector (represented by a Tensor<1,dim>) from the current point. This results in another point and, as discussed in the documentation of this class, the result is then naturally returned as a Point<dim> object rather than as a Tensor<1,dim>.

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

◆ operator-() [3/3]

template<int dim, typename Number = double>
Point< dim, Number > Point< dim, Number >::operator- ( ) const

The opposite vector.

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

◆ operator/()

template<int dim, typename Number = double>
template<typename OtherNumber >
Point< dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > Point< dim, Number >::operator/ ( const  OtherNumber) const

Divide the current point by a factor.

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

◆ operator*()

template<int dim, typename Number = double>
Number Point< dim, Number >::operator* ( const Tensor< 1, dim, Number > &  p) const

Return the scalar product of the vectors representing two points.

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

◆ square()

template<int dim, typename Number = double>
numbers::NumberTraits< Number >::real_type Point< dim, Number >::square ( ) const

Return the scalar product of this point vector with itself, i.e. the square, or the square of the norm. In case of a complex number type it is equivalent to the contraction of this point vector with a complex conjugate of itself.

Note
This function is equivalent to Tensor<rank,dim,Number>::norm_square() which returns the square of the Frobenius norm.
This function can also be used in CUDA device code.

◆ distance()

template<int dim, typename Number = double>
numbers::NumberTraits< Number >::real_type Point< dim, Number >::distance ( const Point< dim, Number > &  p) const

Return the Euclidean distance of this point to the point p, i.e. the \(l_2\) norm of the difference between the vectors representing the two points.

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

◆ distance_square()

template<int dim, typename Number = double>
numbers::NumberTraits< Number >::real_type Point< dim, Number >::distance_square ( const Point< dim, Number > &  p) const

Return the squared Euclidean distance of this point to the point p.

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

◆ serialize()

template<int dim, typename Number = double>
template<class Archive >
void Point< 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.

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

constexpr Tensor< rank_, dim, double >::operator Tensor< 1, dim, Tensor< rank_ - 1, dim, OtherNumber > > ( ) const
constexprinherited

Conversion operator to tensor of tensors.

◆ operator[]() [1/4]

constexpr value_type & Tensor< rank_, dim, double >::operator[] ( const unsigned int  i)
constexprinherited

Read-Write access operator.

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

◆ operator[]() [2/4]

constexpr const value_type & Tensor< rank_, dim, double >::operator[] ( const unsigned int  i) const
constexprinherited

Read-only access operator.

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

◆ operator[]() [3/4]

constexpr const double & Tensor< rank_, dim, double >::operator[] ( const TableIndices< rank_ > &  indices) const
constexprinherited

Read access using TableIndices indices

◆ operator[]() [4/4]

constexpr double & Tensor< rank_, dim, double >::operator[] ( const TableIndices< rank_ > &  indices)
constexprinherited

Read and write access using TableIndices indices

◆ begin_raw() [1/2]

double * Tensor< rank_, dim, double >::begin_raw ( )
inherited

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

◆ begin_raw() [2/2]

const double * Tensor< rank_, dim, double >::begin_raw ( ) const
inherited

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

◆ end_raw() [1/2]

double * Tensor< rank_, dim, double >::end_raw ( )
inherited

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

◆ end_raw() [2/2]

const double * Tensor< rank_, dim, double >::end_raw ( ) const
inherited

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

◆ operator==()

constexpr bool Tensor< rank_, dim, double >::operator== ( const Tensor< rank_, dim, OtherNumber > &  ) const
constexprinherited

Test for equality of two tensors.

◆ operator!=()

constexpr bool Tensor< rank_, dim, double >::operator!= ( const Tensor< rank_, dim, OtherNumber > &  ) const
constexprinherited

Test for inequality of two tensors.

◆ operator+=()

constexpr Tensor & Tensor< rank_, dim, double >::operator+= ( const Tensor< rank_, dim, OtherNumber > &  )
constexprinherited

Add another tensor.

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

◆ operator-=()

constexpr Tensor & Tensor< rank_, dim, double >::operator-= ( const Tensor< rank_, dim, OtherNumber > &  )
constexprinherited

Subtract another tensor.

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

◆ operator*=()

constexpr Tensor & Tensor< rank_, dim, double >::operator*= ( const OtherNumber &  factor)
constexprinherited

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

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

◆ operator/=()

constexpr Tensor & Tensor< rank_, dim, double >::operator/= ( const OtherNumber &  factor)
constexprinherited

Scale the vector by 1/factor.

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

◆ clear()

constexpr void Tensor< rank_, dim, double >::clear ( )
constexprinherited

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

numbers::NumberTraits< double >::real_type Tensor< rank_, dim, double >::norm ( ) const
inherited

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

◆ norm_square()

constexpr numbers::NumberTraits< double >::real_type Tensor< rank_, dim, double >::norm_square ( ) const
constexprinherited

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.

◆ unroll()

void Tensor< rank_, dim, double >::unroll ( Vector< OtherNumber > &  result) const
inherited

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.

◆ component_to_unrolled_index()

static constexpr unsigned int Tensor< rank_, dim, double >::component_to_unrolled_index ( const TableIndices< rank_ > &  indices)
staticconstexprinherited

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

static constexpr TableIndices< rank_ > Tensor< rank_, dim, double >::unrolled_to_component_indices ( const unsigned int  i)
staticconstexprinherited

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

static constexpr std::size_t Tensor< rank_, dim, double >::memory_consumption ( )
staticconstexprinherited

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

◆ unroll_recursion()

void Tensor< rank_, dim, double >::unroll_recursion ( Vector< OtherNumber > &  result,
unsigned int start_index 
) const
privateinherited

Internal helper function for unroll.

Friends And Related Function Documentation

◆ operator*() [1/7]

template<int dim, typename Number = double>
template<typename OtherNumber >
Point< dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > operator* ( const  OtherNumber) const
related

Multiply the current point by a factor.

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

◆ operator<<() [1/3]

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

Output operator for points. Print the elements consecutively, with a space in between.

Definition at line 669 of file point.h.

◆ operator>>()

template<int dim, typename Number >
std::istream & operator>> ( std::istream &  in,
Point< dim, Number > &  p 
)
related

Input operator for points. Inputs the elements consecutively.

Definition at line 687 of file point.h.

◆ operator-() [1/2]

constexpr Tensor< 0, dim, typename ProductType< double , OtherNumber >::type > operator- ( const Tensor< 0, dim, double > &  p,
const Tensor< 0, dim, OtherNumber > &  q 
)
related

Subtract two tensors of rank 0.

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

Definition at line 1966 of file tensor.h.

◆ operator-() [2/2]

constexpr Tensor< rank, dim, typename ProductType< double , OtherNumber >::type > operator- ( const Tensor< rank, dim, double > &  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 CUDA device code.

Definition at line 2132 of file tensor.h.

◆ sum()

Tensor< rank, dim, double > sum ( const Tensor< rank, dim, double > &  local,
const MPI_Comm mpi_communicator 
)
related

Perform an MPI sum of the entries of a tensor.

◆ operator<<() [2/3]

std::ostream & operator<< ( std::ostream &  out,
const Tensor< rank_, dim, double > &  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 1823 of file tensor.h.

◆ operator<<() [3/3]

std::ostream & operator<< ( std::ostream &  out,
const Tensor< 0, dim, double > &  p 
)
related

Output operator for tensors of rank 0. Since such tensors are scalars, we simply print this one value.

Definition at line 1844 of file tensor.h.

◆ operator*() [2/7]

constexpr ProductType< Other, double >::type operator* ( const Other &  object,
const Tensor< 0, dim, double > &  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 CUDA device code.

Definition at line 1872 of file tensor.h.

◆ operator*() [3/7]

constexpr ProductType< double , Other >::type operator* ( const Tensor< 0, dim, double > &  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 CUDA device code.

Definition at line 1892 of file tensor.h.

◆ operator*() [4/7]

constexpr ProductType< double , OtherNumber >::type operator* ( const Tensor< 0, dim, double > &  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 CUDA device code.

Definition at line 1912 of file tensor.h.

◆ operator*() [5/7]

constexpr Tensor< rank, dim, typename ProductType< double , typename EnableIfScalar< OtherNumber >::type >::type > operator* ( const Tensor< rank, dim, double > &  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 CUDA device code.

Definition at line 1991 of file tensor.h.

◆ operator*() [6/7]

constexpr Tensor< rank, dim, typename ProductType< typename EnableIfScalar< double >::type, OtherNumber >::type > operator* ( const double &  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 CUDA device code.

Definition at line 2019 of file tensor.h.

◆ operator*() [7/7]

OtherNumber::type::tensor_type operator* ( const Tensor< rank_1, dim, double > &  src1,
const Tensor< rank_2, dim, OtherNumber > &  src2 
)
related

Definition at line 2232 of file tensor.h.

◆ operator/() [1/2]

constexpr Tensor< 0, dim, typename ProductType< double , typename EnableIfScalar< OtherNumber >::type >::type > operator/ ( const Tensor< 0, dim, double > &  t,
const OtherNumber &  factor 
)
related

Division of a tensor of rank 0 by a scalar number.

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

Definition at line 1933 of file tensor.h.

◆ operator/() [2/2]

constexpr Tensor< rank, dim, typename ProductType< double , typename EnableIfScalar< OtherNumber >::type >::type > operator/ ( const Tensor< rank, dim, double > &  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 CUDA device code.

Definition at line 2090 of file tensor.h.

◆ operator+() [1/2]

constexpr Tensor< 0, dim, typename ProductType< double , OtherNumber >::type > operator+ ( const Tensor< 0, dim, double > &  p,
const Tensor< 0, dim, OtherNumber > &  q 
)
related

Add two tensors of rank 0.

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

Definition at line 1949 of file tensor.h.

◆ operator+() [2/2]

constexpr Tensor< rank, dim, typename ProductType< double , OtherNumber >::type > operator+ ( const Tensor< rank, dim, double > &  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 CUDA device code.

Definition at line 2108 of file tensor.h.

◆ schur_product() [1/2]

constexpr Tensor< 0, dim, typename ProductType< double , OtherNumber >::type > schur_product ( const Tensor< 0, dim, double > &  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 2152 of file tensor.h.

◆ schur_product() [2/2]

constexpr Tensor< rank, dim, typename ProductType< double , OtherNumber >::type > schur_product ( const Tensor< rank, dim, double > &  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 2181 of file tensor.h.

◆ contract()

constexpr Tensor< rank_1+rank_2-2, dim, typenameProductType< double , OtherNumber >::type >::tensor_type contract ( const Tensor< rank_1, dim, double > &  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, this function should be invoked as

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

Definition at line 2288 of file tensor.h.

◆ double_contract()

constexpr Tensor< rank_1+rank_2-4, dim, typenameProductType< double , OtherNumber >::type >::tensor_type double_contract ( const Tensor< rank_1, dim, double > &  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) of a tensor t2, this function should be invoked as

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

Definition at line 2363 of file tensor.h.

◆ scalar_product()

constexpr ProductType< double , OtherNumber >::type scalar_product ( const Tensor< rank, dim, double > &  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} \]

Definition at line 2442 of file tensor.h.

◆ contract3()

constexpr ProductType< T1, typenameProductType< 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.

Definition at line 2479 of file tensor.h.

◆ outer_product()

constexpr Tensor< rank_1+rank_2, dim, typename ProductType< double , OtherNumber >::type > outer_product ( const Tensor< rank_1, dim, double > &  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}.} \]

Definition at line 2508 of file tensor.h.

◆ cross_product_2d()

constexpr Tensor< 1, dim, double > cross_product_2d ( const Tensor< 1, dim, double > &  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).

Definition at line 2539 of file tensor.h.

◆ cross_product_3d()

constexpr Tensor< 1, dim, typename ProductType< Number1, Number2 >::type > cross_product_3d ( const Tensor< 1, dim, Number1 > &  src1,
const Tensor< 1, dim, Number2 > &  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).

Definition at line 2564 of file tensor.h.

◆ determinant() [1/4]

constexpr double determinant ( const Tensor< 2, dim, double > &  t)
related

Compute the determinant of a tensor or rank 2.

Definition at line 2598 of file tensor.h.

◆ determinant() [2/4]

constexpr double determinant ( const Tensor< 2, 1, double > &  t)
related

Specialization for dim==1.

Definition at line 2626 of file tensor.h.

◆ determinant() [3/4]

constexpr double determinant ( const Tensor< 2, 2, double > &  t)
related

Specialization for dim==2.

Definition at line 2638 of file tensor.h.

◆ determinant() [4/4]

constexpr double determinant ( const Tensor< 2, 3, double > &  t)
related

Specialization for dim==3.

Definition at line 2651 of file tensor.h.

◆ trace()

constexpr double trace ( const Tensor< 2, dim, double > &  d)
related

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

Definition at line 2672 of file tensor.h.

◆ invert()

constexpr Tensor< 2, dim, double > invert ( const Tensor< 2, dim, double > &  )
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.

Definition at line 2691 of file tensor.h.

◆ transpose()

constexpr Tensor< 2, dim, double > transpose ( const Tensor< 2, dim, double > &  t)
related

Return the transpose of the given tensor.

Definition at line 2778 of file tensor.h.

◆ adjugate()

constexpr Tensor< 2, dim, double > adjugate ( const Tensor< 2, dim, double > &  t)
related

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

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

Note
This requires that the tensor is invertible.

Definition at line 2809 of file tensor.h.

◆ cofactor()

constexpr Tensor< 2, dim, double > cofactor ( const Tensor< 2, dim, double > &  t)
related

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

\[ \textrm{cof}\mathbf A \dealcoloneq \textrm{det}\mathbf A \; \mathbf{A}^{-T} = \left[ \textrm{adj}\mathbf A \right]^{T} \; . \]

Note
This requires that the tensor is invertible.

Definition at line 2830 of file tensor.h.

◆ project_onto_orthogonal_tensors()

Tensor< 2, dim, double > project_onto_orthogonal_tensors ( const Tensor< 2, dim, double > &  A)
related

Return the nearest orthogonal matrix \(\hat {\mathbf A}=\mathbf U \mathbf{V}^T\) by combining the products of the singular value decomposition (SVD) \({\mathbf A}=\mathbf U \mathbf S \mathbf V^T\) for a given input \({\mathbf A}\), effectively replacing \(\mathbf S\) with the identity matrix.

This is a (nonlinear) projection operation since when applied twice, we have \(\hat{\hat{\mathbf A}}=\hat{\mathbf A}\) as is easy to see. (That is because the SVD of \(\hat {\mathbf A}\) is simply \(\mathbf U \mathbf I \mathbf{V}^T\).) Furthermore, \(\hat {\mathbf A}\) is really an orthogonal matrix because orthogonal matrices have to satisfy \({\hat {\mathbf A}}^T \hat {\mathbf A}={\mathbf I}\), which here implies that

\begin{align*} {\hat {\mathbf A}}^T \hat {\mathbf A} &= \left(\mathbf U \mathbf{V}^T\right)^T\left(\mathbf U \mathbf{V}^T\right) \\ &= \mathbf V \mathbf{U}^T \mathbf U \mathbf{V}^T \\ &= \mathbf V \left(\mathbf{U}^T \mathbf U\right) \mathbf{V}^T \\ &= \mathbf V \mathbf I \mathbf{V}^T \\ &= \mathbf V \mathbf{V}^T \\ &= \mathbf I \end{align*}

due to the fact that the \(\mathbf U\) and \(\mathbf V\) factors that come out of the SVD are themselves orthogonal matrices.

Parameters
AThe tensor for which to find the closest orthogonal tensor.
Template Parameters
NumberThe type used to store the entries of the tensor. Must be either float or double.
Precondition
In order to use this function, this program must be linked with the LAPACK library.
A must not be singular. This is because, conceptually, the problem to be solved here is trying to find a matrix \(\hat{\mathbf A}\) that minimizes some kind of distance from \(\mathbf A\) while satisfying the quadratic constraint \({\hat {\mathbf A}}^T \hat {\mathbf A}={\mathbf I}\). This is not so dissimilar to the kind of problem where one wants to find a vector \(\hat{\mathbf x}\in{\mathbb R}^n\) that minimizes the quadratic objective function \(\|\hat {\mathbf x} - \mathbf x\|^2\) for a given \(\mathbf x\) subject to the constraint \(\|\mathbf x\|^2=1\) – in other words, we are seeking the point \(\hat{\mathbf x}\) on the unit sphere that is closest to \(\mathbf x\). This problem has a solution for all \(\mathbf x\) except if \(\mathbf x=0\). The corresponding condition for the problem we are considering here is that \(\mathbf A\) must not have a zero eigenvalue.

◆ l1_norm()

double l1_norm ( const Tensor< 2, dim, double > &  t)
related

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

◆ linfty_norm()

double linfty_norm ( const Tensor< 2, dim, double > &  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 2939 of file tensor.h.

Member Data Documentation

◆ dimension

constexpr unsigned int Tensor< rank_, dim, double >::dimension
staticconstexprinherited

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

◆ rank

constexpr unsigned int Tensor< rank_, dim, double >::rank
staticconstexprinherited

Publish the rank of this tensor to the outside world.

Definition at line 491 of file tensor.h.

◆ n_independent_components

constexpr unsigned int Tensor< rank_, dim, double >::n_independent_components
staticconstexprinherited

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

◆ values

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

Array of tensors holding the subelements.

Definition at line 814 of file tensor.h.


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