Reference documentation for deal.II version GIT relicensing-437-g81ec864850 2024-04-19 07:30:02+00:00
\(\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\}}\)
Loading...
Searching...
No Matches
bounding_box.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2017 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_base_bounding_box_h
16#define dealii_base_bounding_box_h
17
18
19#include <deal.II/base/config.h>
20
22#include <deal.II/base/point.h>
23
24#include <limits>
25
27
32enum class NeighborType
33{
37 not_neighbors = 0,
38
46
54
68};
69
134template <int spacedim, typename Number = double>
136{
137public:
142 BoundingBox() = default;
143
148
154
159
167
175 template <class Container>
176 BoundingBox(const Container &points);
177
181 std::pair<Point<spacedim, Number>, Point<spacedim, Number>> &
183
187 const std::pair<Point<spacedim, Number>, Point<spacedim, Number>> &
189
193 bool
195
199 bool
201
206 bool
208 const BoundingBox<spacedim, Number> &other_bbox,
209 const double tolerance = std::numeric_limits<Number>::epsilon()) const;
210
216 const BoundingBox<spacedim, Number> &other_bbox,
217 const double tolerance = std::numeric_limits<Number>::epsilon()) const;
218
224 void
225 merge_with(const BoundingBox<spacedim, Number> &other_bbox);
226
233 bool
236 const double tolerance = std::numeric_limits<Number>::epsilon()) const;
237
248 void
249 extend(const Number amount);
250
256 create_extended(const Number amount) const;
257
272 create_extended_relative(const Number relative_amount) const;
273
277 double
278 volume() const;
279
284 center() const;
285
289 Number
290 side_length(const unsigned int direction) const;
291
295 Number
296 lower_bound(const unsigned int direction) const;
297
301 Number
302 upper_bound(const unsigned int direction) const;
303
308 bounds(const unsigned int direction) const;
309
315 vertex(const unsigned int index) const;
316
322 child(const unsigned int index) const;
323
331 BoundingBox<spacedim - 1, Number>
332 cross_section(const unsigned int direction) const;
333
343 real_to_unit(const Point<spacedim, Number> &point) const;
344
354 unit_to_real(const Point<spacedim, Number> &point) const;
362 Number
364 const unsigned int direction) const;
365
371 Number
372 signed_distance(const Point<spacedim, Number> &point) const;
373
379 template <class Archive>
380 void
381 serialize(Archive &ar, const unsigned int version);
382
383private:
384 std::pair<Point<spacedim, Number>, Point<spacedim, Number>> boundary_points;
385};
386
392template <typename Number>
393class BoundingBox<0, Number>
394{
395public:
400
405
409 template <class Container>
410 BoundingBox(const Container &);
411};
412
413
419template <int dim, typename Number = double>
422
423
424namespace internal
425{
458 template <int dim>
459 inline int
460 coordinate_to_one_dim_higher(const int locked_coordinate,
461 const int coordinate_in_dim)
462 {
463 AssertIndexRange(locked_coordinate, dim + 1);
464 AssertIndexRange(coordinate_in_dim, dim);
465 return (locked_coordinate + coordinate_in_dim + 1) % (dim + 1);
466 }
467
468} // namespace internal
469
470/*------------------------ Inline functions: BoundingBox --------------------*/
471
472#ifndef DOXYGEN
473
474
475template <int spacedim, typename Number>
478 : BoundingBox({p, p})
479{}
480
481
482
483template <int spacedim, typename Number>
486 &boundary_points)
487{
488 // We check the Bounding Box is not degenerate
489 for (unsigned int i = 0; i < spacedim; ++i)
490 Assert(boundary_points.first[i] <= boundary_points.second[i],
491 ExcMessage("Bounding Box can't be created: the points' "
492 "order should be bottom left, top right!"));
493
494 this->boundary_points = boundary_points;
495}
496
497
498
499template <int spacedim, typename Number>
500template <class Container>
501inline BoundingBox<spacedim, Number>::BoundingBox(const Container &points)
502{
503 // Use the default constructor in case points is empty instead of setting
504 // things to +oo and -oo
505 if (points.size() > 0)
506 {
507 auto &min = boundary_points.first;
508 auto &max = boundary_points.second;
509 for (unsigned int d = 0; d < spacedim; ++d)
510 {
511 min[d] = std::numeric_limits<Number>::infinity();
512 max[d] = -std::numeric_limits<Number>::infinity();
513 }
514
515 for (const Point<spacedim, Number> &point : points)
516 for (unsigned int d = 0; d < spacedim; ++d)
517 {
518 min[d] = std::min(min[d], point[d]);
519 max[d] = std::max(max[d], point[d]);
520 }
521 }
522}
523
524
525
526template <int spacedim, typename Number>
527inline std::pair<Point<spacedim, Number>, Point<spacedim, Number>> &
529{
530 return this->boundary_points;
531}
532
533
534
535template <int spacedim, typename Number>
536inline const std::pair<Point<spacedim, Number>, Point<spacedim, Number>> &
538{
539 return this->boundary_points;
540}
541
542
543
544template <int spacedim, typename Number>
545inline bool
547 const BoundingBox<spacedim, Number> &box) const
548{
549 return boundary_points == box.boundary_points;
550}
551
552
553
554template <int spacedim, typename Number>
555inline bool
557 const BoundingBox<spacedim, Number> &box) const
558{
559 return boundary_points != box.boundary_points;
560}
561
562
563
564template <int spacedim, typename Number>
565inline void
566BoundingBox<spacedim, Number>::extend(const Number amount)
567{
568 for (unsigned int d = 0; d < spacedim; ++d)
569 {
570 boundary_points.first[d] -= amount;
571 boundary_points.second[d] += amount;
572 Assert(boundary_points.first[d] <= boundary_points.second[d],
573 ExcMessage("Bounding Box can't be shrunk this much: the points' "
574 "order should remain bottom left, top right."));
575 }
576}
577
578
579
580template <int spacedim, typename Number>
582BoundingBox<spacedim, Number>::create_extended(const Number amount) const
583{
584 // create and modify copy
585 auto bb = *this;
586 bb.extend(amount);
587
588 return bb;
589}
590
591
592
593template <int spacedim, typename Number>
596 const Number relative_amount) const
597{
598 // create and modify copy
599 auto bb = *this;
600
601 for (unsigned int d = 0; d < spacedim; ++d)
602 {
603 bb.boundary_points.first[d] -= relative_amount * side_length(d);
604 bb.boundary_points.second[d] += relative_amount * side_length(d);
605 Assert(bb.boundary_points.first[d] <= bb.boundary_points.second[d],
606 ExcMessage("Bounding Box can't be shrunk this much: the points' "
607 "order should remain bottom left, top right."));
608 }
609
610 return bb;
611}
612
613
614
615template <int spacedim, typename Number>
616template <class Archive>
617void
619 const unsigned int /*version*/)
620{
621 ar &boundary_points;
622}
623
624
625
626template <typename Number>
628{
630}
631
632
633
634template <typename Number>
636 const std::pair<Point<0, Number>, Point<0, Number>> &)
637{
639}
640
641
642
643template <typename Number>
644template <class Container>
645inline BoundingBox<0, Number>::BoundingBox(const Container &)
646{
648}
649
650
651
652#endif // DOXYGEN
654
655#endif
NeighborType
BoundingBox< dim, Number > create_unit_bounding_box()
BoundingBox(const std::pair< Point< 0, Number >, Point< 0, Number > > &)
BoundingBox(const Container &)
int coordinate_to_one_dim_higher(const int locked_coordinate, const int coordinate_in_dim)
std::pair< Point< spacedim, Number >, Point< spacedim, Number > > boundary_points
BoundingBox< 1, Number > bounds(const unsigned int direction) const
bool has_overlap_with(const BoundingBox< spacedim, Number > &other_bbox, const double tolerance=std::numeric_limits< Number >::epsilon()) const
BoundingBox(const std::pair< Point< spacedim, Number >, Point< spacedim, Number > > &boundary_points)
Point< spacedim, Number > center() const
BoundingBox()=default
Number lower_bound(const unsigned int direction) const
Number signed_distance(const Point< spacedim, Number > &point, const unsigned int direction) const
bool operator==(const BoundingBox< spacedim, Number > &box) const
void serialize(Archive &ar, const unsigned int version)
void merge_with(const BoundingBox< spacedim, Number > &other_bbox)
std::pair< Point< spacedim, Number >, Point< spacedim, Number > > & get_boundary_points()
bool point_inside(const Point< spacedim, Number > &p, const double tolerance=std::numeric_limits< Number >::epsilon()) const
double volume() const
void extend(const Number amount)
BoundingBox< spacedim, Number > create_extended_relative(const Number relative_amount) const
BoundingBox(const Container &points)
Point< spacedim, Number > real_to_unit(const Point< spacedim, Number > &point) const
const std::pair< Point< spacedim, Number >, Point< spacedim, Number > > & get_boundary_points() const
Number side_length(const unsigned int direction) const
BoundingBox< spacedim, Number > child(const unsigned int index) const
bool operator!=(const BoundingBox< spacedim, Number > &box) const
BoundingBox< spacedim, Number > & operator=(const BoundingBox< spacedim, Number > &t)=default
NeighborType get_neighbor_type(const BoundingBox< spacedim, Number > &other_bbox, const double tolerance=std::numeric_limits< Number >::epsilon()) const
BoundingBox< spacedim, Number > create_extended(const Number amount) const
Point< spacedim, Number > vertex(const unsigned int index) const
BoundingBox(const Point< spacedim, Number > &point)
BoundingBox< spacedim - 1, Number > cross_section(const unsigned int direction) const
BoundingBox(const BoundingBox< spacedim, Number > &box)=default
Number upper_bound(const unsigned int direction) const
Point< spacedim, Number > unit_to_real(const Point< spacedim, Number > &point) const
Definition point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
spacedim const Point< spacedim > & p
Definition grid_tools.h:990
for(unsigned int j=best_vertex+1;j< vertices.size();++j) if(vertices_to_use[j])
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)