Reference documentation for deal.II version 9.6.0
|
#include <deal.II/base/bounding_box.h>
Public Member Functions | |
BoundingBox ()=default | |
BoundingBox (const BoundingBox< spacedim, Number > &box)=default | |
BoundingBox< spacedim, Number > & | operator= (const BoundingBox< spacedim, Number > &t)=default |
BoundingBox (const Point< spacedim, Number > &point) | |
BoundingBox (const std::pair< Point< spacedim, Number >, Point< spacedim, Number > > &boundary_points) | |
template<class Container > | |
BoundingBox (const Container &points) | |
std::pair< Point< spacedim, Number >, Point< spacedim, Number > > & | get_boundary_points () |
const std::pair< Point< spacedim, Number >, Point< spacedim, Number > > & | get_boundary_points () const |
bool | operator== (const BoundingBox< spacedim, Number > &box) const |
bool | operator!= (const BoundingBox< spacedim, Number > &box) const |
bool | has_overlap_with (const BoundingBox< spacedim, Number > &other_bbox, const double tolerance=std::numeric_limits< Number >::epsilon()) const |
NeighborType | get_neighbor_type (const BoundingBox< spacedim, Number > &other_bbox, const double tolerance=std::numeric_limits< Number >::epsilon()) const |
void | merge_with (const BoundingBox< spacedim, Number > &other_bbox) |
bool | point_inside (const Point< spacedim, Number > &p, const double tolerance=std::numeric_limits< Number >::epsilon()) const |
void | extend (const Number amount) |
BoundingBox< spacedim, Number > | create_extended (const Number amount) const |
BoundingBox< spacedim, Number > | create_extended_relative (const Number relative_amount) const |
double | volume () const |
Point< spacedim, Number > | center () const |
Number | side_length (const unsigned int direction) const |
Number | lower_bound (const unsigned int direction) const |
Number | upper_bound (const unsigned int direction) const |
BoundingBox< 1, Number > | bounds (const unsigned int direction) const |
Point< spacedim, Number > | vertex (const unsigned int index) const |
BoundingBox< spacedim, Number > | child (const unsigned int index) const |
BoundingBox< spacedim - 1, Number > | cross_section (const unsigned int direction) const |
Point< spacedim, Number > | real_to_unit (const Point< spacedim, Number > &point) const |
Point< spacedim, Number > | unit_to_real (const Point< spacedim, Number > &point) const |
Number | signed_distance (const Point< spacedim, Number > &point, const unsigned int direction) const |
Number | signed_distance (const Point< spacedim, Number > &point) const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Private Attributes | |
std::pair< Point< spacedim, Number >, Point< spacedim, Number > > | boundary_points |
Related Symbols | |
(Note that these are not member symbols.) | |
template<int dim, typename Number = double> | |
BoundingBox< dim, Number > | create_unit_bounding_box () |
template<int dim> | |
int | coordinate_to_one_dim_higher (const int locked_coordinate, const int coordinate_in_dim) |
A class that represents a box of arbitrary dimension spacedim
and with sides parallel to the coordinate axes, that is, a region
\[ [x_0^L, x_0^U] \times ... \times [x_{spacedim-1}^L, x_{spacedim-1}^U], \]
where \((x_0^L , ..., x_{spacedim-1}^L)\) and \((x_0^U , ..., x_{spacedim-1}^U)\) denote the two vertices (bottom left and top right) which are used to represent the box. The quantities \(x_k^L\) and \(x_k^U\) denote the "lower" and "upper" bounds of values that are within the box for each coordinate direction \(k\).
Geometrically, a bounding box is thus:
Bounding boxes are, for example, useful in parallel distributed meshes to give a general description of the owners of each portion of the mesh. More generally, bounding boxes are often used to roughly describe a region of space in which an object is contained; if a candidate point is not within the bounding box (a test that is cheap to execute), then it is not necessary to perform an expensive test whether the candidate point is in fact inside the object itself. Bounding boxes are therefore often used as a first, cheap rejection test before more detailed checks. As such, bounding boxes serve many of the same purposes as the convex hull, for which it is also relatively straightforward to compute whether a point is inside or outside, though not quite as cheap as for the bounding box.
Taking the cross section of a BoundingBox<spacedim> orthogonal to a given direction gives a box in one dimension lower: BoundingBox<spacedim - 1>. In 3d, the 2 coordinates of the cross section of BoundingBox<3> can be ordered in 2 different ways. That is, if we take the cross section orthogonal to the y direction we could either order a 3d-coordinate into a 2d-coordinate as \((x,z)\) or as \((z,x)\). This class uses the second convention, corresponding to the coordinates being ordered cyclicly \(x \rightarrow y \rightarrow z \rightarrow x \rightarrow ... \) To be precise, if we take a cross section:
Orthogonal to | Cross section coordinates ordered as |
---|---|
x | (y, z) |
y | (z, x) |
z | (x, y) |
This is according to the convention set by the function coordinate_to_one_dim_higher
.
Definition at line 135 of file bounding_box.h.
|
default |
Standard constructor. Creates an object that corresponds to an empty box, i.e. a degenerate box with both points being the origin.
|
default |
Standard copy constructor operator.
BoundingBox< spacedim, Number >::BoundingBox | ( | const Point< spacedim, Number > & | point | ) |
Standard constructor for an empty box around a point point
.
BoundingBox< spacedim, Number >::BoundingBox | ( | const std::pair< Point< spacedim, Number >, Point< spacedim, Number > > & | boundary_points | ) |
Standard constructor for non-empty boxes: it uses a pair of points which describe the box: one for the bottom and one for the top corner.
BoundingBox< spacedim, Number >::BoundingBox | ( | const Container & | points | ) |
Construct the bounding box that encloses all the points in the given container.
The constructor supports any Container that provides begin() and end() iterators to Point<spacedim, Number> elements.
|
default |
Standard copy assignment operator.
std::pair< Point< spacedim, Number >, Point< spacedim, Number > > & BoundingBox< spacedim, Number >::get_boundary_points | ( | ) |
Return a reference to the boundary_points
const std::pair< Point< spacedim, Number >, Point< spacedim, Number > > & BoundingBox< spacedim, Number >::get_boundary_points | ( | ) | const |
Return a const reference to the boundary_points
bool BoundingBox< spacedim, Number >::operator== | ( | const BoundingBox< spacedim, Number > & | box | ) | const |
Test for equality.
bool BoundingBox< spacedim, Number >::operator!= | ( | const BoundingBox< spacedim, Number > & | box | ) | const |
Test for inequality.
bool BoundingBox< spacedim, Number >::has_overlap_with | ( | const BoundingBox< spacedim, Number > & | other_bbox, |
const double | tolerance = std::numeric_limits<Number>::epsilon() ) const |
Check if the current object and other_bbox
are neighbors, i.e. if the boxes have dimension spacedim, check if their intersection is non empty.
Definition at line 61 of file bounding_box.cc.
NeighborType BoundingBox< spacedim, Number >::get_neighbor_type | ( | const BoundingBox< spacedim, Number > & | other_bbox, |
const double | tolerance = std::numeric_limits<Number>::epsilon() ) const |
Check which NeighborType other_bbox
is to the current object.
Definition at line 79 of file bounding_box.cc.
void BoundingBox< spacedim, Number >::merge_with | ( | const BoundingBox< spacedim, Number > & | other_bbox | ) |
Enlarge the current object so that it contains other_bbox
. If the current object already contains other_bbox
then it is not changed by this function.
Definition at line 45 of file bounding_box.cc.
bool BoundingBox< spacedim, Number >::point_inside | ( | const Point< spacedim, Number > & | p, |
const double | tolerance = std::numeric_limits<Number>::epsilon() ) const |
Return true if the point is inside the Bounding Box, false otherwise. The parameter tolerance
is a factor by which the bounding box is enlarged relative to the dimensions of the bounding box in order to determine in a numerically robust way whether the point is inside.
Definition at line 25 of file bounding_box.cc.
void BoundingBox< spacedim, Number >::extend | ( | const Number | amount | ) |
Increase (or decrease) the size of the bounding box by the given amount. After calling this method, the lower left corner of the bounding box will have each coordinate decreased by amount
, and the upper right corner of the bounding box will have each coordinate increased by amount
.
If you call this method with a negative number, and one of the axes of the original bounding box is smaller than amount/2, the method will trigger an assertion.
BoundingBox< spacedim, Number > BoundingBox< spacedim, Number >::create_extended | ( | const Number | amount | ) | const |
The same as above with the difference that a new BoundingBox instance is created without changing the current object.
BoundingBox< spacedim, Number > BoundingBox< spacedim, Number >::create_extended_relative | ( | const Number | relative_amount | ) | const |
Increase (or decrease) each side of the bounding box by the given relative_amount
.
After calling this method, the lower left corner of the bounding box will have each coordinate decreased by relative_amount
* side_length(direction), and the upper right corner of the bounding box will have each coordinate increased by relative_amount
* side_length(direction).
If you call this method with a negative number, and one of the axes of the original bounding box is smaller than relative_amount * side_length(direction) / 2, the method will trigger an assertion.
double BoundingBox< spacedim, Number >::volume | ( | ) | const |
Compute the volume (i.e. the dim-dimensional measure) of the BoundingBox.
Definition at line 159 of file bounding_box.cc.
Point< spacedim, Number > BoundingBox< spacedim, Number >::center | ( | ) | const |
Returns the point in the center of the box.
Definition at line 193 of file bounding_box.cc.
Number BoundingBox< spacedim, Number >::side_length | ( | const unsigned int | direction | ) | const |
Returns the side length of the box in direction
.
Definition at line 221 of file bounding_box.cc.
Number BoundingBox< spacedim, Number >::lower_bound | ( | const unsigned int | direction | ) | const |
Return the lower bound of the box in direction
.
Definition at line 171 of file bounding_box.cc.
Number BoundingBox< spacedim, Number >::upper_bound | ( | const unsigned int | direction | ) | const |
Return the upper bound of the box in direction
.
Definition at line 182 of file bounding_box.cc.
BoundingBox< 1, Number > BoundingBox< spacedim, Number >::bounds | ( | const unsigned int | direction | ) | const |
Return the bounds of the box in direction
, as a one-dimensional box.
Definition at line 206 of file bounding_box.cc.
Point< spacedim, Number > BoundingBox< spacedim, Number >::vertex | ( | const unsigned int | index | ) | const |
Returns the indexth vertex of the box. Vertex is meant in the same way as for a cell, so that index
\(\in [0, 2^{\text{dim}} - 1]\).
Definition at line 232 of file bounding_box.cc.
BoundingBox< spacedim, Number > BoundingBox< spacedim, Number >::child | ( | const unsigned int | index | ) | const |
Returns the indexth child of the box. Child is meant in the same way as for a cell.
Definition at line 250 of file bounding_box.cc.
BoundingBox< spacedim - 1, Number > BoundingBox< spacedim, Number >::cross_section | ( | const unsigned int | direction | ) | const |
Returns the cross section of the box orthogonal to direction
. This is a box in one dimension lower.
BoundingBox<0>
is not implemented. Definition at line 286 of file bounding_box.cc.
Point< spacedim, Number > BoundingBox< spacedim, Number >::real_to_unit | ( | const Point< spacedim, Number > & | point | ) | const |
Apply the affine transformation that transforms this BoundingBox to a unit BoundingBox object.
If \(B\) is this bounding box, and \(\hat{B}\) is the unit bounding box, compute the affine mapping that satisfies \(G(B) = \hat{B}\) and apply it to point
.
Definition at line 311 of file bounding_box.cc.
Point< spacedim, Number > BoundingBox< spacedim, Number >::unit_to_real | ( | const Point< spacedim, Number > & | point | ) | const |
Apply the affine transformation that transforms the unit BoundingBox object to this object.
If \(B\) is this bounding box, and \(\hat{B}\) is the unit bounding box, compute the affine mapping that satisfies \(F(\hat{B}) = B\) and apply it to point
.
Definition at line 326 of file bounding_box.cc.
Number BoundingBox< spacedim, Number >::signed_distance | ( | const Point< spacedim, Number > & | point, |
const unsigned int | direction ) const |
Returns the signed distance from a point
orthogonal to the bounds of the box in direction
. The signed distance is negative for points inside the interval described by the bounds of the rectangle in the respective direction, zero for points on the interval boundary and positive for points outside.
Definition at line 340 of file bounding_box.cc.
Number BoundingBox< spacedim, Number >::signed_distance | ( | const Point< spacedim, Number > & | point | ) | const |
Returns the signed distance from a point
to the bounds of the box. The signed distance is negative for points inside the rectangle, zero for points on the rectangle and positive for points outside the rectangle.
Definition at line 359 of file bounding_box.cc.
void BoundingBox< spacedim, Number >::serialize | ( | Archive & | ar, |
const unsigned int | version ) |
Write or read the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.
|
related |
Returns the unit box \([0,1]^\text{dim}\).
Definition at line 395 of file bounding_box.cc.
|
related |
This function defines a convention for how coordinates in dim dimensions should translate to the coordinates in dim + 1 dimensions, when one of the coordinates in dim + 1 dimensions is locked to a given value.
The convention is the following: Starting from the locked coordinate we store the lower dimensional coordinates consecutively and wrap around when going over the dimension. This relationship is, in 2d,
locked in 2D | 1d coordinate | 2d coordinate |
---|---|---|
x0 | (a) | (x0, a) |
x1 | (a) | (a , x1) |
and, in 3d,
locked in 3D | 2d coordinates | 3d coordinates |
---|---|---|
x0 | (a, b) | (x0, a, b) |
x1 | (a, b) | ( b, x1, a) |
x2 | (a, b) | ( a, b, x2) |
Given a locked coordinate, this function maps a coordinate index in dim dimension to a coordinate index in dim + 1 dimensions.
locked_coordinate | should be in the range [0, dim+1). |
coordinate_in_dim | should be in the range [0, dim). |
Definition at line 460 of file bounding_box.h.
|
private |
Definition at line 384 of file bounding_box.h.