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\}}\)
bounding_box.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2017 - 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_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>
25
27
32enum class NeighborType
33{
37 not_neighbors = 0,
38
45
52
64};
65
120template <int spacedim, typename Number = double>
122{
123public:
128 BoundingBox() = default;
129
137
145 template <class Container>
146 BoundingBox(const Container &points);
147
151 std::pair<Point<spacedim, Number>, Point<spacedim, Number>> &
153
157 const std::pair<Point<spacedim, Number>, Point<spacedim, Number>> &
159
163 bool
165
169 bool
171
179 get_neighbor_type(const BoundingBox<spacedim, Number> &other_bbox) const;
180
186 void
187 merge_with(const BoundingBox<spacedim, Number> &other_bbox);
188
195 bool
198 const double tolerance = std::numeric_limits<Number>::epsilon()) const;
199
210 void
211 extend(const Number &amount);
212
216 double
217 volume() const;
218
223 center() const;
224
228 Number
229 side_length(const unsigned int direction) const;
230
234 Number
235 lower_bound(const unsigned int direction) const;
236
240 Number
241 upper_bound(const unsigned int direction) const;
242
247 bounds(const unsigned int direction) const;
248
254 vertex(const unsigned int index) const;
255
261 child(const unsigned int index) const;
262
270 BoundingBox<spacedim - 1, Number>
271 cross_section(const unsigned int direction) const;
272
283
294
300 template <class Archive>
301 void
302 serialize(Archive &ar, const unsigned int version);
303
304private:
305 std::pair<Point<spacedim, Number>, Point<spacedim, Number>> boundary_points;
306};
307
313template <typename Number>
314class BoundingBox<0, Number>
315{
316public:
321
326
330 template <class Container>
331 BoundingBox(const Container &);
332};
333
334
340template <int dim, typename Number = double>
343
344
345namespace internal
346{
379 template <int dim>
380 inline int
381 coordinate_to_one_dim_higher(const int locked_coordinate,
382 const int coordinate_in_dim)
383 {
384 AssertIndexRange(locked_coordinate, dim + 1);
385 AssertIndexRange(coordinate_in_dim, dim);
386 return (locked_coordinate + coordinate_in_dim + 1) % (dim + 1);
387 }
388
389} // namespace internal
390
391/*------------------------ Inline functions: BoundingBox --------------------*/
392
393#ifndef DOXYGEN
394
395
396template <int spacedim, typename Number>
399 &boundary_points)
400{
401 // We check the Bounding Box is not degenerate
402 for (unsigned int i = 0; i < spacedim; ++i)
403 Assert(boundary_points.first[i] <= boundary_points.second[i],
404 ExcMessage("Bounding Box can't be created: the points' "
405 "order should be bottom left, top right!"));
406
407 this->boundary_points = boundary_points;
408}
409
410
411
412template <int spacedim, typename Number>
413template <class Container>
414inline BoundingBox<spacedim, Number>::BoundingBox(const Container &points)
415{
416 // Use the default constructor in case points is empty instead of setting
417 // things to +oo and -oo
418 if (points.size() > 0)
419 {
420 auto &min = boundary_points.first;
421 auto &max = boundary_points.second;
422 std::fill(min.begin_raw(),
423 min.end_raw(),
424 std::numeric_limits<Number>::infinity());
425 std::fill(max.begin_raw(),
426 max.end_raw(),
427 -std::numeric_limits<Number>::infinity());
428
429 for (const Point<spacedim, Number> &point : points)
430 for (unsigned int d = 0; d < spacedim; ++d)
431 {
432 min[d] = std::min(min[d], point[d]);
433 max[d] = std::max(max[d], point[d]);
434 }
435 }
436}
437
438
439
440template <int spacedim, typename Number>
441inline std::pair<Point<spacedim, Number>, Point<spacedim, Number>> &
443{
444 return this->boundary_points;
445}
446
447
448
449template <int spacedim, typename Number>
450inline const std::pair<Point<spacedim, Number>, Point<spacedim, Number>> &
452{
453 return this->boundary_points;
454}
455
456
457
458template <int spacedim, typename Number>
459inline bool
462{
463 return boundary_points == box.boundary_points;
464}
465
466
467
468template <int spacedim, typename Number>
469inline bool
472{
473 return boundary_points != box.boundary_points;
474}
475
476
477
478template <int spacedim, typename Number>
479inline void
480BoundingBox<spacedim, Number>::extend(const Number &amount)
481{
482 for (unsigned int d = 0; d < spacedim; ++d)
483 {
484 boundary_points.first[d] -= amount;
485 boundary_points.second[d] += amount;
486 Assert(boundary_points.first[d] <= boundary_points.second[d],
487 ExcMessage("Bounding Box can't be shrunk this much: the points' "
488 "order should remain bottom left, top right."));
489 }
490}
491
492
493template <int spacedim, typename Number>
494template <class Archive>
495void
497 const unsigned int /*version*/)
498{
499 ar &boundary_points;
500}
501
502
503
504template <typename Number>
506{
508}
509
510
511
512template <typename Number>
514 const std::pair<Point<0, Number>, Point<0, Number>> &)
515{
517}
518
519
520
521template <typename Number>
522template <class Container>
523inline BoundingBox<0, Number>::BoundingBox(const Container &)
524{
526}
527
528
529
530#endif // DOXYGEN
532
533#endif
NeighborType
Definition: bounding_box.h:33
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)
Definition: bounding_box.h:381
std::pair< Point< spacedim, Number >, Point< spacedim, Number > > boundary_points
Definition: bounding_box.h:305
BoundingBox< 1, Number > bounds(const unsigned int direction) const
BoundingBox(const std::pair< Point< spacedim, Number >, Point< spacedim, Number > > &boundary_points)
BoundingBox< dim, Number > create_unit_bounding_box()
NeighborType get_neighbor_type(const BoundingBox< spacedim, Number > &other_bbox) const
Definition: bounding_box.cc:60
Point< spacedim, Number > center() const
BoundingBox()=default
Number lower_bound(const unsigned int direction) const
void extend(const Number &amount)
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)
Definition: bounding_box.cc:44
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
Definition: bounding_box.cc:23
double volume() 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
Point< spacedim, Number > vertex(const unsigned int index) const
BoundingBox< spacedim - 1, Number > cross_section(const unsigned int direction) const
Number upper_bound(const unsigned int direction) const
Point< spacedim, Number > unit_to_real(const Point< spacedim, Number > &point) const
Definition: point.h:111
VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &x, const ::VectorizedArray< Number, width > &y)
VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &x, const ::VectorizedArray< Number, width > &y)
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
#define AssertIndexRange(index, range)
Definition: exceptions.h:1690
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1575
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition: utilities.cc:188
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > epsilon(const Tensor< 2, dim, Number > &Grad_u)
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)