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\}}\)
Public Member Functions | Private Attributes | Related Functions | List of all members
BoundingBox< spacedim, Number > Class Template Reference

#include <deal.II/base/bounding_box.h>

Inheritance diagram for BoundingBox< spacedim, Number >:
[legend]

Public Member Functions

 BoundingBox ()=default
 
 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
 
NeighborType get_neighbor_type (const BoundingBox< spacedim, Number > &other_bbox) const
 
void merge_with (const BoundingBox< spacedim, Number > &other_bbox)
 
bool point_inside (const Point< spacedim, Number > &p) const
 
void extend (const Number &amount)
 
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
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Private Attributes

std::pair< Point< spacedim, Number >, Point< spacedim, Number > > boundary_points
 

Related Functions

(Note that these are not member functions.)

template<int dim, typename Number = double>
BoundingBox< dim, Number > create_unit_bounding_box ()
 
template<int dim>
unsigned int coordinate_to_one_dim_higher (const unsigned int locked_coordinate, const unsigned int coordiante_in_dim)
 

Detailed Description

template<int spacedim, typename Number = double>
class BoundingBox< spacedim, Number >

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.

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.

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.

Author
Giovanni Alzetta, 2017, Simon Sticko 2020.

Definition at line 128 of file bounding_box.h.

Constructor & Destructor Documentation

◆ BoundingBox() [1/3]

template<int spacedim, typename Number = double>
BoundingBox< spacedim, Number >::BoundingBox ( )
default

Standard constructor. Creates an object that corresponds to an empty box, i.e. a degenerate box with both points being the origin.

◆ BoundingBox() [2/3]

template<int spacedim, typename Number = double>
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() [3/3]

template<int spacedim, typename Number = double>
template<class Container >
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.

Member Function Documentation

◆ get_boundary_points() [1/2]

template<int spacedim, typename Number >
const std::pair< Point< spacedim, Number >, Point< spacedim, Number > > & BoundingBox< spacedim, Number >::get_boundary_points

Return a reference to the boundary_points

Definition at line 150 of file bounding_box.cc.

◆ get_boundary_points() [2/2]

template<int spacedim, typename Number = double>
const std::pair<Point<spacedim, Number>, Point<spacedim, Number> >& BoundingBox< spacedim, Number >::get_boundary_points ( ) const

Return a const reference to the boundary_points

◆ get_neighbor_type()

template<int spacedim, typename Number >
NeighborType BoundingBox< spacedim, Number >::get_neighbor_type ( const BoundingBox< spacedim, Number > &  other_bbox) 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.

Return an enumerator of type NeighborType.

Definition at line 59 of file bounding_box.cc.

◆ merge_with()

template<int spacedim, typename Number >
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 43 of file bounding_box.cc.

◆ point_inside()

template<int spacedim, typename Number >
bool BoundingBox< spacedim, Number >::point_inside ( const Point< spacedim, Number > &  p) const

Return true if the point is inside the Bounding Box, false otherwise.

Definition at line 22 of file bounding_box.cc.

◆ extend()

template<int spacedim, typename Number = double>
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.

◆ volume()

template<int spacedim, typename Number >
double BoundingBox< spacedim, Number >::volume

Compute the volume (i.e. the dim-dimensional measure) of the BoundingBox.

Definition at line 168 of file bounding_box.cc.

◆ center()

template<int spacedim, typename Number >
Point< spacedim, Number > BoundingBox< spacedim, Number >::center

Returns the point in the center of the box.

Definition at line 202 of file bounding_box.cc.

◆ side_length()

template<int spacedim, typename Number >
Number BoundingBox< spacedim, Number >::side_length ( const unsigned int  direction) const

Returns the side length of the box in direction.

Definition at line 230 of file bounding_box.cc.

◆ lower_bound()

template<int spacedim, typename Number >
Number BoundingBox< spacedim, Number >::lower_bound ( const unsigned int  direction) const

Return the lower bound of the box in direction.

Definition at line 180 of file bounding_box.cc.

◆ upper_bound()

template<int spacedim, typename Number >
Number BoundingBox< spacedim, Number >::upper_bound ( const unsigned int  direction) const

Return the upper bound of the box in direction.

Definition at line 191 of file bounding_box.cc.

◆ bounds()

template<int spacedim, typename Number >
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 215 of file bounding_box.cc.

◆ vertex()

template<int spacedim, typename Number >
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 241 of file bounding_box.cc.

◆ child()

template<int spacedim, typename Number >
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 259 of file bounding_box.cc.

◆ cross_section()

template<int spacedim, typename Number >
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.

Note
Calling this method in 1D will result in an exception since BoundingBox<0> is not implemented.

Definition at line 295 of file bounding_box.cc.

◆ serialize()

template<int spacedim, typename Number = double>
template<class Archive >
void BoundingBox< spacedim, Number >::serialize ( Archive &  ar,
const unsigned int  version 
)

Boost serialization function

Friends And Related Function Documentation

◆ create_unit_bounding_box()

template<int dim, typename Number = double>
BoundingBox< dim, Number > create_unit_bounding_box ( )
related

Returns the unit box \([0,1]^\text{dim}\).

Definition at line 320 of file bounding_box.cc.

◆ coordinate_to_one_dim_higher()

template<int dim>
unsigned int coordinate_to_one_dim_higher ( const unsigned int  locked_coordinate,
const unsigned int  coordiante_in_dim 
)
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.

Parameters
locked_coordinateshould be in the range [0, dim+1).
coordiante_in_dimshould be in the range [0, dim).
Returns
A coordinate index in the range [0, dim+1)

Definition at line 347 of file bounding_box.h.

Member Data Documentation

◆ boundary_points

template<int spacedim, typename Number = double>
std::pair<Point<spacedim, Number>, Point<spacedim, Number> > BoundingBox< spacedim, Number >::boundary_points
private

Definition at line 271 of file bounding_box.h.


The documentation for this class was generated from the following files:
LAPACKSupport::V
static const char V
Definition: lapack_support.h:175