24#include <boost/geometry/core/cs.hpp>
25#include <boost/geometry/geometries/point.hpp>
108template <
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 std::enable_if_t<(dim == dummy_dim) && (dummy_dim != 0),
int> = 0>
169 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>
339 template <
class Archive>
350template <
int dim,
typename Number>
357template <
int dim,
typename Number>
360 const
Tensor<1, dim, Number> &t)
361 :
Tensor<1, dim, Number>(t)
366template <
int dim,
typename Number>
372 "You can only initialize Point<1> objects using the constructor "
373 "that takes only one argument. Point<dim> objects with dim!=1 "
374 "require initialization with the constructor that takes 'dim' "
392template <
int dim,
typename Number>
399 "You can only initialize Point<2> objects using the constructor "
400 "that takes two arguments. Point<dim> objects with dim!=2 "
401 "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>
422 "You can only initialize Point<3> objects using the constructor "
423 "that takes three arguments. Point<dim> objects with dim!=3 "
424 "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>
441template <
std::
size_t dummy_dim,
442 std::enable_if_t<(dim == dummy_dim) && (dummy_dim != 0),
int>>
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>
464 Point<dim, Number>
Point<dim, Number>::unit_vector(
unsigned int i)
472template <
int dim,
typename Number>
475 const
unsigned int index)
const
478 return this->values[
index];
483template <
int dim,
typename Number>
486 const
unsigned int index)
489 return this->values[
index];
494template <
int dim,
typename Number>
496template <typename OtherNumber>
498 &
Point<dim, Number>::operator=(const
Tensor<1, dim, OtherNumber> &p)
506template <
int dim,
typename Number>
509 const
Tensor<1, dim, Number> &p)
const
518template <
int dim,
typename Number>
521 Tensor<1, dim, Number>
Point<dim, Number>::operator-(
522 const
Point<dim, Number> &p)
const
529template <
int dim,
typename Number>
532 const
Tensor<1, dim, Number> &p)
const
541template <
int dim,
typename Number>
547 for (
unsigned int i = 0; i < dim; ++i)
548 result.
values[i] = -this->values[i];
554template <
int dim,
typename Number>
556template <typename OtherNumber>
560 type>
Point<dim, Number>::operator*(const OtherNumber factor)
const
563 for (
unsigned int i = 0; i < dim; ++i)
564 tmp[i] = this->
operator[](i) * factor;
570template <
int dim,
typename Number>
572template <typename OtherNumber>
576 type>
Point<dim, Number>::operator/(const OtherNumber factor)
const
583 ::operator/(base_object, factor));
588template <
int dim,
typename Number>
591 const
Tensor<1, dim, Number> &p)
const
593 Number res = Number();
594 for (
unsigned int i = 0; i < dim; ++i)
595 res += this->
operator[](i) * p[i];
600template <
int dim,
typename Number>
603 Point<dim, Number>::square()
const
605 return this->norm_square();
610template <
int dim,
typename Number>
613 Point<dim, Number>::distance(const
Point<dim, Number> &p)
const
620template <
int dim,
typename Number>
623 Point<dim, Number>::distance_square(const
Point<dim, Number> &p)
const
626 for (
unsigned int i = 0; i < dim; ++i)
628 const Number diff =
static_cast<Number
>(this->values[i]) - p[i];
637template <
int dim,
typename Number>
639template <class Archive>
640inline
void Point<dim, Number>::serialize(Archive &ar, const
unsigned int)
660template <
int dim,
typename Number,
typename OtherNumber>
679template <
int dim,
typename Number>
682operator<<(
std::ostream &out, const
Point<dim, Number> &p)
684 for (
unsigned int i = 0; i < dim - 1; ++i)
699template <
int dim,
typename Number>
704 for (
unsigned int i = 0; i < dim; ++i)
718template <
typename Number>
std::ostream & operator<<(std::ostream &out, const Point< dim, Number > &p)
constexpr Point< dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > operator*(const OtherNumber) const
constexpr Point< dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > operator/(const OtherNumber) const
constexpr Tensor< 1, dim, Number > operator-(const Point< dim, Number > &) const
constexpr Point(const boost::geometry::model::point< Number, dummy_dim, boost::geometry::cs::cartesian > &boost_pt)
constexpr Point< dim, Number > & operator=(const Tensor< 1, dim, OtherNumber > &p)
constexpr numbers::NumberTraits< Number >::real_type distance_square(const Point< dim, Number > &p) const
constexpr Number & operator()(const unsigned int index)
constexpr Point(const Number x, const Number y)
constexpr Point< dim, Number > operator+(const Tensor< 1, dim, Number > &) const
constexpr Point< dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > operator*(const OtherNumber factor, const Point< dim, Number > &p)
constexpr Point(const Number x, const Number y, const Number z)
static constexpr Point< dim, Number > unit_vector(const unsigned int i)
constexpr Point< dim, Number > operator-() const
constexpr Point(const Number x)
void serialize(Archive &ar, const unsigned int version)
constexpr Point(const Tensor< 1, dim, Number > &)
numbers::NumberTraits< Number >::real_type distance(const Point< dim, Number > &p) const
constexpr numbers::NumberTraits< Number >::real_type square() const
constexpr Number operator()(const unsigned int index) const
constexpr Point< dim, Number > operator-(const Tensor< 1, dim, Number > &) const
constexpr Number operator*(const Tensor< 1, dim, Number > &p) const
std::conditional_t< rank_==1, std::array< Number, dim >, std::array< Tensor< rank_ - 1, dim, Number >, dim > > values
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)