16 #ifndef dealii__point_h 17 #define dealii__point_h 20 #include <deal.II/base/config.h> 21 #include <deal.II/base/exceptions.h> 22 #include <deal.II/base/tensor.h> 25 DEAL_II_NAMESPACE_OPEN
88 template <
int dim,
typename Number =
double>
109 explicit Point (
const Number x);
118 Point (
const Number x,
128 Point (
const Number x,
142 Number
operator () (
const unsigned int index)
const;
195 template <
typename OtherNumber>
202 template <
typename OtherNumber>
238 template <
class Archive>
239 void serialize(Archive &ar,
const unsigned int version);
246 template <
int dim,
typename Number>
253 template <
int dim,
typename Number>
262 template <
int dim,
typename Number>
267 ExcMessage (
"You can only initialize Point<1> objects using the constructor " 268 "that takes only one argument. Point<dim> objects with dim!=1 " 269 "require initialization with the constructor that takes 'dim' " 288 template <
int dim,
typename Number>
294 ExcMessage (
"You can only initialize Point<2> objects using the constructor " 295 "that takes two arguments. Point<dim> objects with dim!=2 " 296 "require initialization with the constructor that takes 'dim' " 315 template <
int dim,
typename Number>
322 ExcMessage (
"You can only initialize Point<3> objects using the constructor " 323 "that takes three arguments. Point<dim> objects with dim!=3 " 324 "require initialization with the constructor that takes 'dim' " 344 template <
int dim,
typename Number>
355 template <
int dim,
typename Number>
361 return this->values[index];
366 template <
int dim,
typename Number>
372 return this->values[index];
377 template <
int dim,
typename Number>
389 template <
int dim,
typename Number>
399 template <
int dim,
typename Number>
411 template <
int dim,
typename Number>
417 for (
unsigned int i=0; i<dim; ++i)
418 result.values[i] = -this->values[i];
424 template <
int dim,
typename Number>
425 template<
typename OtherNumber>
431 for (
unsigned int i=0; i<dim; ++i)
432 tmp[i] = this->
operator[](i) * factor;
438 template <
int dim,
typename Number>
439 template<
typename OtherNumber>
445 for (
unsigned int i=0; i<dim; ++i)
446 tmp[i] = this->
operator[](i) / factor;
452 template <
int dim,
typename Number>
457 Number res = Number();
458 for (
unsigned int i=0; i<dim; ++i)
459 res += this->
operator[](i) * p[i];
464 template <
int dim,
typename Number>
469 return this->norm_square();
474 template <
int dim,
typename Number>
479 Number
sum = Number();
480 for (
unsigned int i=0; i<dim; ++i)
482 const Number diff=this->values[i]-p(i);
486 return std::sqrt(
sum);
491 template <
int dim,
typename Number>
492 template <
class Archive>
514 template <
int dim,
typename Number,
typename OtherNumber>
530 template <
int dim,
typename Number>
535 for (
unsigned int i=0; i<dim-1; ++i)
549 template <
int dim,
typename Number>
554 for (
unsigned int i=0; i<dim; ++i)
568 template <
typename Number>
579 DEAL_II_NAMESPACE_CLOSE
Point< dim, Number > operator-() const
Point< dim, Number > operator+(const Tensor< 1, dim, Number > &) const
#define AssertIndexRange(index, range)
Number operator()(const unsigned int index) const
static ::ExceptionBase & ExcMessage(std::string arg1)
static real_type abs_square(const number &x)
#define Assert(cond, exc)
static Point< dim, Number > unit_vector(const unsigned int i)
void serialize(Archive &ar, const unsigned int version)
numbers::NumberTraits< Number >::real_type distance(const Point< dim, Number > &p) const
SymmetricTensor< rank, dim, Number > sum(const SymmetricTensor< rank, dim, Number > &local, const MPI_Comm &mpi_communicator)
OutputOperator< VectorType > & operator<<(OutputOperator< VectorType > &out, unsigned int step)
numbers::NumberTraits< Number >::real_type square() const
Point< dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > operator*(const OtherNumber) const
std::istream & operator>>(std::istream &in, Point< dim, Number > &p)
Point< dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > operator/(const OtherNumber) const