Reference documentation for deal.II version 9.6.0
|
#include <deal.II/base/point.h>
Public Types | |
using | value_type |
using | array_type |
using | tensor_type |
Public Member Functions | |
constexpr | Point () |
constexpr | Point (const Tensor< 1, dim, Number > &) |
constexpr | Point (const Number x) |
constexpr | Point (const Number x, const Number y) |
constexpr | Point (const Number x, const Number y, const Number z) |
template<std::size_t dummy_dim, std::enable_if_t<(dim==dummy_dim) &&(dummy_dim !=0), int > = 0> | |
constexpr | Point (const boost::geometry::model::point< Number, dummy_dim, boost::geometry::cs::cartesian > &boost_pt) |
constexpr Number | operator() (const unsigned int index) const |
constexpr Number & | operator() (const unsigned int index) |
template<typename OtherNumber > | |
constexpr 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_type & | operator[] (const unsigned int i) |
constexpr const value_type & | operator[] (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 Tensor & | operator+= (const Tensor< rank_, dim, OtherNumber > &) |
constexpr Tensor & | operator-= (const Tensor< rank_, dim, OtherNumber > &) |
constexpr Tensor & | operator*= (const OtherNumber &factor) |
constexpr Tensor & | operator/= (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 |
void | unroll (const Iterator begin, const Iterator end) const |
Addition and subtraction of points. | |
constexpr Point< dim, Number > | operator+ (const Tensor< 1, dim, Number > &) const |
constexpr Tensor< 1, dim, Number > | operator- (const Point< dim, Number > &) const |
constexpr Point< dim, Number > | operator- (const Tensor< 1, dim, Number > &) const |
constexpr Point< dim, Number > | operator- () const |
Static Public Member Functions | |
static constexpr 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 Attributes | |
std::conditional_t< rank_==1, std::array< double, dim >, std::array< Tensor< rank_ - 1, dim, double >, dim > > | values |
Related Symbols | |
(Note that these are not member symbols.) | |
template<int dim, typename Number , typename OtherNumber > | |
constexpr Point< dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > | operator* (const OtherNumber factor, const Point< dim, Number > &p) |
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 DEAL_II_HOST_DEVICE_ALWAYS_INLINE 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 DEAL_II_HOST_DEVICE_ALWAYS_INLINE 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 | |
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 > | |
constexpr Point< dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > | operator/ (const OtherNumber) const |
constexpr Number | operator* (const Tensor< 1, dim, Number > &p) const |
constexpr numbers::NumberTraits< Number >::real_type | square () const |
numbers::NumberTraits< Number >::real_type | distance (const Point< dim, Number > &p) const |
constexpr numbers::NumberTraits< Number >::real_type | distance_square (const Point< dim, Number > &p) const |
template<typename OtherNumber > | |
constexpr Point< dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > | operator* (const OtherNumber) const |
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).
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
.
dim | An 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. |
Number | The 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. |
|
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 rank-1 tensors.
|
inherited |
Declare an array type which can be used to initialize an object of this type statically. For rank-1 tensors, this array is simply an array of length dim
of scalars of type Number
. For higher-rank tensors, it is an array of length dim
of the array_type
of the next lower-rank tensor.
This class is occasionally instantiated for dim == 0
. C++ does not allow the creation of zero-length arrays. As a consequence, if dim==0
, then all arrays that should have length dim
are instead declared as having length 1 to avoid compiler errors.
|
inherited |
Internal type declaration that is used to specialize the return type of operator[]() for Tensor<1,dim,Number>
Standard constructor. Creates an object that corresponds to the origin, i.e., all coordinates are set to zero.
|
explicitconstexpr |
Convert a tensor to a point.
|
explicitconstexpr |
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.
|
constexpr |
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).
|
constexpr |
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).
|
constexpr |
Convert a boost::geometry::point to a Point.
|
staticconstexpr |
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 i
th coordinate.
|
constexpr |
Read access to the index
th coordinate.
|
constexpr |
Read and write access to the index
th coordinate.
|
constexpr |
Assignment operator from Tensor<1, dim, Number> with different underlying scalar type. This obviously requires that the OtherNumber
type is convertible to Number
.
|
constexpr |
Add an offset given as Tensor<1,dim,Number> to a point.
|
constexpr |
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>.
|
constexpr |
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>.
|
constexpr |
The opposite vector.
|
constexpr |
Divide the current point by a factor.
|
constexpr |
Return the scalar product of the vectors representing two points.
|
constexpr |
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.
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.
|
constexpr |
Return the squared Euclidean distance of this
point to the point p
.
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.
|
constexprinherited |
Conversion operator to tensor of tensors.
|
constexprinherited |
Read-Write access operator.
|
constexprinherited |
Read-only access operator.
|
constexprinherited |
Read access using TableIndices indices
|
constexprinherited |
Read and write access using TableIndices indices
|
inherited |
Return a pointer to the first element of the underlying storage.
|
inherited |
Return a const pointer to the first element of the underlying storage.
|
inherited |
Return a pointer to the element past the end of the underlying storage.
|
inherited |
Return a pointer to the element past the end of the underlying storage.
|
constexprinherited |
Test for equality of two tensors.
|
constexprinherited |
Test for inequality of two tensors.
|
constexprinherited |
Add another tensor.
|
constexprinherited |
Subtract another tensor.
Scale the tensor by factor
, i.e. multiply all components by factor
.
Scale the vector by 1/factor
.
|
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 an empty 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.
|
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.
|
constexprinherited |
Return the square of the Frobenius-norm of a tensor, i.e. the sum of the absolute squares of all entries.
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.
|
inherited |
Fill a range with all tensor elements.
This function unrolls all tensor entries into a single, linearly numbered sequence. The order of the elements is the one given by component_to_unrolled_index().
The template type Number must be convertible to the type of *begin
.
|
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.
|
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.
|
staticconstexprinherited |
Determine an estimate for the memory consumption (in bytes) of this object.
|
related |
Multiply the current point by a factor.
|
related |
|
related |
Output operator for points. Print the elements consecutively, with a space in between.
|
related |
Input operator for points. Inputs the elements consecutively.
|
related |
|
related |
Perform an MPI sum of the entries of a tensor.
|
related |
|
related |
|
related |
|
related |
|
related |
|
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.
|
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.
|
related |
|
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.
|
related |
|
related |
|
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} \]
rank | The rank of both tensors. |
|
related |
The dot product (single contraction) for tensors. This function return a tensor of rank \((\text{rank}_1 + \text{rank}_2 - 2)\) that is the contraction of the last index of a tensor src1
of rank rank_1
with the first index of a tensor src2
of rank rank_2:
\[ \text{result}_{i_1,\ldots,i_{r1},j_1,\ldots,j_{r2}} = \sum_{k} \text{left}_{i_1,\ldots,i_{r1}, k} \text{right}_{k, j_1,\ldots,j_{r2}} \]
operator*()
performs a double contraction. The origin of the difference in how operator*()
is implemented between Tensor and SymmetricTensor is that for the former, the product between two Tensor objects of same rank and dimension results in another Tensor object – that it, operator*()
corresponds to the multiplicative group action within the group of tensors. On the other hand, there is no corresponding multiplicative group action with the set of symmetric tensors because, in general, the product of two symmetric tensors is a nonsymmetric tensor. As a consequence, for a mathematician, it is clear that operator*()
for symmetric tensors must have a different meaning: namely the dot or scalar product that maps two symmetric tensors of rank 2 to a scalar. This corresponds to the double-dot (colon) operator whose meaning is then extended to the product of any two even-ranked symmetric tensors.rank_1==rank_2==1
, then a scalar number is returned as an unwrapped number type. Return the \(l_1\) norm of the given rank-2 tensor, where \(\|\mathbf T\|_1 = \max_j \sum_i |T_{ij}|\) (maximum of the sums over columns).
|
related |
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.
Number of independent components of a tensor of current rank. This is dim times the number of independent components of each sub-tensor, which equates to dim^rank
.
This number can be used to loop over all of the entries of a tensor, using the unrolled_to_component_indices() function: