boundBox Class Reference

A bounding box defined in terms of min/max extrema points. More...

Inheritance diagram for boundBox:
[legend]
Collaboration diagram for boundBox:
[legend]

Public Types

enum  directionBit : direction { XDIR = 0x1, YDIR = 0x2, ZDIR = 0x4 }
 Bits used for (x/y/z) direction encoding. More...
 

Public Member Functions

 boundBox ()
 Construct without any points - an inverted bounding box. More...
 
 boundBox (const point &pt)
 Construct a bounding box containing a single initial point. More...
 
 boundBox (const point &min, const point &max)
 Construct from components. More...
 
 boundBox (const UList< point > &points, bool doReduce=true)
 Construct as the bounding box of the given points. More...
 
 boundBox (const tmp< pointField > &tpoints, bool doReduce=true)
 Construct as the bounding box of the given temporary pointField. More...
 
 boundBox (const UList< point > &points, const labelUList &indices, bool doReduce=true)
 Construct bounding box as an indirect subset of the points. More...
 
template<unsigned N>
 boundBox (const UList< point > &points, const FixedList< label, N > &indices, bool doReduce=true)
 Construct bounding box as an indirect subset of the points. More...
 
 boundBox (Istream &is)
 Construct from Istream. More...
 
bool empty () const
 Bounding box is inverted, contains no points. More...
 
bool valid () const
 Bounding box is non-inverted. 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 centre () const
 The centre (midpoint) of the bounding box. More...
 
point midpoint () const
 The midpoint (centre) of the bounding box. Identical to centre() 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...
 
label nDim () const
 Count the number of positive, non-zero dimensions. More...
 
tmp< pointFieldpoints () const
 Corner points in an order corresponding to a 'hex' cell. More...
 
tmp< pointFieldfaceCentres () const
 Face midpoints. More...
 
point faceCentre (const direction facei) const
 Face centre of given face index. More...
 
void clear ()
 Clear bounding box and make it an inverted box. More...
 
void add (const boundBox &bb)
 Extend to include the second box. More...
 
void add (const point &pt)
 Extend to include the point. More...
 
void add (const UList< point > &points)
 Extend to include the points. More...
 
void add (const tmp< pointField > &tpoints)
 Extend to include the points from the temporary point field. More...
 
template<unsigned N>
void add (const FixedList< point, N > &points)
 Extend to include the points. More...
 
template<unsigned N>
void add (const UList< point > &points, const FixedList< label, N > &indices)
 Extend to include a (subsetted) point field. More...
 
template<class IntContainer >
void add (const UList< point > &points, const IntContainer &indices)
 Extend to include a (subsetted) point field. More...
 
void inflate (const scalar s)
 Inflate box by factor*mag(span) in all dimensions. More...
 
void reduce ()
 Parallel reduction of min/max values. More...
 
bool intersect (const boundBox &bb)
 Intersection (union) with the second box. More...
 
bool intersects (const plane &pln) const
 Does plane intersect this bounding box. More...
 
bool overlaps (const boundBox &bb) const
 Overlaps/touches boundingBox? More...
 
bool overlaps (const point &centre, const scalar radiusSqr) const
 Overlaps boundingSphere (centre + sqr(radius))? More...
 
bool contains (const point &pt) const
 Contains point? (inside or on edge) More...
 
bool contains (const boundBox &bb) const
 Fully contains other boundingBox? More...
 
bool containsInside (const point &pt) const
 Contains point? (inside only) More...
 
bool contains (const UList< point > &points) const
 Contains all points? (inside or on edge) More...
 
template<unsigned N>
bool contains (const UList< point > &points, const FixedList< label, N > &indices) const
 Contains all of the (subsetted) points? (inside or on edge) More...
 
template<class IntContainer >
bool contains (const UList< point > &points, const IntContainer &indices) const
 Contains all of the (subsetted) points? (inside or on edge) More...
 
bool containsAny (const UList< point > &points) const
 Contains any of the points? (inside or on edge) More...
 
template<unsigned N>
bool containsAny (const UList< point > &points, const FixedList< label, N > &indices) const
 Contains any of the (subsetted) points? (inside or on edge) More...
 
template<class IntContainer >
bool containsAny (const UList< point > &points, const IntContainer &indices) const
 Contains any of the (subsetted) points? (inside or on edge) More...
 
point nearest (const point &pt) const
 Return the nearest point on the boundBox to the supplied point. More...
 
void operator+= (const boundBox &bb)
 Extend box to include the second box, as per the add() method. More...
 

Static Public Attributes

static const boundBox greatBox
 A large boundBox: min/max == -/+ ROOTVGREAT. More...
 
static const boundBox invertedBox
 A large inverted boundBox: min/max == +/- ROOTVGREAT. More...
 
static const faceList faces
 Faces to point addressing, as per a 'hex' cell. More...
 
static const FixedList< vector, 6 > faceNormals
 The unit normal per face. More...
 

Friends

Istreamoperator>> (Istream &is, boundBox &bb)
 
Ostreamoperator<< (Ostream &os, const boundBox &bb)
 

Detailed Description

A bounding box defined in terms of min/max extrema points.

Note
When a bounding box is created without any points, it creates an inverted bounding box. Points can be added later and the bounding box will grow to include them.

Definition at line 63 of file boundBox.H.

Member Enumeration Documentation

◆ directionBit

Bits used for (x/y/z) direction encoding.

Enumerator
XDIR 

1: x-direction (vector component 0)

YDIR 

2: y-direction (vector component 1)

ZDIR 

4: z-direction (vector component 2)

Definition at line 75 of file boundBox.H.

Constructor & Destructor Documentation

◆ boundBox() [1/8]

boundBox ( )
inline

Construct without any points - an inverted bounding box.

Definition at line 33 of file boundBoxI.H.

◆ boundBox() [2/8]

boundBox ( const point pt)
inlineexplicit

Construct a bounding box containing a single initial point.

Definition at line 40 of file boundBoxI.H.

◆ boundBox() [3/8]

boundBox ( const point min,
const point max 
)
inline

Construct from components.

Definition at line 47 of file boundBoxI.H.

◆ boundBox() [4/8]

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

Construct as the bounding box of the given points.

Does parallel communication (doReduce = true)

Definition at line 72 of file boundBox.C.

References boundBox::add(), boundBox::points(), and boundBox::reduce().

Here is the call graph for this function:

◆ boundBox() [5/8]

boundBox ( const tmp< pointField > &  tpoints,
bool  doReduce = true 
)
explicit

Construct as the bounding box of the given temporary pointField.

Does parallel communication (doReduce = true)

Definition at line 85 of file boundBox.C.

References boundBox::add(), and boundBox::reduce().

Here is the call graph for this function:

◆ boundBox() [6/8]

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

Construct bounding box as an indirect subset of the points.

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

Definition at line 99 of file boundBox.C.

References Foam::add(), points, and Foam::reduce().

Here is the call graph for this function:

◆ boundBox() [7/8]

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

Construct bounding box as an indirect subset of the points.

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

Definition at line 36 of file boundBoxTemplates.C.

References Foam::add(), points, and reduce().

Here is the call graph for this function:

◆ boundBox() [8/8]

boundBox ( Istream is)
inline

Construct from Istream.

Definition at line 54 of file boundBoxI.H.

References Foam::operator>>().

Here is the call graph for this function:

Member Function Documentation

◆ empty()

bool empty ( ) const
inline

Bounding box is inverted, contains no points.

Definition at line 62 of file boundBoxI.H.

References VectorSpace< Vector< scalar >, scalar, 3 >::nComponents.

Referenced by fieldExtents::calcFieldExtents(), and cuttingSurfaceBase::cellSelection().

Here is the caller graph for this function:

◆ valid()

bool valid ( ) const
inline

Bounding box is non-inverted.

Definition at line 76 of file boundBoxI.H.

References VectorSpace< Vector< scalar >, scalar, 3 >::nComponents.

Referenced by isoSurfaceBase::blockCells(), cuttingPlane::checkOverlap(), cuttingSurfaceBase::checkOverlap(), searchableSphere::overlaps(), and searchableBox::searchableBox().

Here is the caller graph for this function:

◆ min() [1/2]

const Foam::point & min ( ) const
inline

Minimum describing the bounding box.

Definition at line 91 of file boundBoxI.H.

Referenced by PDRblock::blockMeshDict(), eddy::bounds(), fieldExtents::calcFieldExtents(), faceAreaWeightAMI2D::calculate(), inverseDistance::cellBb(), polyMesh::cellTree(), voxelMeshSearch::centre(), boundBox::contains(), AABBTree< Type >::createBoxes(), triSurfaceMesh::edgeTree(), extendedEdgeMesh::edgeTree(), extendedEdgeMesh::edgeTreesByType(), treeBoundBox::extend(), inverseDistance::fill(), voxelMeshSearch::fill(), patchProbes::findElements(), mappedPatchBase::findLocalSamples(), dynamicTreeDataPoint::findNearest(), dynamicIndexedOctree< Foam::dynamicTreeDataPoint >::findNearest(), indexedOctree< Foam::treeDataPrimitivePatch< PatchType > >::findNearest(), meshSearch::findNearestBoundaryFace(), boundaryMesh::getNearest(), voxelMeshSearch::index(), inverseDistance::index3(), voxelMeshSearch::index3(), waveModel::initialiseGeometry(), waveMakerPointPatchVectorField::initialiseGeometry(), trackingInverseDistance::markBoundaries(), inverseDistance::markBoundaries(), trackingInverseDistance::markDonors(), inverseDistance::markDonors(), trackingInverseDistance::markPatchesAsHoles(), inverseDistance::markPatchesAsHoles(), meshToMeshMethod::maskCells(), projectVertex::operator point(), Foam::operator==(), inverseDistance::overlaps(), voxelMeshSearch::overlaps(), searchableSphere::overlaps(), extendedEdgeMesh::pointTree(), inverseDistance::position(), energySpectrum::read(), Foam::readBoxDim(), refinementFeatures::regionEdgeTrees(), searchableExtrudedCircle::searchableExtrudedCircle(), searchableSurfaceCollection::searchableSurfaceCollection(), treeBoundBox::subBbox(), triSurfaceSearch::tree(), triSurfaceRegionSearch::treeByRegion(), topoSet::writeDebug(), and voxelMeshSearch::writeGrid().

◆ max() [1/2]

const Foam::point & max ( ) const
inline

Maximum describing the bounding box.

Definition at line 97 of file boundBoxI.H.

Referenced by PDRblock::blockMeshDict(), eddy::bounds(), fieldExtents::calcFieldExtents(), faceAreaWeightAMI2D::calculate(), inverseDistance::cellBb(), polyMesh::cellTree(), boundBox::contains(), triSurfaceMesh::edgeTree(), extendedEdgeMesh::edgeTree(), extendedEdgeMesh::edgeTreesByType(), treeBoundBox::extend(), inverseDistance::fill(), voxelMeshSearch::fill(), patchProbes::findElements(), dynamicTreeDataPoint::findNearest(), dynamicIndexedOctree< Foam::dynamicTreeDataPoint >::findNearest(), indexedOctree< Foam::treeDataPrimitivePatch< PatchType > >::findNearest(), meshSearch::findNearestBoundaryFace(), advancingFrontAMI::findTargetFace(), boundaryMesh::getNearest(), waveModel::initialiseGeometry(), waveMakerPointPatchVectorField::initialiseGeometry(), Kmesh::Kmesh(), trackingInverseDistance::markDonors(), inverseDistance::markDonors(), trackingInverseDistance::markPatchesAsHoles(), inverseDistance::markPatchesAsHoles(), meshToMeshMethod::maskCells(), projectVertex::operator point(), Foam::operator==(), inverseDistance::overlaps(), voxelMeshSearch::overlaps(), searchableSphere::overlaps(), extendedEdgeMesh::pointTree(), energySpectrum::read(), Foam::readBoxDim(), refinementFeatures::regionEdgeTrees(), searchableExtrudedCircle::searchableExtrudedCircle(), searchableSurfaceCollection::searchableSurfaceCollection(), treeBoundBox::subBbox(), triSurfaceSearch::tree(), triSurfaceRegionSearch::treeByRegion(), and topoSet::writeDebug().

Here is the caller graph for this function:

◆ min() [2/2]

Foam::point & min ( )
inline

Minimum describing the bounding box, non-const access.

Definition at line 103 of file boundBoxI.H.

◆ max() [2/2]

Foam::point & max ( )
inline

Maximum describing the bounding box, non-const access.

Definition at line 109 of file boundBoxI.H.

◆ centre()

Foam::point centre ( ) const
inline

The centre (midpoint) of the bounding box.

Definition at line 115 of file boundBoxI.H.

Referenced by PDRblock::blockMeshDict(), boundBox::faceCentre(), and advancingFrontAMI::findTargetFace().

Here is the caller graph for this function:

◆ midpoint()

Foam::point midpoint ( ) const
inline

The midpoint (centre) of the bounding box. Identical to centre()

Definition at line 121 of file boundBoxI.H.

◆ span()

Foam::vector span ( ) const
inline

The bounding box span (from minimum to maximum)

Definition at line 127 of file boundBoxI.H.

Referenced by PDRblock::blockMeshDict(), voxelMeshSearch::centre(), sensitivitySurface::computeRadius(), AABBTree< Type >::createBoxes(), treeBoundBox::extend(), voxelMeshSearch::index(), inverseDistance::index3(), voxelMeshSearch::index3(), Kmesh::Kmesh(), projectVertex::operator point(), porosityModel::porosityModel(), inverseDistance::position(), points0MotionSolver::updateMesh(), and voxelMeshSearch::writeGrid().

Here is the caller graph for this function:

◆ mag()

Foam::scalar mag ( ) const
inline

The magnitude of the bounding box span.

Definition at line 133 of file boundBoxI.H.

References Foam::mag().

Referenced by nearWallFields::calcAddressing(), refinementFeatures::checkSizes(), searchableSurfaces::checkSizes(), isoSurfaceCell::isoSurfaceCell(), isoSurfaceTopo::isoSurfaceTopo(), streamLineParticle::move(), and triSurfaceTools::writeCloseness().

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 139 of file boundBoxI.H.

References Foam::cmptProduct().

Here is the call graph for this function:

◆ minDim()

Foam::scalar minDim ( ) const
inline

Smallest length/height/width dimension.

Definition at line 145 of file boundBoxI.H.

References Foam::cmptMin().

Here is the call graph for this function:

◆ maxDim()

Foam::scalar maxDim ( ) const
inline

Largest length/height/width dimension.

Definition at line 151 of file boundBoxI.H.

References Foam::cmptMax().

Here is the call graph for this function:

◆ avgDim()

Foam::scalar avgDim ( ) const
inline

Average length/height/width dimension.

Definition at line 157 of file boundBoxI.H.

References Foam::cmptAv().

Referenced by boundaryMesh::getNearest().

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

◆ nDim()

Foam::label nDim ( ) const
inline

Count the number of positive, non-zero dimensions.

Returns
-1 if any dimensions are negative, 0 = 0D (point), 1 = 1D (line aligned with an axis), 2 = 2D (plane aligned with an axis), 3 = 3D (box)

Definition at line 163 of file boundBoxI.H.

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

Here is the call graph for this function:

◆ points()

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

Corner points in an order corresponding to a 'hex' cell.

Definition at line 118 of file boundBox.C.

References tmp< T >::New().

Referenced by PDRblock::blockMeshDict(), boundBox::boundBox(), searchableRotatedBox::overlaps(), and faceAreaWeightAMI2D::writeNoMatch().

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

◆ faceCentres()

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

Face midpoints.

Definition at line 136 of file boundBox.C.

References forAll, and tmp< T >::New().

Here is the call graph for this function:

◆ faceCentre()

Foam::point faceCentre ( const direction  facei) const

Face centre of given face index.

Definition at line 150 of file boundBox.C.

References Foam::abort(), boundBox::centre(), Foam::FatalError, FatalErrorInFunction, Foam::max(), Foam::min(), Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().

Here is the call graph for this function:

◆ clear()

void clear ( )
inline

Clear bounding box and make it an inverted box.

Definition at line 184 of file boundBoxI.H.

Referenced by cuttingSurfaceBase::cellSelection(), and Foam::simpleGeometricFilter().

Here is the caller graph for this function:

◆ add() [1/7]

void add ( const boundBox bb)
inline

Extend to include the second box.

Definition at line 191 of file boundBoxI.H.

References Foam::max(), and Foam::min().

Referenced by boundBox::boundBox(), triPoints::bounds(), tetPoints::bounds(), searchableSurfacesQueries::bounds(), PatchTools::calcBounds(), fieldExtents::calcFieldExtents(), cuttingSurfaceBase::cellSelection(), patchProbes::findElements(), meshToMeshMethod::intersect(), meshToMeshMethod::interVol(), meshToMeshMethod::interVolAndCentroid(), porosityModel::porosityModel(), energySpectrum::read(), Foam::simpleGeometricFilter(), and triSurface::writeStats().

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

◆ add() [2/7]

void add ( const point pt)
inline

Extend to include the point.

Definition at line 198 of file boundBoxI.H.

References Foam::max(), and Foam::min().

Here is the call graph for this function:

◆ add() [3/7]

void add ( const UList< point > &  points)
inline

Extend to include the points.

Definition at line 205 of file boundBoxI.H.

References Foam::add(), p, and points.

Here is the call graph for this function:

◆ add() [4/7]

void add ( const tmp< pointField > &  tpoints)
inline

Extend to include the points from the temporary point field.

Definition at line 214 of file boundBoxI.H.

References Foam::add(), and tmp< T >::clear().

Here is the call graph for this function:

◆ add() [5/7]

void add ( const FixedList< point, N > &  points)

Extend to include the points.

Definition at line 57 of file boundBoxTemplates.C.

References Foam::add(), p, and points.

Here is the call graph for this function:

◆ add() [6/7]

void add ( const UList< point > &  points,
const FixedList< label, N > &  indices 
)

Extend to include a (subsetted) point field.

The indices could be from edge/triFace etc.

Definition at line 70 of file boundBoxTemplates.C.

References Foam::add(), and points.

Here is the call graph for this function:

◆ add() [7/7]

void add ( const UList< point > &  points,
const IntContainer &  indices 
)

Extend to include a (subsetted) point field.

Template Parameters
IntContainerA container with an iterator that dereferences to an label

Definition at line 95 of file boundBoxTemplates.C.

References Foam::add(), and points.

Here is the call graph for this function:

◆ inflate()

void inflate ( const scalar  s)

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

Definition at line 175 of file boundBox.C.

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

Referenced by AABBTree< Type >::AABBTree(), box::box(), AMIInterpolation::createTree(), meshToMeshMethod::maskCells(), surfaceFeatures::nearestSurfEdge(), primitiveMesh::pointInCellBB(), Foam::simpleGeometricFilter(), and propellerInfo::updateSampleDiskCells().

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

◆ reduce()

void reduce ( )

Parallel reduction of min/max values.

Definition at line 184 of file boundBox.C.

References Foam::reduce().

Referenced by boundBox::boundBox(), fieldExtents::calcFieldExtents(), cuttingSurfaceBase::cellSelection(), distributedTriSurfaceMesh::distributedTriSurfaceMesh(), porosityModel::porosityModel(), and distributedTriSurfaceMesh::writeStats().

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

◆ intersect()

bool intersect ( const boundBox bb)

Intersection (union) with the second box.

The return value is true if the intersection is non-empty.

Definition at line 191 of file boundBox.C.

References Foam::max(), and Foam::min().

Referenced by sampledMeshedSurface::update().

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

◆ intersects()

bool intersects ( const plane pln) const

Does plane intersect this bounding box.

There is an intersection if the plane segments the corner points

Note
the results are unreliable when plane coincides almost exactly with a box face

Definition at line 200 of file boundBox.C.

References plane::FRONT, p, points, and plane::sideOfPlane().

Referenced by cuttingPlane::checkOverlap().

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

◆ overlaps() [1/2]

bool overlaps ( const boundBox bb) const
inline

Overlaps/touches boundingBox?

Definition at line 221 of file boundBoxI.H.

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

Referenced by cuttingSurfaceBase::checkOverlap(), trackingInverseDistance::markBoundaries(), inverseDistance::markBoundaries(), meshToMeshMethod::maskCells(), faceAreaWeightAMI2D::overlappingTgtFaces(), searchableBox::overlaps(), searchableRotatedBox::overlaps(), searchablePlate::overlaps(), searchableSphere::overlaps(), indexedOctree< Foam::treeDataPrimitivePatch< PatchType > >::overlaps(), and interRegionOption::setMapper().

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

◆ overlaps() [2/2]

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

Overlaps boundingSphere (centre + sqr(radius))?

Definition at line 233 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 271 of file boundBoxI.H.

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

Referenced by isoSurfaceBase::blockCells(), cuttingSurfaceBase::cellSelection(), trackingInverseDistance::markBoundaries(), primitiveMesh::pointInCellBB(), and Foam::simpleGeometricFilter().

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 282 of file boundBoxI.H.

References 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 288 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 points? (inside or on edge)

Definition at line 230 of file boundBox.C.

References p, and points.

◆ contains() [4/5]

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

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

Definition at line 120 of file boundBoxTemplates.C.

References points.

◆ contains() [5/5]

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

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

Template Parameters
IntContainerA container with an iterator that dereferences to an label

Definition at line 149 of file boundBoxTemplates.C.

References points.

◆ containsAny() [1/3]

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

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

Definition at line 249 of file boundBox.C.

References p, and points.

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

Here is the caller graph for this function:

◆ containsAny() [2/3]

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

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

Definition at line 178 of file boundBoxTemplates.C.

References points.

◆ containsAny() [3/3]

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

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

Template Parameters
IntContainerA container with an iterator that dereferences to an label

Definition at line 211 of file boundBoxTemplates.C.

References points.

◆ 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 268 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:

◆ operator+=()

void operator+= ( const boundBox bb)
inline

Extend box to include the second box, as per the add() method.

Definition at line 301 of file boundBoxI.H.

References Foam::add().

Here is the call graph for this function:

Friends And Related Function Documentation

◆ operator>>

Istream& operator>> ( Istream is,
boundBox bb 
)
friend

◆ operator<<

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

Member Data Documentation

◆ greatBox

const Foam::boundBox greatBox
static

A large boundBox: min/max == -/+ ROOTVGREAT.

Definition at line 83 of file boundBox.H.

◆ invertedBox

◆ faces

const Foam::faceList faces
static

Faces to point addressing, as per a 'hex' cell.

Definition at line 89 of file boundBox.H.

◆ faceNormals

const Foam::FixedList< Foam::vector, 6 > faceNormals
static

The unit normal per face.

Definition at line 92 of file boundBox.H.

Referenced by searchableBox::getNormal().


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