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>
336 template <
class Archive>
347template <
int dim,
typename Number>
354template <
int dim,
typename Number>
357 :
Tensor<1, dim, Number>(t)
362template <
int dim,
typename Number>
366# ifndef __CUDA_ARCH__
369 "You can only initialize Point<1> objects using the constructor "
370 "that takes only one argument. Point<dim> objects with dim!=1 "
371 "require initialization with the constructor that takes 'dim' "
390template <
int dim,
typename Number>
394# ifndef __CUDA_ARCH__
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' "
406 constexpr unsigned int y_index = (dim < 2) ? 0 : 1;
408 this->
values[y_index] = y;
413template <
int dim,
typename Number>
417# ifndef __CUDA_ARCH__
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' "
429 constexpr unsigned int y_index = (dim < 2) ? 0 : 1;
430 constexpr unsigned int z_index = (dim < 3) ? 0 : 2;
432 this->
values[y_index] = y;
433 this->
values[z_index] = z;
438template <
int dim,
typename Number>
440 std::size_t dummy_dim,
441 typename std::enable_if<(dim == dummy_dim) && (dummy_dim != 0),
int>::type>
443 const boost::geometry::model::
444 point<Number, dummy_dim, boost::geometry::cs::cartesian> &boost_pt)
447 this->
values[0] = boost::geometry::get<0>(boost_pt);
448 constexpr unsigned int y_index = (dim < 2) ? 0 : 1;
449 constexpr unsigned int z_index = (dim < 3) ? 0 : 2;
452 this->
values[y_index] = boost::geometry::get<y_index>(boost_pt);
455 this->
values[z_index] = boost::geometry::get<z_index>(boost_pt);
460template <
int dim,
typename Number>
470template <
int dim,
typename Number>
474# ifndef __CUDA_ARCH__
477 return this->
values[index];
482template <
int dim,
typename Number>
486# ifndef __CUDA_ARCH__
489 return this->
values[index];
494template <
int dim,
typename Number>
495template <
typename OtherNumber>
505template <
int dim,
typename Number>
516template <
int dim,
typename Number>
525template <
int dim,
typename Number>
536template <
int dim,
typename Number>
541 for (
unsigned int i = 0; i < dim; ++i)
542 result.
values[i] = -this->values[i];
548template <
int dim,
typename Number>
549template <
typename OtherNumber>
557 for (
unsigned int i = 0; i < dim; ++i)
558 tmp[i] = this->
operator[](i) * factor;
564template <
int dim,
typename Number>
565template <
typename OtherNumber>
577 ::operator/(base_object, factor));
582template <
int dim,
typename Number>
586 Number res = Number();
587 for (
unsigned int i = 0; i < dim; ++i)
588 res += this->
operator[](i) * p[i];
593template <
int dim,
typename Number>
602template <
int dim,
typename Number>
611template <
int dim,
typename Number>
616 for (
unsigned int i = 0; i < dim; ++i)
618 const Number diff =
static_cast<Number
>(this->
values[i]) - p(i);
627template <
int dim,
typename Number>
628template <
class Archive>
651template <
int dim,
typename Number,
typename OtherNumber>
668template <
int dim,
typename Number>
672 for (
unsigned int i = 0; i < dim - 1; ++i)
685template <
int dim,
typename Number>
689 for (
unsigned int i = 0; i < dim; ++i)
703template <
typename Number>
OutputOperator< VectorType > & operator<<(OutputOperator< VectorType > &out, unsigned int step)
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
constexpr numbers::NumberTraits< double >::real_type norm_square() const
Tensor< rank_ - 1, dim, double > values[(dim !=0) ? dim :1]
Tensor< rank, dim, double > sum(const Tensor< rank, dim, double > &local, const MPI_Comm &mpi_communicator)
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)
::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)