boundBox Class Reference

A bounding box defined in terms of the points at its extremities. More...

Inheritance diagram for boundBox:
Collaboration diagram for boundBox:

Public Member Functions

 boundBox ()
 Construct null, setting points to zero. More...
 
 boundBox (const point &min, const point &max)
 Construct from components. More...
 
 boundBox (const UList< point > &, const bool doReduce=true)
 Construct as the bounding box of the given points. More...
 
 boundBox (const tmp< pointField > &, const bool doReduce=true)
 Construct as the bounding box of the given temporary pointField. More...
 
 boundBox (const UList< point > &, const labelUList &indices, const bool doReduce=true)
 Construct bounding box as subset of the pointField. More...
 
template<unsigned Size>
 boundBox (const UList< point > &, const FixedList< label, Size > &indices, const bool doReduce=true)
 Construct bounding box as subset of the pointField. More...
 
 boundBox (Istream &)
 Construct from Istream. More...
 
const pointmin () const
 Minimum describing the bounding box. More...
 
const pointmax () const
 Maximum describing the bounding box. More...
 
pointmin ()
 Minimum describing the bounding box, non-const access. More...
 
pointmax ()
 Maximum describing the bounding box, non-const access. More...
 
point midpoint () const
 The midpoint of the bounding box. More...
 
vector span () const
 The bounding box span (from minimum to maximum) More...
 
scalar mag () const
 The magnitude of the bounding box span. More...
 
scalar volume () const
 The volume of the bound box. More...
 
scalar minDim () const
 Smallest length/height/width dimension. More...
 
scalar maxDim () const
 Largest length/height/width dimension. More...
 
scalar avgDim () const
 Average length/height/width dimension. More...
 
tmp< pointFieldpoints () const
 Return corner points in an order corresponding to a 'hex' cell. More...
 
void inflate (const scalar s)
 Inflate box by factor*mag(span) in all dimensions. More...
 
bool overlaps (const boundBox &) const
 Overlaps/touches boundingBox? More...
 
bool overlaps (const point &, const scalar radiusSqr) const
 Overlaps boundingSphere (centre + sqr(radius))? More...
 
bool contains (const point &) const
 Contains point? (inside or on edge) More...
 
bool contains (const boundBox &) const
 Fully contains other boundingBox? More...
 
bool containsInside (const point &) const
 Contains point? (inside only) More...
 
bool contains (const UList< point > &) const
 Contains all of the points? (inside or on edge) More...
 
bool contains (const UList< point > &, const labelUList &indices) const
 Contains all of the points? (inside or on edge) More...
 
template<unsigned Size>
bool contains (const UList< point > &, const FixedList< label, Size > &indices) const
 Contains all of the points? (inside or on edge) More...
 
bool containsAny (const UList< point > &) const
 Contains any of the points? (inside or on edge) More...
 
bool containsAny (const UList< point > &, const labelUList &indices) const
 Contains any of the points? (inside or on edge) More...
 
template<unsigned Size>
bool containsAny (const UList< point > &, const FixedList< label, Size > &indices) const
 Contains any of the points? (inside or on edge) More...
 
point nearest (const point &) const
 Return the nearest point on the boundBox to the supplied point. More...
 

Static Public Member Functions

static faceList faces ()
 Return faces with correct point order. More...
 

Static Public Attributes

static const scalar great
 The great value used for greatBox and invertedBox. More...
 
static const boundBox greatBox
 A very large boundBox: min/max == -/+ vGreat. More...
 
static const boundBox invertedBox
 A very large inverted boundBox: min/max == +/- vGreat. More...
 

Friends

bool operator== (const boundBox &, const boundBox &)
 
bool operator!= (const boundBox &, const boundBox &)
 
Istreamoperator>> (Istream &, boundBox &)
 
Ostreamoperator<< (Ostream &, const boundBox &)
 

Detailed Description

A bounding box defined in terms of the points at its extremities.

Definition at line 58 of file boundBox.H.

Constructor & Destructor Documentation

◆ boundBox() [1/7]

boundBox ( )
inline

Construct null, setting points to zero.

Definition at line 32 of file boundBoxI.H.

Referenced by boundBox::boundBox().

Here is the caller graph for this function:

◆ boundBox() [2/7]

boundBox ( const point min,
const point max 
)
inline

Construct from components.

Definition at line 39 of file boundBoxI.H.

◆ boundBox() [3/7]

boundBox ( const UList< point > &  points,
const bool  doReduce = true 
)

Construct as the bounding box of the given points.

Does parallel communication (doReduce = true)

Definition at line 88 of file boundBox.C.

◆ boundBox() [4/7]

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 97 of file boundBox.C.

References boundBox::boundBox(), tmp< T >::clear(), and boundBox::points().

Here is the call graph for this function:

◆ boundBox() [5/7]

boundBox ( const UList< point > &  points,
const labelUList indices,
const bool  doReduce = true 
)

Construct bounding box as subset of the pointField.

The indices could be from cell/face etc. Does parallel communication (doReduce = true)

Definition at line 108 of file boundBox.C.

References UList< T >::empty(), Foam::max(), Foam::min(), UPstream::parRun(), and Foam::reduce().

Here is the call graph for this function:

◆ boundBox() [6/7]

boundBox ( const UList< point > &  points,
const FixedList< label, Size > &  indices,
const bool  doReduce = true 
)

Construct bounding box as subset of the pointField.

The indices could be from edge/triFace etc. Does parallel communication (doReduce = true)

Definition at line 35 of file boundBoxTemplates.C.

References boundBox::contains(), UList< T >::empty(), Foam::max(), Foam::min(), and Foam::reduce().

Here is the call graph for this function:

◆ boundBox() [7/7]

boundBox ( Istream is)
inline

Construct from Istream.

Definition at line 46 of file boundBoxI.H.

References boundBox::operator>>.

Member Function Documentation

◆ min() [1/2]

◆ max() [1/2]

◆ min() [2/2]

Foam::point & min ( )
inline

Minimum describing the bounding box, non-const access.

Definition at line 66 of file boundBoxI.H.

◆ max() [2/2]

Foam::point & max ( )
inline

Maximum describing the bounding box, non-const access.

Definition at line 72 of file boundBoxI.H.

◆ midpoint()

Foam::point midpoint ( ) const
inline

The midpoint of the bounding box.

Definition at line 78 of file boundBoxI.H.

Referenced by dynamicIndexedOctree< Foam::dynamicTreeDataPoint >::overlaps(), indexedOctree< Foam::treeDataFace >::overlaps(), treeBoundBox::searchOrder(), treeBoundBox::subBbox(), and treeBoundBox::subOctant().

Here is the caller graph for this function:

◆ span()

◆ mag()

Foam::scalar mag ( ) const
inline

The magnitude of the bounding box span.

Definition at line 90 of file boundBoxI.H.

References Foam::mag().

Referenced by searchableSurfaces::checkSizes(), triSurfaceMesh::extractCloseness(), triSurfaceMesh::extractPointCloseness(), patchProbes::findElements(), globalMeshData::geometricSharedPoints(), boundBox::inflate(), streamlinesParticle::move(), and globalMeshData::updateMesh().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ volume()

Foam::scalar volume ( ) const
inline

The volume of the bound box.

Definition at line 96 of file boundBoxI.H.

References Foam::cmptProduct(), and boundBox::span().

Here is the call graph for this function:

◆ minDim()

Foam::scalar minDim ( ) const
inline

Smallest length/height/width dimension.

Definition at line 102 of file boundBoxI.H.

References Foam::cmptMin(), and boundBox::span().

Here is the call graph for this function:

◆ maxDim()

Foam::scalar maxDim ( ) const
inline

Largest length/height/width dimension.

Definition at line 108 of file boundBoxI.H.

References Foam::cmptMax(), and boundBox::span().

Here is the call graph for this function:

◆ avgDim()

Foam::scalar avgDim ( ) const
inline

Average length/height/width dimension.

Definition at line 114 of file boundBoxI.H.

References Foam::cmptAv(), and boundBox::span().

Referenced by repatchMesh::getNearest(), and treeBoundBox::typDim().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ points()

Foam::tmp< Foam::pointField > points ( ) const

Return corner points in an order corresponding to a 'hex' cell.

Definition at line 149 of file boundBox.C.

References tmp< T >::ref(), Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().

Referenced by boundBox::boundBox().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ faces()

Foam::faceList faces ( )
static

Return faces with correct point order.

Definition at line 167 of file boundBox.C.

References forAll, and List< T >::setSize().

Referenced by searchableBox::boundingSpheres(), and searchableBox::coordinates().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ inflate()

void inflate ( const scalar  s)

Inflate box by factor*mag(span) in all dimensions.

Definition at line 210 of file boundBox.C.

References boundBox::mag(), and VectorSpace< Vector< scalar >, scalar, 3 >::one.

Referenced by meshToMesh::mapAndOpTgtToSrc(), meshToMeshMethod::maskCells(), surfaceFeatures::nearestSurfEdge(), and AMIMethod::resetTree().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ overlaps() [1/2]

◆ overlaps() [2/2]

bool overlaps ( const point centre,
const scalar  radiusSqr 
) const
inline

Overlaps boundingSphere (centre + sqr(radius))?

Definition at line 132 of file boundBoxI.H.

References Foam::mag(), and VectorSpace< Vector< scalar >, scalar, 3 >::nComponents.

Here is the call graph for this function:

◆ contains() [1/5]

bool contains ( const point pt) const
inline

Contains point? (inside or on edge)

Definition at line 170 of file boundBoxI.H.

References Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().

Referenced by boundBox::boundBox(), boundBox::contains(), boundBox::containsAny(), primitiveMesh::pointInCellBB(), and Foam::selectBox().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ contains() [2/5]

bool contains ( const boundBox bb) const
inline

Fully contains other boundingBox?

Definition at line 182 of file boundBoxI.H.

References boundBox::contains(), boundBox::max(), and boundBox::min().

Here is the call graph for this function:

◆ containsInside()

bool containsInside ( const point pt) const
inline

Contains point? (inside only)

Definition at line 188 of file boundBoxI.H.

References Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().

Here is the call graph for this function:

◆ contains() [3/5]

bool contains ( const UList< point > &  points) const

Contains all of the points? (inside or on edge)

Definition at line 219 of file boundBox.C.

References boundBox::contains(), UList< T >::empty(), and forAll.

Here is the call graph for this function:

◆ contains() [4/5]

bool contains ( const UList< point > &  points,
const labelUList indices 
) const

Contains all of the points? (inside or on edge)

Definition at line 239 of file boundBox.C.

References boundBox::contains(), UList< T >::empty(), and forAll.

Here is the call graph for this function:

◆ contains() [5/5]

bool contains ( const UList< point > &  points,
const FixedList< label, Size > &  indices 
) const

Contains all of the points? (inside or on edge)

Definition at line 79 of file boundBoxTemplates.C.

References boundBox::containsAny(), UList< T >::empty(), and forAll.

Here is the call graph for this function:

◆ containsAny() [1/3]

bool containsAny ( const UList< point > &  points) const

Contains any of the points? (inside or on edge)

Definition at line 261 of file boundBox.C.

References boundBox::contains(), UList< T >::empty(), and forAll.

Referenced by boundBox::contains(), treeDataFace::overlaps(), treeDataPrimitivePatch< PatchType >::overlaps(), and triSurfaceTools::triangulate().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ containsAny() [2/3]

bool containsAny ( const UList< point > &  points,
const labelUList indices 
) const

Contains any of the points? (inside or on edge)

Definition at line 281 of file boundBox.C.

References boundBox::contains(), UList< T >::empty(), and forAll.

Here is the call graph for this function:

◆ containsAny() [3/3]

bool containsAny ( const UList< point > &  points,
const FixedList< label, Size > &  indices 
) const

Contains any of the points? (inside or on edge)

Definition at line 104 of file boundBoxTemplates.C.

References UList< T >::empty(), and forAll.

Here is the call graph for this function:

◆ nearest()

Foam::point nearest ( const point pt) const

Return the nearest point on the boundBox to the supplied point.

If point is inside the boundBox then the point is returned unchanged.

Definition at line 303 of file boundBox.C.

References Foam::max(), Foam::min(), Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().

Here is the call graph for this function:

Friends And Related Function Documentation

◆ operator==

bool operator== ( const boundBox ,
const boundBox  
)
friend

◆ operator!=

bool operator!= ( const boundBox ,
const boundBox  
)
friend

◆ operator>>

Istream& operator>> ( Istream ,
boundBox  
)
friend

Referenced by boundBox::boundBox().

◆ operator<<

Ostream& operator<< ( Ostream ,
const boundBox  
)
friend

Member Data Documentation

◆ great

const Foam::scalar great
static

The great value used for greatBox and invertedBox.

Definition at line 76 of file boundBox.H.

◆ greatBox

const Foam::boundBox greatBox
static

A very large boundBox: min/max == -/+ vGreat.

Definition at line 79 of file boundBox.H.

◆ invertedBox

const Foam::boundBox invertedBox
static

A very large inverted boundBox: min/max == +/- vGreat.

Definition at line 82 of file boundBox.H.

Referenced by NamedEnum< compressibleField, 8 >::names(), and triSurface::writeStats().


The documentation for this class was generated from the following files: