16 #ifndef dealii_point_h 17 #define dealii_point_h 20 #include <deal.II/base/config.h> 22 #include <deal.II/base/exceptions.h> 23 #include <deal.II/base/tensor.h> 25 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
26 #include <boost/geometry.hpp> 27 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
31 DEAL_II_NAMESPACE_OPEN
109 template <
int dim,
typename Number =
double>
130 explicit Point(
const Number x);
139 Point(
const Number x,
const Number y);
148 Point(
const Number x,
const Number y,
const Number z);
153 template <std::size_t dummy_dim,
154 typename std::enable_if<(dim == dummy_dim) && (dummy_dim != 0),
156 Point(
const boost::geometry::model::
157 point<Number, dummy_dim, boost::geometry::cs::cartesian> &boost_pt);
229 template <
typename OtherNumber>
238 template <
typename OtherNumber>
285 template <
class Archive>
287 serialize(Archive &ar,
const unsigned int version);
296 template <
int dim,
typename Number>
302 template <
int dim,
typename Number>
304 :
Tensor<1, dim, Number>(t)
309 template <
int dim,
typename Number>
314 "You can only initialize Point<1> objects using the constructor " 315 "that takes only one argument. Point<dim> objects with dim!=1 " 316 "require initialization with the constructor that takes 'dim' " 334 template <
int dim,
typename Number>
339 "You can only initialize Point<2> objects using the constructor " 340 "that takes two arguments. Point<dim> objects with dim!=2 " 341 "require initialization with the constructor that takes 'dim' " 346 constexpr
unsigned int y_index = (dim < 2) ? 0 : 1;
348 this->values[y_index] = y;
353 template <
int dim,
typename Number>
358 "You can only initialize Point<3> objects using the constructor " 359 "that takes three arguments. Point<dim> objects with dim!=3 " 360 "require initialization with the constructor that takes 'dim' " 366 constexpr
unsigned int y_index = (dim < 2) ? 0 : 1;
367 constexpr
unsigned int z_index = (dim < 3) ? 0 : 2;
369 this->values[y_index] = y;
370 this->values[z_index] = z;
375 template <
int dim,
typename Number>
377 std::size_t dummy_dim,
378 typename std::enable_if<(dim == dummy_dim) && (dummy_dim != 0), int>::type>
380 const boost::geometry::model::
381 point<Number, dummy_dim, boost::geometry::cs::cartesian> &boost_pt)
384 this->values[0] = boost::geometry::get<0>(boost_pt);
385 constexpr
unsigned int y_index = (dim < 2) ? 0 : 1;
386 constexpr
unsigned int z_index = (dim < 3) ? 0 : 2;
389 this->values[y_index] = boost::geometry::get<y_index>(boost_pt);
392 this->values[z_index] = boost::geometry::get<z_index>(boost_pt);
397 template <
int dim,
typename Number>
407 template <
int dim,
typename Number>
412 return this->values[index];
417 template <
int dim,
typename Number>
422 return this->values[index];
427 template <
int dim,
typename Number>
438 template <
int dim,
typename Number>
447 template <
int dim,
typename Number>
458 template <
int dim,
typename Number>
463 for (
unsigned int i = 0; i < dim; ++i)
464 result.values[i] = -this->values[i];
470 template <
int dim,
typename Number>
471 template <
typename OtherNumber>
479 for (
unsigned int i = 0; i < dim; ++i)
480 tmp[i] = this->
operator[](i) * factor;
486 template <
int dim,
typename Number>
487 template <
typename OtherNumber>
495 for (
unsigned int i = 0; i < dim; ++i)
496 tmp[i] = this->
operator[](i) / factor;
502 template <
int dim,
typename Number>
506 Number res = Number();
507 for (
unsigned int i = 0; i < dim; ++i)
508 res += this->
operator[](i) * p[i];
513 template <
int dim,
typename Number>
517 return this->norm_square();
522 template <
int dim,
typename Number>
526 return std::sqrt(distance_square(p));
531 template <
int dim,
typename Number>
536 for (
unsigned int i = 0; i < dim; ++i)
538 const Number diff =
static_cast<Number
>(this->values[i]) - p(i);
547 template <
int dim,
typename Number>
548 template <
class Archive>
569 template <
int dim,
typename Number,
typename OtherNumber>
586 template <
int dim,
typename Number>
587 inline std::ostream &
588 operator<<(std::ostream &out, const Point<dim, Number> &p)
590 for (
unsigned int i = 0; i < dim - 1; ++i)
603 template <
int dim,
typename Number>
604 inline std::istream &
607 for (
unsigned int i = 0; i < dim; ++i)
621 template <
typename Number>
622 inline std::ostream &
623 operator<<(std::ostream &out, const Point<1, Number> &p)
631 DEAL_II_NAMESPACE_CLOSE
Point< dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > operator*(const OtherNumber factor, const Point< dim, Number > &p)
Point< dim, Number > operator+(const Tensor< 1, dim, Number > &) const
static constexpr std::enable_if< std::is_same< Dummy, number >::value &&is_cuda_compatible< Dummy >::value, real_type >::type abs_square(const number &x)
#define AssertIndexRange(index, range)
static Point< dim, Number > unit_vector(const unsigned int i)
Number operator()(const unsigned int index) const
static ::ExceptionBase & ExcMessage(std::string arg1)
Point< dim, Number > operator-() const
#define Assert(cond, exc)
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)
Point< dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > operator/(const OtherNumber) const
numbers::NumberTraits< Number >::real_type square() const
static ::ExceptionBase & ExcNotImplemented()
numbers::NumberTraits< Number >::real_type distance_square(const Point< dim, Number > &p) 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