25#include <boost/geometry/core/cs.hpp>
26#include <boost/geometry/geometries/point.hpp>
109template <
int dim,
typename Number =
double>
150 Point(
const Number x,
const Number y);
162 Point(
const Number x,
const Number y,
const Number z);
167 template <std::size_t dummy_dim,
168 std::enable_if_t<(dim == dummy_dim) && (dummy_dim != 0),
int> = 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 const
Tensor<1, dim, Number> &t)
359 :
Tensor<1, dim, Number>(t)
364template <
int dim,
typename Number>
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' "
390template <
int dim,
typename Number>
397 "You can only initialize Point<2> objects using the constructor "
398 "that takes two arguments. Point<dim> objects with dim!=2 "
399 "require initialization with the constructor that takes 'dim' "
405 constexpr unsigned int y_index = (dim < 2) ? 0 : 1;
407 this->
values[y_index] = y;
412template <
int dim,
typename Number>
420 "You can only initialize Point<3> objects using the constructor "
421 "that takes three arguments. Point<dim> objects with dim!=3 "
422 "require initialization with the constructor that takes 'dim' "
428 constexpr unsigned int y_index = (dim < 2) ? 0 : 1;
429 constexpr unsigned int z_index = (dim < 3) ? 0 : 2;
431 this->
values[y_index] = y;
432 this->
values[z_index] = z;
437template <
int dim,
typename Number>
439template <
std::
size_t dummy_dim,
440 std::enable_if_t<(dim == dummy_dim) && (dummy_dim != 0),
int>>
442 const
boost::geometry::model::
443 point<Number, dummy_dim,
boost::geometry::cs::cartesian> &boost_pt)
446 this->
values[0] = boost::geometry::get<0>(boost_pt);
447 constexpr unsigned int y_index = (dim < 2) ? 0 : 1;
448 constexpr unsigned int z_index = (dim < 3) ? 0 : 2;
451 this->
values[y_index] = boost::geometry::get<y_index>(boost_pt);
454 this->
values[z_index] = boost::geometry::get<z_index>(boost_pt);
459template <
int dim,
typename Number>
470template <
int dim,
typename Number>
473 const
unsigned int index)
const
481template <
int dim,
typename Number>
484 const
unsigned int index)
492template <
int dim,
typename Number>
494template <typename OtherNumber>
496 const
Tensor<1, dim, OtherNumber> &p)
504template <
int dim,
typename Number>
507 const
Tensor<1, dim, Number> &p)
const
516template <
int dim,
typename Number>
519 const
Point<dim, Number> &p)
const
526template <
int dim,
typename Number>
529 const
Tensor<1, dim, Number> &p)
const
538template <
int dim,
typename Number>
544 for (
unsigned int i = 0; i < dim; ++i)
545 result.
values[i] = -this->values[i];
551template <
int dim,
typename Number>
553template <typename OtherNumber>
557 type>
Point<dim, Number>::operator*(const OtherNumber factor)
const
560 for (
unsigned int i = 0; i < dim; ++i)
561 tmp[i] = this->
operator[](i) * factor;
567template <
int dim,
typename Number>
569template <typename OtherNumber>
573 type>
Point<dim, Number>::operator/(const OtherNumber factor)
const
580 ::operator/(base_object, factor));
585template <
int dim,
typename Number>
588 const
Tensor<1, dim, Number> &p)
const
590 Number res = Number();
591 for (
unsigned int i = 0; i < dim; ++i)
592 res += this->
operator[](i) * p[i];
597template <
int dim,
typename Number>
600 Point<dim, Number>::square()
const
602 return this->norm_square();
607template <
int dim,
typename Number>
610 Point<dim, Number>::distance(const
Point<dim, Number> &p)
const
617template <
int dim,
typename Number>
620 Point<dim, Number>::distance_square(const
Point<dim, Number> &p)
const
623 for (
unsigned int i = 0; i < dim; ++i)
625 const Number diff =
static_cast<Number
>(this->
values[i]) - p(i);
634template <
int dim,
typename Number>
636template <class Archive>
637inline
void Point<dim, Number>::serialize(Archive &ar, const
unsigned int)
657template <
int dim,
typename Number,
typename OtherNumber>
676template <
int dim,
typename Number>
679operator<<(
std::ostream &out, const
Point<dim, Number> &p)
681 for (
unsigned int i = 0; i < dim - 1; ++i)
696template <
int dim,
typename Number>
701 for (
unsigned int i = 0; i < dim; ++i)
715template <
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
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
Point< dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > operator*(const OtherNumber factor, const Point< dim, Number > &p)
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_CXX20_REQUIRES(condition)
#define DEAL_II_NAMESPACE_CLOSE
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_HOST_DEVICE
static constexpr DEAL_II_HOST_DEVICE_ALWAYS_INLINE const T & value(const T &t)
static constexpr real_type abs_square(const number &x)