A bounding box defined in terms of the points at its extremities. More...
#include <OpenFOAM/boundBox.H>
A bounding box defined in terms of the points at its extremities.
Definition at line 54 of file boundBox.H.
Public Member Functions | |
boundBox () | |
Construct null, setting points to zero.
| |
boundBox (const point &min, const point &max) | |
Construct from components.
| |
boundBox (const pointField &, const bool doReduce=true) | |
Construct as the bounding box of the given pointField.
| |
boundBox (const tmp< pointField > &, const bool doReduce=true) | |
Construct as the bounding box of the given temporary pointField.
| |
boundBox (Istream &) | |
Construct from Istream.
| |
const point & | min () const |
Minimum describing the bounding box.
| |
const point & | max () const |
Maximum describing the bounding box.
| |
point & | min () |
Minimum describing the bounding box, non-const access.
| |
point & | max () |
Maximum describing the bounding box, non-const access.
| |
point | midpoint () const |
The midpoint of the bounding box.
| |
vector | span () const |
The bounding box span (from minimum to maximum)
| |
scalar | mag () const |
The magnitude of the bounding box span.
| |
scalar | minDim () const |
Smallest length/height/width dimension.
| |
scalar | maxDim () const |
Largest length/height/width dimension.
| |
scalar | avgDim () const |
Average length/height/width dimension.
| |
bool | overlaps (const boundBox &bb) const |
Overlaps/touches boundingBox?
| |
bool | contains (const point &pt) const |
Contains point? (inside or on edge)
| |
bool | containsInside (const point &pt) const |
Contains point? (inside only)
| |
Static Public Attributes | |
static const scalar | great |
The great value used for greatBox and invertedBox.
| |
static const boundBox | greatBox |
A very large boundBox: min/max == -/+ VGREAT.
| |
static const boundBox | invertedBox |
A very large inverted boundBox: min/max == +/- VGREAT.
| |
Friends | |
bool | operator== (const boundBox &a, const boundBox &b) |
bool | operator!= (const boundBox &a, const boundBox &b) |
Istream & | operator>> (Istream &is, boundBox &) |
Ostream & | operator<< (Ostream &os, const boundBox &) |
boundBox | ( | ) | [inline]
|
Construct null, setting points to zero.
Definition at line 84 of file boundBox.H.
Construct from components.
Definition at line 91 of file boundBox.H.
boundBox | ( | const pointField & | points, |
const bool | doReduce = true
|
||
) |
Construct as the bounding box of the given pointField.
Does parallel communication (doReduce = true)
Definition at line 87 of file boundBox.C.
boundBox | ( | const tmp< pointField > & | points, |
const bool | doReduce = true
|
||
) |
Construct as the bounding box of the given temporary pointField.
Does parallel communication (doReduce = true)
Definition at line 96 of file boundBox.C.
References tmp< T >::clear(), and points.
const point& min | ( | ) | const [inline]
|
Minimum describing the bounding box.
Definition at line 114 of file boundBox.H.
Referenced by triSurfaceMesh::calcBounds(), Foam::meshTools::constrainToMeshCentre(), treeBoundBox::contains(), triSurfaceMesh::edgeTree(), treeBoundBox::extend(), treeLeaf< Type >::findNearest(), treeDataPoint::findNearest(), treeDataEdge::findNearest(), indexedOctree< Type >::findNearest(), octreeDataTriSurface::findTightest(), octreeDataFace::findTightest(), octreeDataEdges::findTightest(), octreeDataCell::findTightest(), octreeDataFaceList::findTightest(), boundaryMesh::getNearest(), triangleFuncs::intersectBb(), treeDataTriSurface::intersects(), treeDataTriSurface::overlaps(), treeBoundBox::subBbox(), triSurfaceMesh::tree(), treeBoundBox::treeBoundBox(), triSurfaceSearch::triSurfaceSearch(), topoSet::writeDebug(), triSurface::writeStats(), and distributedTriSurfaceMesh::writeStats().
const point& max | ( | ) | const [inline]
|
Maximum describing the bounding box.
Definition at line 120 of file boundBox.H.
Referenced by triSurfaceMesh::calcBounds(), Foam::meshTools::constrainToMeshCentre(), treeBoundBox::contains(), triSurfaceMesh::edgeTree(), treeBoundBox::extend(), treeLeaf< Type >::findNearest(), treeDataPoint::findNearest(), treeDataEdge::findNearest(), indexedOctree< Type >::findNearest(), octreeDataTriSurface::findTightest(), octreeDataFace::findTightest(), octreeDataEdges::findTightest(), octreeDataCell::findTightest(), octreeDataFaceList::findTightest(), boundaryMesh::getNearest(), triangleFuncs::intersectBb(), treeDataTriSurface::intersects(), Kmesh::Kmesh(), treeDataTriSurface::overlaps(), treeBoundBox::subBbox(), triSurfaceMesh::tree(), treeBoundBox::treeBoundBox(), topoSet::writeDebug(), triSurface::writeStats(), and distributedTriSurfaceMesh::writeStats().
point& min | ( | ) | [inline]
|
Minimum describing the bounding box, non-const access.
Definition at line 126 of file boundBox.H.
point& max | ( | ) | [inline]
|
Maximum describing the bounding box, non-const access.
Definition at line 132 of file boundBox.H.
point midpoint | ( | ) | const [inline]
|
The midpoint of the bounding box.
Definition at line 138 of file boundBox.H.
Referenced by treeNode< Type >::findBox(), treeNode< Type >::setSubNodeType(), and treeNode< Type >::writeOBJ().
vector span | ( | ) | const [inline]
|
The bounding box span (from minimum to maximum)
Definition at line 144 of file boundBox.H.
Referenced by boundBox::avgDim(), treeBoundBox::extend(), treeDataTriSurface::getVolumeType(), Kmesh::Kmesh(), boundBox::maxDim(), boundBox::minDim(), and displacementFvMotionSolver::updateMesh().
scalar mag | ( | ) | const [inline]
|
The magnitude of the bounding box span.
Definition at line 150 of file boundBox.H.
Referenced by patchProbes::findElements(), meshSearch::findNearestBoundaryFace(), and globalMeshData::updateMesh().
scalar minDim | ( | ) | const [inline]
|
Smallest length/height/width dimension.
Definition at line 156 of file boundBox.H.
References Foam::cmptMin(), and boundBox::span().
scalar maxDim | ( | ) | const [inline]
|
Largest length/height/width dimension.
Definition at line 162 of file boundBox.H.
References Foam::cmptMax(), and boundBox::span().
scalar avgDim | ( | ) | const [inline]
|
Average length/height/width dimension.
Definition at line 168 of file boundBox.H.
References Foam::cmptAv(), and boundBox::span().
Referenced by boundaryMesh::getNearest().
bool overlaps | ( | const boundBox & | bb ) | const [inline]
|
Overlaps/touches boundingBox?
Definition at line 177 of file boundBox.H.
References Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().
Referenced by treeBoundBox::overlaps().
bool contains | ( | const point & | pt ) | const [inline]
|
Contains point? (inside or on edge)
Reimplemented in treeBoundBox.
Definition at line 188 of file boundBox.H.
References Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().
Referenced by treeBoundBox::contains().
bool containsInside | ( | const point & | pt ) | const [inline]
|
Contains point? (inside only)
Definition at line 199 of file boundBox.H.
References Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().
Definition at line 212 of file boundBox.H.
Definition at line 217 of file boundBox.H.
const Foam::scalar great [static]
|
The great value used for greatBox and invertedBox.
Reimplemented in treeBoundBox.
Definition at line 72 of file boundBox.H.
const Foam::boundBox greatBox [static]
|
A very large boundBox: min/max == -/+ VGREAT.
Reimplemented in treeBoundBox.
Definition at line 75 of file boundBox.H.
const Foam::boundBox invertedBox [static]
|
A very large inverted boundBox: min/max == +/- VGREAT.
Reimplemented in treeBoundBox.
Definition at line 78 of file boundBox.H.
Referenced by triSurfaceMesh::calcBounds(), and triSurface::writeStats().