Reference documentation for deal.II version 9.2.0
\(\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\}}\)
bounding_box.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2017 - 2020 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_base_bounding_box_h
17 #define dealii_base_bounding_box_h
18 
19 
20 #include <deal.II/base/config.h>
21 
23 #include <deal.II/base/point.h>
24 #include <deal.II/base/utilities.h>
25 
27 #include <boost/geometry/algorithms/envelope.hpp>
28 #include <boost/geometry/geometries/multi_point.hpp>
30 
32 
37 enum class NeighborType
38 {
42  not_neighbors = 0,
43 
49  simple_neighbors = 1,
50 
57 
69 };
70 
127 template <int spacedim, typename Number = double>
129 {
130 public:
135  BoundingBox() = default;
136 
143  &boundary_points);
144 
152  template <class Container>
153  BoundingBox(const Container &points);
154 
158  std::pair<Point<spacedim, Number>, Point<spacedim, Number>> &
160 
164  const std::pair<Point<spacedim, Number>, Point<spacedim, Number>> &
165  get_boundary_points() const;
166 
174  get_neighbor_type(const BoundingBox<spacedim, Number> &other_bbox) const;
175 
181  void
182  merge_with(const BoundingBox<spacedim, Number> &other_bbox);
183 
187  bool
188  point_inside(const Point<spacedim, Number> &p) const;
189 
200  void
201  extend(const Number &amount);
202 
206  double
207  volume() const;
208 
213  center() const;
214 
218  Number
219  side_length(const unsigned int direction) const;
220 
224  Number
225  lower_bound(const unsigned int direction) const;
226 
230  Number
231  upper_bound(const unsigned int direction) const;
232 
237  bounds(const unsigned int direction) const;
238 
244  vertex(const unsigned int index) const;
245 
251  child(const unsigned int index) const;
252 
260  BoundingBox<spacedim - 1, Number>
261  cross_section(const unsigned int direction) const;
262 
266  template <class Archive>
267  void
268  serialize(Archive &ar, const unsigned int version);
269 
270 private:
271  std::pair<Point<spacedim, Number>, Point<spacedim, Number>> boundary_points;
272 };
273 
279 template <typename Number>
280 class BoundingBox<0, Number>
281 {
282 public:
286  BoundingBox();
287 
291  BoundingBox(const std::pair<Point<0, Number>, Point<0, Number>> &);
292 
296  template <class Container>
297  BoundingBox(const Container &);
298 };
299 
300 
306 template <int dim, typename Number = double>
309 
310 
311 namespace internal
312 {
345  template <int dim>
346  inline unsigned int
347  coordinate_to_one_dim_higher(const unsigned int locked_coordinate,
348  const unsigned int coordiante_in_dim)
349  {
350  AssertIndexRange(locked_coordinate, dim + 1);
351  AssertIndexRange(coordiante_in_dim, dim);
352  return (locked_coordinate + coordiante_in_dim + 1) % (dim + 1);
353  }
354 
355 } // namespace internal
356 
357 /*------------------------ Inline functions: BoundingBox --------------------*/
358 
359 #ifndef DOXYGEN
360 
361 
362 template <int spacedim, typename Number>
365  &boundary_points)
366 {
367  // We check the Bounding Box is not degenerate
368  for (unsigned int i = 0; i < spacedim; ++i)
369  Assert(boundary_points.first[i] <= boundary_points.second[i],
370  ExcMessage("Bounding Box can't be created: the points' "
371  "order should be bottom left, top right!"));
372 
373  this->boundary_points = boundary_points;
374 }
375 
376 
377 
378 template <int spacedim, typename Number>
379 template <class Container>
380 inline BoundingBox<spacedim, Number>::BoundingBox(const Container &points)
381 {
382  boost::geometry::envelope(
383  boost::geometry::model::multi_point<Point<spacedim, Number>>(points.begin(),
384  points.end()),
385  *this);
386 }
387 
388 
389 
390 template <int spacedim, typename Number>
391 inline void
392 BoundingBox<spacedim, Number>::extend(const Number &amount)
393 {
394  for (unsigned int d = 0; d < spacedim; ++d)
395  {
396  boundary_points.first[d] -= amount;
397  boundary_points.second[d] += amount;
398  Assert(boundary_points.first[d] <= boundary_points.second[d],
399  ExcMessage("Bounding Box can't be shrunk this much: the points' "
400  "order should remain bottom left, top right."));
401  }
402 }
403 
404 
405 template <int spacedim, typename Number>
406 template <class Archive>
407 void
409  const unsigned int /*version*/)
410 {
411  ar &boundary_points;
412 }
413 
414 
415 
416 template <typename Number>
418 {
419  AssertThrow(false, ExcImpossibleInDim(0));
420 }
421 
422 
423 
424 template <typename Number>
426  const std::pair<Point<0, Number>, Point<0, Number>> &)
427 {
428  AssertThrow(false, ExcImpossibleInDim(0));
429 }
430 
431 
432 
433 template <typename Number>
434 template <class Container>
435 inline BoundingBox<0, Number>::BoundingBox(const Container &)
436 {
437  AssertThrow(false, ExcImpossibleInDim(0));
438 }
439 
440 
441 
442 #endif // DOXYGEN
444 
445 #endif
BoundingBox::volume
double volume() const
Definition: bounding_box.cc:168
BoundingBox::lower_bound
Number lower_bound(const unsigned int direction) const
Definition: bounding_box.cc:180
utilities.h
BoundingBox::child
BoundingBox< spacedim, Number > child(const unsigned int index) const
Definition: bounding_box.cc:259
AssertIndexRange
#define AssertIndexRange(index, range)
Definition: exceptions.h:1649
NeighborType::attached_neighbors
@ attached_neighbors
NeighborType
NeighborType
Definition: bounding_box.h:37
NeighborType::mergeable_neighbors
@ mergeable_neighbors
BoundingBox::center
Point< spacedim, Number > center() const
Definition: bounding_box.cc:202
Physics::Elasticity::Kinematics::d
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
BoundingBox::extend
void extend(const Number &amount)
BoundingBox
Definition: bounding_box.h:128
DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
Definition: config.h:409
NeighborType::not_neighbors
@ not_neighbors
BoundingBox::merge_with
void merge_with(const BoundingBox< spacedim, Number > &other_bbox)
Definition: bounding_box.cc:43
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
BoundingBox::upper_bound
Number upper_bound(const unsigned int direction) const
Definition: bounding_box.cc:191
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
BoundingBox::BoundingBox
BoundingBox()=default
BoundingBox::cross_section
BoundingBox< spacedim - 1, Number > cross_section(const unsigned int direction) const
Definition: bounding_box.cc:295
exceptions.h
BoundingBox::create_unit_bounding_box
BoundingBox< dim, Number > create_unit_bounding_box()
Definition: bounding_box.cc:320
BoundingBox::point_inside
bool point_inside(const Point< spacedim, Number > &p) const
Definition: bounding_box.cc:22
BoundingBox::boundary_points
std::pair< Point< spacedim, Number >, Point< spacedim, Number > > boundary_points
Definition: bounding_box.h:271
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
BoundingBox::get_neighbor_type
NeighborType get_neighbor_type(const BoundingBox< spacedim, Number > &other_bbox) const
Definition: bounding_box.cc:59
StandardExceptions::ExcImpossibleInDim
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
Point
Definition: point.h:111
config.h
internal
Definition: aligned_vector.h:369
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
BoundingBox::side_length
Number side_length(const unsigned int direction) const
Definition: bounding_box.cc:230
BoundingBox::bounds
BoundingBox< 1, Number > bounds(const unsigned int direction) const
Definition: bounding_box.cc:215
BoundingBox::get_boundary_points
std::pair< Point< spacedim, Number >, Point< spacedim, Number > > & get_boundary_points()
Definition: bounding_box.cc:150
AssertThrow
#define AssertThrow(cond, exc)
Definition: exceptions.h:1531
BoundingBox::vertex
Point< spacedim, Number > vertex(const unsigned int index) const
Definition: bounding_box.cc:241
BoundingBox::coordinate_to_one_dim_higher
unsigned int coordinate_to_one_dim_higher(const unsigned int locked_coordinate, const unsigned int coordiante_in_dim)
Definition: bounding_box.h:347
NeighborType::simple_neighbors
@ simple_neighbors
DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
Definition: config.h:372
BoundingBox::serialize
void serialize(Archive &ar, const unsigned int version)
point.h