|
| 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) |
|
Number | operator() (const unsigned int index) const |
|
Number & | operator() (const unsigned int index) |
|
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 |
|
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 |
|
template<class Archive > |
void | serialize (Archive &ar, const unsigned int version) |
|
| Tensor () |
|
| Tensor (const array_type &initializer) |
|
| Tensor (const Tensor< rank_, dim, OtherNumber > &initializer) |
|
| Tensor (const Tensor< 1, dim, Tensor< rank_-1, dim, OtherNumber > > &initializer) |
|
| operator Tensor< 1, dim, Tensor< rank_-1, dim, OtherNumber > > () const |
|
value_type & | operator[] (const unsigned int i) |
|
const value_type & | operator[] (const unsigned int i) const |
|
const Number & | operator[] (const TableIndices< rank_ > &indices) const |
|
Number & | operator[] (const TableIndices< rank_ > &indices) |
|
Tensor & | operator= (const Tensor< rank_, dim, OtherNumber > &rhs) |
|
Tensor< rank_, dim, Number > & | operator= (const Number d) |
|
bool | operator== (const Tensor< rank_, dim, OtherNumber > &) const |
|
bool | operator!= (const Tensor< rank_, dim, OtherNumber > &) const |
|
Tensor< rank_, dim, Number > & | operator+= (const Tensor< rank_, dim, OtherNumber > &) |
|
Tensor< rank_, dim, Number > & | operator-= (const Tensor< rank_, dim, OtherNumber > &) |
|
Tensor< rank_, dim, Number > & | operator*= (const OtherNumber factor) |
|
Tensor< rank_, dim, Number > & | operator/= (const OtherNumber factor) |
|
Tensor< rank_, dim, Number > | operator- () const |
|
void | clear () |
|
numbers::NumberTraits< Number >::real_type | norm () const |
|
numbers::NumberTraits< Number >::real_type | norm_square () const |
|
void | unroll (Vector< OtherNumber > &result) const |
|
void | serialize (Archive &ar, const unsigned int version) |
|
|
(Note that these are not member functions.)
|
template<typename OtherNumber > |
Point< dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > | operator* (const OtherNumber) const |
|
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, Number > | sum (const Tensor< rank, dim, Number > &local, const MPI_Comm &mpi_communicator) |
|
Tensor< 0, dim, typename ProductType< Number, OtherNumber >::type > | operator- (const Tensor< 0, dim, Number > &p, const Tensor< 0, dim, OtherNumber > &q) |
|
Tensor< rank, dim, typename ProductType< Number, OtherNumber >::type > | operator- (const Tensor< rank, dim, Number > &p, const Tensor< rank, dim, OtherNumber > &q) |
|
ProductType< Other, Number >::type | operator* (const Other object, const Tensor< 0, dim, Number > &t) |
|
ProductType< Number, Other >::type | operator* (const Tensor< 0, dim, Number > &t, const Other object) |
|
ProductType< Number, OtherNumber >::type | operator* (const Tensor< 0, dim, Number > &src1, const Tensor< 0, dim, OtherNumber > &src2) |
|
Tensor< rank, dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > | operator* (const Tensor< rank, dim, Number > &t, const OtherNumber factor) |
|
Tensor< rank, dim, typename ProductType< typename EnableIfScalar< Number >::type, OtherNumber >::type > | operator* (const Number factor, const Tensor< rank, dim, OtherNumber > &t) |
|
Tensor< 0, dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > | operator/ (const Tensor< 0, dim, Number > &t, const OtherNumber factor) |
|
Tensor< rank, dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > | operator/ (const Tensor< rank, dim, Number > &t, const OtherNumber factor) |
|
Tensor< 0, dim, typename ProductType< Number, OtherNumber >::type > | operator+ (const Tensor< 0, dim, Number > &p, const Tensor< 0, dim, OtherNumber > &q) |
|
Tensor< rank, dim, typename ProductType< Number, OtherNumber >::type > | operator+ (const Tensor< rank, dim, Number > &p, const Tensor< rank, dim, OtherNumber > &q) |
|
std::ostream & | operator<< (std::ostream &out, const Tensor< rank_, dim, Number > &p) |
|
std::ostream & | operator<< (std::ostream &out, const Tensor< 0, dim, Number > &p) |
|
DEAL_II_ALWAYS_INLINE Tensor< rank_1+rank_2 - 2, dim, typename ProductType< Number, OtherNumber >::type >::tensor_type | operator* (const Tensor< rank_1, dim, Number > &src1, const Tensor< rank_2, dim, OtherNumber > &src2) |
|
Tensor< rank_1+rank_2 - 2, dim, typename ProductType< Number, OtherNumber >::type >::tensor_type | contract (const Tensor< rank_1, dim, Number > &src1, const Tensor< rank_2, dim, OtherNumber > &src2) |
|
Tensor< rank_1+rank_2 - 4, dim, typename ProductType< Number, OtherNumber >::type >::tensor_type | double_contract (const Tensor< rank_1, dim, Number > &src1, const Tensor< rank_2, dim, OtherNumber > &src2) |
|
ProductType< Number, OtherNumber >::type | scalar_product (const Tensor< rank, dim, Number > &left, const Tensor< rank, dim, OtherNumber > &right) |
|
ProductType< T1, typename ProductType< T2, T3 >::type >::type | contract3 (const TensorT1< rank_1, dim, T1 > &left, const TensorT2< rank_1+rank_2, dim, T2 > &middle, const TensorT3< rank_2, dim, T3 > &right) |
|
Tensor< rank_1+rank_2, dim, typename ProductType< Number, OtherNumber >::type > | outer_product (const Tensor< rank_1, dim, Number > &src1, const Tensor< rank_2, dim, OtherNumber > &src2) |
|
void | contract (Tensor< 2, dim, Number > &dest, const Tensor< 2, dim, Number > &src1, const unsigned int index1, const Tensor< 2, dim, Number > &src2, const unsigned int index3) 1 |
|
void | contract (Tensor< 2, dim, Number > &dest, const Tensor< 3, dim, Number > &src1, const unsigned int index1, const Tensor< 1, dim, Number > &src2) 1 |
|
void | contract (Tensor< 3, dim, Number > &dest, const Tensor< 3, dim, Number > &src1, const unsigned int index1, const Tensor< 2, dim, Number > &src2, const unsigned int index2) 1 |
|
void | contract (Tensor< rank_1+rank_2 - 2, dim, Number > &dest, const Tensor< rank_1, dim, Number > &src1, const Tensor< rank_2, dim, Number > &src2) 1 |
|
ProductType< Number, OtherNumber >::type | contract (const Tensor< 1, dim, Number > &src1, const Tensor< 1, dim, OtherNumber > &src2) 1 |
|
Number | double_contract (const Tensor< 2, dim, Number > &src1, const Tensor< 2, dim, Number > &src2) 1 |
|
void | double_contract (Tensor< 2, dim, Number > &dest, const Tensor< 4, dim, Number > &src1, const Tensor< 2, dim, Number > &src2) 1 |
|
void | outer_product (Tensor< rank_1+rank_2, dim, Number > &dst, const Tensor< rank_1, dim, Number > &src1, const Tensor< rank_2, dim, Number > &src2) 1 |
|
void | outer_product (Tensor< 1, dim, Number > &dst, const Number src1, const Tensor< 1, dim, Number > &src2) 1 |
|
void | outer_product (Tensor< 1, dim, Number > &dst, const Tensor< 1, dim, Number > src1, const Number src2) 1 |
|
Number | determinant (const Tensor< rank, 1, Number > &t) 1 |
|
Number | determinant (const Tensor< 1, 1, Number > &t) 1 |
|
void | cross_product (Tensor< 1, dim, Number > &dst, const Tensor< 1, dim, Number > &src) 1 |
|
void | cross_product (Tensor< 1, dim, Number > &dst, const Tensor< 1, dim, Number > &src1, const Tensor< 1, dim, Number > &src2) 1 |
|
Tensor< 1, dim, Number > | cross_product_2d (const Tensor< 1, dim, Number > &src) |
|
Tensor< 1, dim, Number > | cross_product_3d (const Tensor< 1, dim, Number > &src1, const Tensor< 1, dim, Number > &src2) |
|
Number | determinant (const Tensor< 2, dim, Number > &t) |
|
Number | determinant (const Tensor< 2, 1, Number > &t) |
|
Number | trace (const Tensor< 2, dim, Number > &d) |
|
Tensor< 2, dim, Number > | invert (const Tensor< 2, dim, Number > &) |
|
Tensor< 2, dim, Number > | transpose (const Tensor< 2, dim, Number > &t) |
|
Tensor< 2, dim, Number > | adjugate (const Tensor< 2, dim, Number > &t) |
|
Tensor< 2, dim, Number > | cofactor (const Tensor< 2, dim, Number > &t) |
|
double | l1_norm (const Tensor< 2, dim, Number > &t) |
|
double | linfty_norm (const Tensor< 2, dim, Number > &t) |
|
template<int dim, typename Number = double>
class Point< dim, Number >
A class that represents a point in a space with arbitrary dimension dim
.
Objects of this class are used to represent points, i.e., vectors anchored at the origin of a Cartesian vector space. 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(double x)
and double f(double x, double y)
, you should use double f(Point<dim> &p)
instead as it allows writing dimension independent code.
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_cxx11::array
class. Alternatively, as in the case of vector-valued functions, you can use objects of type Vector or std::vector
.
- Template Parameters
-
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. |
- Author
- Wolfgang Bangerth, 1997
Definition at line 89 of file point.h.