26#include <boost/geometry/core/cs.hpp>
27#include <boost/geometry/geometries/point.hpp>
109template <
int dim,
typename Number =
double>
149 Point(
const Number x,
const Number y);
161 Point(
const Number x,
const Number y,
const Number z);
166 template <std::size_t dummy_dim,
167 typename std::enable_if<(dim == dummy_dim) && (dummy_dim != 0),
169 Point(
const boost::geometry::model::
170 point<Number, dummy_dim, boost::geometry::cs::cartesian> &boost_pt);
203 template <
typename OtherNumber>
267 template <
typename OtherNumber>
279 template <
typename OtherNumber>
337 template <
class Archive>
348template <
int dim,
typename Number>
355template <
int dim,
typename Number>
358 :
Tensor<1, dim, Number>(t)
363template <
int dim,
typename Number>
367# ifndef __CUDA_ARCH__
370 "You can only initialize Point<1> objects using the constructor "
371 "that takes only one argument. Point<dim> objects with dim!=1 "
372 "require initialization with the constructor that takes 'dim' "
391template <
int dim,
typename Number>
395# ifndef __CUDA_ARCH__
398 "You can only initialize Point<2> objects using the constructor "
399 "that takes two arguments. Point<dim> objects with dim!=2 "
400 "require initialization with the constructor that takes 'dim' "
407 constexpr unsigned int y_index = (dim < 2) ? 0 : 1;
409 this->
values[y_index] = y;
414template <
int dim,
typename Number>
418# ifndef __CUDA_ARCH__
421 "You can only initialize Point<3> objects using the constructor "
422 "that takes three arguments. Point<dim> objects with dim!=3 "
423 "require initialization with the constructor that takes 'dim' "
430 constexpr unsigned int y_index = (dim < 2) ? 0 : 1;
431 constexpr unsigned int z_index = (dim < 3) ? 0 : 2;
433 this->
values[y_index] = y;
434 this->
values[z_index] = z;
439template <
int dim,
typename Number>
441 std::size_t dummy_dim,
442 typename std::enable_if<(dim == dummy_dim) && (dummy_dim != 0),
int>::type>
444 const boost::geometry::model::
445 point<Number, dummy_dim, boost::geometry::cs::cartesian> &boost_pt)
448 this->
values[0] = boost::geometry::get<0>(boost_pt);
449 constexpr unsigned int y_index = (dim < 2) ? 0 : 1;
450 constexpr unsigned int z_index = (dim < 3) ? 0 : 2;
453 this->
values[y_index] = boost::geometry::get<y_index>(boost_pt);
456 this->
values[z_index] = boost::geometry::get<z_index>(boost_pt);
461template <
int dim,
typename Number>
471template <
int dim,
typename Number>
475# ifndef __CUDA_ARCH__
483template <
int dim,
typename Number>
487# ifndef __CUDA_ARCH__
495template <
int dim,
typename Number>
496template <
typename OtherNumber>
506template <
int dim,
typename Number>
517template <
int dim,
typename Number>
526template <
int dim,
typename Number>
537template <
int dim,
typename Number>
542 for (
unsigned int i = 0; i < dim; ++i)
543 result.
values[i] = -this->values[i];
549template <
int dim,
typename Number>
550template <
typename OtherNumber>
558 for (
unsigned int i = 0; i < dim; ++i)
559 tmp[i] = this->
operator[](i) * factor;
565template <
int dim,
typename Number>
566template <
typename OtherNumber>
578 ::operator/(base_object, factor));
583template <
int dim,
typename Number>
587 Number res = Number();
588 for (
unsigned int i = 0; i < dim; ++i)
589 res += this->
operator[](i) * p[i];
594template <
int dim,
typename Number>
598 return this->norm_square();
603template <
int dim,
typename Number>
612template <
int dim,
typename Number>
617 for (
unsigned int i = 0; i < dim; ++i)
619 const Number diff =
static_cast<Number
>(this->
values[i]) - p(i);
628template <
int dim,
typename Number>
629template <
class Archive>
652template <
int dim,
typename Number,
typename OtherNumber>
669template <
int dim,
typename Number>
673 for (
unsigned int i = 0; i < dim - 1; ++i)
686template <
int dim,
typename Number>
690 for (
unsigned int i = 0; i < dim; ++i)
704template <
typename Number>
std::ostream & operator<<(std::ostream &out, const Point< dim, Number > &p)
Tensor< 1, dim, Number > operator-(const Point< dim, Number > &) const
Number operator*(const Tensor< 1, dim, Number > &p) const
Point(const boost::geometry::model::point< Number, dummy_dim, boost::geometry::cs::cartesian > &boost_pt)
Point(const Tensor< 1, dim, Number > &)
Point(const Number x, const Number y, const Number z)
Point< dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > operator*(const OtherNumber) const
Point< dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > operator/(const OtherNumber) const
Point< dim, Number > operator+(const Tensor< 1, dim, Number > &) const
Number & operator()(const unsigned int index)
void serialize(Archive &ar, const unsigned int version)
Point< dim, Number > operator-() const
numbers::NumberTraits< Number >::real_type distance(const Point< dim, Number > &p) const
static Point< dim, Number > unit_vector(const unsigned int i)
numbers::NumberTraits< Number >::real_type square() const
std::istream & operator>>(std::istream &in, Point< dim, Number > &p)
Point< dim, Number > & operator=(const Tensor< 1, dim, OtherNumber > &p)
numbers::NumberTraits< Number >::real_type distance_square(const Point< dim, Number > &p) const
Point< dim, Number > operator-(const Tensor< 1, dim, Number > &) const
Point(const Number x, const Number y)
Number operator()(const unsigned int index) const
Tensor< rank_ - 1, dim, Number > values[(dim !=0) ? dim :1]
constexpr Tensor & operator=(const Tensor< rank_, dim, OtherNumber > &rhs)
#define DEAL_II_ALWAYS_INLINE
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcMessage(std::string arg1)
T sum(const T &t, const MPI_Comm &mpi_communicator)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
#define DEAL_II_CUDA_HOST_DEV
Point< dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > operator*(const OtherNumber factor, const Point< dim, Number > &p)
static constexpr const T & value(const T &t)
static constexpr std::enable_if< std::is_same< Dummy, number >::value &&is_cuda_compatible< Dummy >::value, real_type >::type abs_square(const number &x)