Reference documentation for deal.II version 9.3.3
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
point.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 1998 - 2021 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_point_h
17#define dealii_point_h
18
19
20#include <deal.II/base/config.h>
21
23#include <deal.II/base/tensor.h>
24
26#include <boost/geometry/core/cs.hpp>
27#include <boost/geometry/geometries/point.hpp>
29
30#include <cmath>
31
33
109template <int dim, typename Number = double>
110class Point : public Tensor<1, dim, Number>
111{
112public:
121
125 explicit DEAL_II_CUDA_HOST_DEV
127
136 explicit DEAL_II_CUDA_HOST_DEV
137 Point(const Number x);
138
149 Point(const Number x, const Number y);
150
161 Point(const Number x, const Number y, const Number z);
162
166 template <std::size_t dummy_dim,
167 typename std::enable_if<(dim == dummy_dim) && (dummy_dim != 0),
168 int>::type = 0>
169 Point(const boost::geometry::model::
170 point<Number, dummy_dim, boost::geometry::cs::cartesian> &boost_pt);
171
180 unit_vector(const unsigned int i);
181
188 operator()(const unsigned int index) const;
189
195 DEAL_II_CUDA_HOST_DEV Number &
196 operator()(const unsigned int index);
197
203 template <typename OtherNumber>
206
219
231
242
249 operator-() const;
250
267 template <typename OtherNumber>
269 dim,
270 typename ProductType<Number,
271 typename EnableIfScalar<OtherNumber>::type>::type>
272 operator*(const OtherNumber) const;
273
279 template <typename OtherNumber>
281 dim,
282 typename ProductType<Number,
283 typename EnableIfScalar<OtherNumber>::type>::type>
284 operator/(const OtherNumber) const;
285
292
306 square() const;
307
317
326
336 template <class Archive>
337 void
338 serialize(Archive &ar, const unsigned int version);
339};
340
341/*--------------------------- Inline functions: Point -----------------------*/
342
343#ifndef DOXYGEN
344
345// At least clang-3.7 requires us to have a user-defined constructor
346// and we can't use 'Point<dim,Number>::Point () = default' here.
347template <int dim, typename Number>
350{}
351
352
353
354template <int dim, typename Number>
357 : Tensor<1, dim, Number>(t)
358{}
359
360
361
362template <int dim, typename Number>
364Point<dim, Number>::Point(const Number x)
365{
366# ifndef __CUDA_ARCH__
367 Assert(dim == 1,
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' "
372 "arguments."));
373# endif
374
375 // we can only get here if we pass the assertion. use the switch anyway so
376 // as to avoid compiler warnings about uninitialized elements or writing
377 // beyond the end of the 'values' array
378 switch (dim)
379 {
380 case 1:
381 this->values[0] = x;
382 break;
383
384 default:;
385 }
386}
387
388
389
390template <int dim, typename Number>
392Point<dim, Number>::Point(const Number x, const Number y)
393{
394# ifndef __CUDA_ARCH__
395 Assert(dim == 2,
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' "
400 "arguments."));
401# endif
402
403 // we can only get here if we pass the assertion. use the indirection anyway
404 // so as to avoid compiler warnings about uninitialized elements or writing
405 // beyond the end of the 'values' array
406 constexpr unsigned int y_index = (dim < 2) ? 0 : 1;
407 this->values[0] = x;
408 this->values[y_index] = y;
409}
410
411
412
413template <int dim, typename Number>
415Point<dim, Number>::Point(const Number x, const Number y, const Number z)
416{
417# ifndef __CUDA_ARCH__
418 Assert(dim == 3,
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' "
423 "arguments."));
424# endif
425
426 // we can only get here if we pass the assertion. use the indirection anyway
427 // so as to avoid compiler warnings about uninitialized elements or writing
428 // beyond the end of the 'values' array
429 constexpr unsigned int y_index = (dim < 2) ? 0 : 1;
430 constexpr unsigned int z_index = (dim < 3) ? 0 : 2;
431 this->values[0] = x;
432 this->values[y_index] = y;
433 this->values[z_index] = z;
434}
435
436
437
438template <int dim, typename Number>
439template <
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)
445{
446 Assert(dim <= 3, ExcNotImplemented());
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;
450
451 if (dim >= 2)
452 this->values[y_index] = boost::geometry::get<y_index>(boost_pt);
453
454 if (dim >= 3)
455 this->values[z_index] = boost::geometry::get<z_index>(boost_pt);
456}
457
458
459
460template <int dim, typename Number>
462 Point<dim, Number>::unit_vector(unsigned int i)
463{
465 p[i] = 1.;
466 return p;
467}
468
469
470template <int dim, typename Number>
471inline DEAL_II_CUDA_HOST_DEV Number
472Point<dim, Number>::operator()(const unsigned int index) const
473{
474# ifndef __CUDA_ARCH__
475 AssertIndexRange(static_cast<int>(index), dim);
476# endif
477 return this->values[index];
478}
479
480
481
482template <int dim, typename Number>
483inline DEAL_II_CUDA_HOST_DEV Number &
484Point<dim, Number>::operator()(const unsigned int index)
485{
486# ifndef __CUDA_ARCH__
487 AssertIndexRange(static_cast<int>(index), dim);
488# endif
489 return this->values[index];
490}
491
492
493
494template <int dim, typename Number>
495template <typename OtherNumber>
498{
500 return *this;
501}
502
503
504
505template <int dim, typename Number>
508{
509 Point<dim, Number> tmp = *this;
510 tmp += p;
511 return tmp;
512}
513
514
515
516template <int dim, typename Number>
519{
520 return (Tensor<1, dim, Number>(*this) -= p);
521}
522
523
524
525template <int dim, typename Number>
528{
529 Point<dim, Number> tmp = *this;
530 tmp -= p;
531 return tmp;
532}
533
534
535
536template <int dim, typename Number>
539{
540 Point<dim, Number> result;
541 for (unsigned int i = 0; i < dim; ++i)
542 result.values[i] = -this->values[i];
543 return result;
544}
545
546
547
548template <int dim, typename Number>
549template <typename OtherNumber>
551 Point<dim,
552 typename ProductType<Number,
553 typename EnableIfScalar<OtherNumber>::type>::type>
554 Point<dim, Number>::operator*(const OtherNumber factor) const
555{
557 for (unsigned int i = 0; i < dim; ++i)
558 tmp[i] = this->operator[](i) * factor;
559 return tmp;
560}
561
562
563
564template <int dim, typename Number>
565template <typename OtherNumber>
567 Point<dim,
568 typename ProductType<Number,
569 typename EnableIfScalar<OtherNumber>::type>::type>
570 Point<dim, Number>::operator/(const OtherNumber factor) const
571{
572 const Tensor<1, dim, Number> &base_object = *this;
573 return Point<
574 dim,
575 typename ProductType<Number,
576 typename EnableIfScalar<OtherNumber>::type>::type>(
577 ::operator/(base_object, factor));
578}
579
580
581
582template <int dim, typename Number>
585{
586 Number res = Number();
587 for (unsigned int i = 0; i < dim; ++i)
588 res += this->operator[](i) * p[i];
589 return res;
590}
591
592
593template <int dim, typename Number>
596{
597 return this->norm_square();
598}
599
600
601
602template <int dim, typename Number>
605{
606 return std::sqrt(distance_square(p));
607}
608
609
610
611template <int dim, typename Number>
614{
616 for (unsigned int i = 0; i < dim; ++i)
617 {
618 const Number diff = static_cast<Number>(this->values[i]) - p(i);
620 }
621
622 return sum;
623}
624
625
626
627template <int dim, typename Number>
628template <class Archive>
629inline void
630Point<dim, Number>::serialize(Archive &ar, const unsigned int)
631{
632 // forward to serialization
633 // function in the base class
634 ar &static_cast<Tensor<1, dim, Number> &>(*this);
635}
636
637#endif // DOXYGEN
638
639
640/*--------------------------- Global functions: Point -----------------------*/
641
642
651template <int dim, typename Number, typename OtherNumber>
653 Point<dim,
654 typename ProductType<Number,
655 typename EnableIfScalar<OtherNumber>::type>::type>
656 operator*(const OtherNumber factor, const Point<dim, Number> &p)
657{
658 return p * factor;
659}
660
661
662
668template <int dim, typename Number>
669inline std::ostream &
670operator<<(std::ostream &out, const Point<dim, Number> &p)
671{
672 for (unsigned int i = 0; i < dim - 1; ++i)
673 out << p[i] << ' ';
674 out << p[dim - 1];
675
676 return out;
677}
678
679
680
685template <int dim, typename Number>
686inline std::istream &
687operator>>(std::istream &in, Point<dim, Number> &p)
688{
689 for (unsigned int i = 0; i < dim; ++i)
690 in >> p[i];
691
692 return in;
693}
694
695
696#ifndef DOXYGEN
697
703template <typename Number>
704inline std::ostream &
705operator<<(std::ostream &out, const Point<1, Number> &p)
706{
707 out << p[0];
708
709 return out;
710}
711
712#endif // DOXYGEN
714
715#endif
OutputOperator< VectorType > & operator<<(OutputOperator< VectorType > &out, unsigned int step)
Definition: operator.h:165
Definition: point.h:111
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
Point(const Number x)
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)
Definition: point.h:687
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
Definition: tensor.h:472
constexpr numbers::NumberTraits< double >::real_type norm_square() const
Tensor< rank_ - 1, dim, double > values[(dim !=0) ? dim :1]
Definition: tensor.h:814
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
Definition: config.h:100
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
Definition: config.h:416
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
Definition: config.h:454
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
Definition: exceptions.h:1465
#define AssertIndexRange(index, range)
Definition: exceptions.h:1690
static ::ExceptionBase & ExcMessage(std::string arg1)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
#define DEAL_II_CUDA_HOST_DEV
Definition: numbers.h:34
Point< dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > operator*(const OtherNumber factor, const Point< dim, Number > &p)
Definition: point.h:656
static constexpr const T & value(const T &t)
Definition: numbers.h:693
static constexpr std::enable_if< std::is_same< Dummy, number >::value &&is_cuda_compatible< Dummy >::value, real_type >::type abs_square(const number &x)
Definition: numbers.h:577