Loading...
Searching...
No Matches
boundBox Class Reference

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

#include <boundBox.H>

Inheritance diagram for boundBox:
Collaboration diagram for boundBox:

Public Types

enum  directionBit : direction { XDIR = 1 , YDIR = 2 , ZDIR = 4 }
 Bits used for (x/y/z) directional encoding. More...

Public Member Functions

 boundBox ()
 Default construct: an inverted bounding box.
 boundBox (const boundBox &)=default
 Copy construct.
boundBoxoperator= (const boundBox &)=default
 Copy assignment.
 boundBox (const boundBox &bb, bool doReduce)
 Copy construct with specified global reduction.
 boundBox (const point &p)
 Construct a bounding box containing a single initial point.
 boundBox (Foam::zero_one)
 Construct a 0/1 unit bounding box.
 boundBox (const point &min, const point &max)
 Construct from bound box min/max points.
 boundBox (const Pair< point > &bb)
 Construct from bound box min/max points.
 boundBox (const UList< point > &points, bool doReduce=true)
 Construct as the bounding box of the given points.
 boundBox (const tmp< pointField > &tpoints, bool doReduce=true)
 Construct as the bounding box of the given temporary pointField.
 boundBox (const UList< point > &points, const labelUList &indices, bool doReduce=true)
 Construct bounding box as an indirect subset of the points.
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.
 boundBox (Istream &is)
 Construct from Istream.
bool empty () const
 Bounding box is inverted, contains no points.
bool good () const
 Bounding box is non-inverted.
bool valid () const
 Bounding box is non-inverted - same as good().
const pointmin () const noexcept
 Minimum describing the bounding box.
const pointmax () const noexcept
 Maximum describing the bounding box.
pointmin () noexcept
 Minimum describing the bounding box, non-const access.
pointmax () noexcept
 Maximum describing the bounding box, non-const access.
point centre () const
 The centre (midpoint) of the bounding box.
vector span () const
 The bounding box span (from minimum to maximum).
scalar mag () const
 The magnitude/length of the bounding box diagonal.
scalar magSqr () const
 The magnitude/length squared of bounding box diagonal.
scalar volume () const
 The volume of the bound box.
scalar minDim () const
 Smallest length/height/width dimension.
scalar maxDim () const
 Largest length/height/width dimension.
scalar avgDim () const
 Average length/height/width dimension.
direction minDir () const
 Direction (X/Y/Z) of the smallest span (for empty box: 0).
direction maxDir () const
 Direction (X/Y/Z) of the largest span (for empty box: 0).
int nDim () const
 Count the number of positive, non-zero dimensions.
template<direction CornerNumber>
point hexCorner () const
 Return corner point [0..7] corresponding to a 'hex' cell.
tmp< pointFieldhexCorners () const
 Corner points in an order corresponding to a 'hex' cell.
tmp< pointFieldpoints () const
 Corner points in an order corresponding to a 'hex' cell.
tmp< pointFieldfaceCentres () const
 Face midpoints.
point faceCentre (const direction facei) const
 Face centre of given face index.
void reset ()
 Reset to an inverted box.
void reset (Foam::zero_one)
 Reset to a 0/1 unit bounding box.
void reset (const point &pt)
 Reset min/max to be identical to the specified point.
void reset (const point &min, const point &max)
 Reset min/max to specified values.
void clear ()
 Same as reset() - reset to an inverted box.
void add (const boundBox &bb)
 Extend to include the second box.
void add (const point &pt)
 Extend to include the point.
void add (const point &pt0, const point &pt1)
 Extend to include two additional points.
void add (const Pair< point > &points)
 Extend to include two additional points.
void add (const UList< point > &points)
 Extend to include the points.
void add (const tmp< pointField > &tpoints)
 Extend to include the points from the temporary point field.
template<unsigned N>
void add (const FixedList< point, N > &points)
 Extend to include the points.
template<unsigned N>
void add (const UList< point > &points, const FixedList< label, N > &indices)
 Extend to include a (subsetted) point field.
template<class IntContainer>
void add (const UList< point > &points, const IntContainer &indices)
 Extend to include a (subsetted) point field.
void grow (const scalar delta)
 Expand box by adjusting min/max by specified amount in each dimension.
void grow (const vector &delta)
 Expand box by adjusting min/max by specified amounts.
void inflate (const scalar factor)
 Expand box by factor*mag(span) in all dimensions.
void inflate (Random &rndGen, const scalar factor)
 Expand box slightly by expanding all dimensions with factor*span*(random 0-1) and guarantees factor*mag(span) minimum width in any direction.
void inflate (Random &r, const scalar factor, const scalar delta)
 As per two parameter version but with additional grow() by given amount in each dimension.
void reduce ()
 Inplace parallel reduction of min/max values, using UPstream::worldComm.
void reduce (int communicator)
 Inplace parallel reduction of min/max values, using the specified communicator.
bool intersects (const plane &pln) const
 Does plane intersect this bounding box.
bool intersects (const triPointRef &tri) const
 Does triangle intersect this bounding box or is contained within this bounding box.
bool overlaps (const boundBox &bb) const
 Overlaps/touches boundingBox?
bool overlaps (const point &centre, const scalar radiusSqr) const
 Overlaps boundingSphere (centre + sqr(radius))?
bool contains (const point &pt) const
 Contains point? (inside or on edge).
bool contains (const boundBox &bb) const
 Fully contains other boundingBox?
bool containsInside (const point &pt) const
 Contains point? (inside only).
bool contains (const UList< point > &points) const
 Contains all points? (inside or on edge).
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).
template<class IntContainer>
bool contains (const UList< point > &points, const IntContainer &indices) const
 Contains all of the (subsetted) points? (inside or on edge).
bool containsAny (const UList< point > &points) const
 Contains any of the points? (inside or on edge).
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).
template<class IntContainer>
bool containsAny (const UList< point > &points, const IntContainer &indices) const
 Contains any of the (subsetted) points? (inside or on edge).
point nearest (const point &p) const
 Return the nearest point on the boundBox to the supplied point.
void operator&= (const boundBox &bb)
 Restrict min/max to union with other box.
void operator+= (const boundBox &bb)
 Extend box to include the second box, as per the add() method.
bool intersect (const boundBox &bb)
 Deprecated(2022-10) - use 'operator&=' to avoid confusion with other intersects() methods.
point midpoint () const
 Identical to centre().
template<>
Foam::point hexCorner () const
template<>
Foam::point hexCorner () const
template<>
Foam::point hexCorner () const
template<>
Foam::point hexCorner () const
template<>
Foam::point hexCorner () const
template<>
Foam::point hexCorner () const
template<>
Foam::point hexCorner () const
template<>
Foam::point hexCorner () const
template<Foam::direction CornerNumber>
Foam::point hexCorner () const

Static Public Member Functions

static const boundBoxnull () noexcept
 The null boundBox is the same as an inverted box.
static const Foam::faceListhexFaces ()
 The boundBox faces as a hexCell, using hexCorner points. Same as hexCell::modelFaces().
static constexpr label nPoints () noexcept
 Number of points for boundBox and HEX.
static constexpr label nEdges () noexcept
 Number of edges for boundBox and HEX.
static constexpr label nFaces () noexcept
 Number of faces for boundBox and HEX.
static boundBox returnReduce (const boundBox &bb)
 Perform a reduction on a copy and return the result, using UPstream::worldComm.
static boundBox returnReduce (const boundBox &bb, int communicator)
 Perform a reduction on a copy and return the result, using the specified communicator.

Static Public Attributes

static const boundBox greatBox
 A large boundBox: min/max == -/+ ROOTVGREAT.
static const boundBox invertedBox
 A large inverted boundBox: min/max == +/- ROOTVGREAT.
static const FixedList< vector, 6 > faceNormals
 The unit normal per face.

Static Protected Member Functions

static bool box_box_overlaps (const point &minA, const point &maxA, const point &minB, const point &maxB)
 Test for overlap of box and box (inclusive check).
static bool box_sphere_overlaps (const point &corner0, const point &corner1, const point &centre, const scalar radiusSqr)
 Test for overlap of box and boundingSphere (centre + sqr(radius)).

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.
See also
Foam::treeBoundBox
Source files

Definition at line 70 of file boundBox.H.

Member Enumeration Documentation

◆ directionBit

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

Enumerator
XDIR 

1: x-direction. Same as (1 << vector::X)

YDIR 

2: y-direction. Same as (1 << vector::Y)

ZDIR 

4: z-direction. Same as (1 << vector::Z)

Definition at line 116 of file boundBox.H.

Constructor & Destructor Documentation

◆ boundBox() [1/12]

◆ boundBox() [2/12]

boundBox ( const boundBox & )
default

Copy construct.

References boundBox().

Here is the call graph for this function:

◆ boundBox() [3/12]

boundBox ( const boundBox & bb,
bool doReduce )

Copy construct with specified global reduction.

References boundBox(), max(), min(), and p.

Here is the call graph for this function:

◆ boundBox() [4/12]

boundBox ( const point & p)
inlineexplicit

Construct a bounding box containing a single initial point.

Definition at line 115 of file boundBoxI.H.

References p.

◆ boundBox() [5/12]

boundBox ( Foam::zero_one )
inlineexplicit

Construct a 0/1 unit bounding box.

Definition at line 108 of file boundBoxI.H.

◆ boundBox() [6/12]

boundBox ( const point & min,
const point & max )
inline

Construct from bound box min/max points.

Definition at line 122 of file boundBoxI.H.

References max(), and min().

Here is the call graph for this function:

◆ boundBox() [7/12]

boundBox ( const Pair< point > & bb)
inlineexplicit

Construct from bound box min/max points.

TBD: Construct from bound box min/max points inline explicit boundBox(const MinMax<point>& bb);

Definition at line 129 of file boundBoxI.H.

◆ boundBox() [8/12]

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

Construct as the bounding box of the given points.

Does parallel communication (doReduce = true)

References points().

Here is the call graph for this function:

◆ boundBox() [9/12]

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

Construct as the bounding box of the given temporary pointField.

Does parallel communication (doReduce = true)

◆ boundBox() [10/12]

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)

References points().

Here is the call graph for this function:

◆ boundBox() [11/12]

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.

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

References points().

Here is the call graph for this function:

◆ boundBox() [12/12]

boundBox ( Istream & is)
inlineexplicit

Construct from Istream.

Definition at line 136 of file boundBoxI.H.

References Foam::operator>>().

Here is the call graph for this function:

Member Function Documentation

◆ box_box_overlaps()

bool box_box_overlaps ( const point & minA,
const point & maxA,
const point & minB,
const point & maxB )
inlinestaticprotected

Test for overlap of box and box (inclusive check).

Definition at line 387 of file boundBoxI.H.

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

Referenced by overlaps().

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

◆ box_sphere_overlaps()

bool box_sphere_overlaps ( const point & corner0,
const point & corner1,
const point & centre,
const scalar radiusSqr )
inlinestaticprotected

Test for overlap of box and boundingSphere (centre + sqr(radius)).

Note: ordering of corners is irrelevant

Definition at line 402 of file boundBoxI.H.

References centre(), Foam::magSqr(), Foam::min(), and VectorSpace< Form, Cmpt, Ncmpts >::nComponents.

Referenced by overlaps(), and treeBoundBox::subOverlaps().

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

◆ null()

const boundBox & null ( )
inlinestaticnoexcept

The null boundBox is the same as an inverted box.

Definition at line 144 of file boundBox.H.

References boundBox(), invertedBox, and Foam::noexcept.

Referenced by cellCentreSet::TypeName(), and sampledCuttingSurface::TypeName().

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

◆ hexFaces()

const Foam::faceList & hexFaces ( )
static

The boundBox faces as a hexCell, using hexCorner points. Same as hexCell::modelFaces().

◆ nPoints()

constexpr label nPoints ( )
inlinestaticconstexprnoexcept

Number of points for boundBox and HEX.

Definition at line 158 of file boundBox.H.

References Foam::noexcept.

◆ nEdges()

constexpr label nEdges ( )
inlinestaticconstexprnoexcept

Number of edges for boundBox and HEX.

Definition at line 163 of file boundBox.H.

References Foam::noexcept.

◆ nFaces()

constexpr label nFaces ( )
inlinestaticconstexprnoexcept

Number of faces for boundBox and HEX.

Definition at line 168 of file boundBox.H.

References Foam::noexcept.

◆ operator=()

boundBox & operator= ( const boundBox & )
default

Copy assignment.

References boundBox().

Here is the call graph for this function:

◆ empty()

bool empty ( ) const
inline

Bounding box is inverted, contains no points.

Definition at line 144 of file boundBoxI.H.

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

Here is the caller graph for this function:

◆ good()

bool good ( ) const
inline

Bounding box is non-inverted.

Definition at line 156 of file boundBoxI.H.

References empty().

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

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

◆ valid()

bool valid ( ) const
inline

Bounding box is non-inverted - same as good().

Definition at line 283 of file boundBox.H.

References good().

Here is the call graph for this function:

◆ min() [1/2]

◆ max() [1/2]

◆ min() [2/2]

Foam::point & min ( )
inlinenoexcept

Minimum describing the bounding box, non-const access.

Definition at line 174 of file boundBoxI.H.

References Foam::noexcept.

◆ max() [2/2]

Foam::point & max ( )
inlinenoexcept

Maximum describing the bounding box, non-const access.

Definition at line 180 of file boundBoxI.H.

References Foam::noexcept.

◆ centre()

◆ span()

◆ mag()

Foam::scalar mag ( ) const
inline

The magnitude/length of the bounding box diagonal.

Definition at line 198 of file boundBoxI.H.

References Vector< Cmpt >::dist().

Referenced by refinementFeatures::checkSizes(), searchableSurfaces::checkSizes(), patchFieldProbe::findElements(), streamLineParticle::move(), wallBoundedStreamLineParticle::move(), and triSurfaceTools::writeCloseness().

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

◆ magSqr()

Foam::scalar magSqr ( ) const
inline

The magnitude/length squared of bounding box diagonal.

Definition at line 204 of file boundBoxI.H.

Referenced by mappedPatchBase::findLocalSamples(), advancingFrontAMI::findTargetFace(), and projectVertex::operator point().

Here is the caller graph for this function:

◆ volume()

Foam::scalar volume ( ) const
inline

The volume of the bound box.

Definition at line 210 of file boundBoxI.H.

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

Here is the call graph for this function:

◆ minDim()

Foam::scalar minDim ( ) const
inline

Smallest length/height/width dimension.

Definition at line 216 of file boundBoxI.H.

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

Here is the call graph for this function:

◆ maxDim()

Foam::scalar maxDim ( ) const
inline

Largest length/height/width dimension.

Definition at line 222 of file boundBoxI.H.

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

Here is the call graph for this function:

◆ avgDim()

Foam::scalar avgDim ( ) const
inline

Average length/height/width dimension.

Definition at line 228 of file boundBoxI.H.

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

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

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

◆ minDir()

Foam::direction minDir ( ) const
inline

Direction (X/Y/Z) of the smallest span (for empty box: 0).

Definition at line 234 of file boundBoxI.H.

References VectorSpace< Form, Cmpt, Ncmpts >::nComponents.

◆ maxDir()

Foam::direction maxDir ( ) const
inline

Direction (X/Y/Z) of the largest span (for empty box: 0).

Definition at line 254 of file boundBoxI.H.

References VectorSpace< Form, Cmpt, Ncmpts >::nComponents.

Referenced by AABBTree< Type >::createBoxes().

Here is the caller graph for this function:

◆ nDim()

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

References VectorSpace< Form, Cmpt, Ncmpts >::nComponents.

◆ hexCorner() [1/10]

template<direction CornerNumber>
point hexCorner ( ) const
inline

Return corner point [0..7] corresponding to a 'hex' cell.

References hexCorner().

Referenced by hexCorner().

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

◆ hexCorners()

tmp< pointField > hexCorners ( ) const

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

References hexCorners().

Referenced by hexCorners(), points(), and ensightCells::writeBox().

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

◆ points()

tmp< pointField > points ( ) const
inline

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

Definition at line 381 of file boundBox.H.

References hexCorners(), and points().

Referenced by add(), add(), add(), add(), add(), PDRblock::blockMeshDict(), boundBox(), boundBox(), boundBox(), contains(), contains(), contains(), containsAny(), containsAny(), containsAny(), searchableRotatedBox::overlaps(), points(), and faceAreaWeightAMI2D::writeNoMatch().

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

◆ faceCentres()

tmp< pointField > faceCentres ( ) const

Face midpoints.

◆ faceCentre()

point faceCentre ( const direction facei) const

Face centre of given face index.

References max(), and min().

Here is the call graph for this function:

◆ reset() [1/4]

void reset ( )
inline

◆ reset() [2/4]

void reset ( Foam::zero_one )
inline

Reset to a 0/1 unit bounding box.

Definition at line 302 of file boundBoxI.H.

References VectorSpace< Form, Cmpt, Ncmpts >::one, and VectorSpace< Form, Cmpt, Ncmpts >::zero.

◆ reset() [3/4]

void reset ( const point & pt)
inline

Reset min/max to be identical to the specified point.

Definition at line 309 of file boundBoxI.H.

◆ reset() [4/4]

void reset ( const point & min,
const point & max )
inline

Reset min/max to specified values.

Definition at line 316 of file boundBoxI.H.

References max(), and min().

Here is the call graph for this function:

◆ clear()

void clear ( )
inline

Same as reset() - reset to an inverted box.

Definition at line 419 of file boundBox.H.

References reset().

Referenced by cuttingSurfaceBase::cellSelection().

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

◆ add() [1/9]

◆ add() [2/9]

void add ( const point & pt)
inline

Extend to include the point.

Definition at line 330 of file boundBoxI.H.

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

Here is the call graph for this function:

◆ add() [3/9]

void add ( const point & pt0,
const point & pt1 )
inline

Extend to include two additional points.

Definition at line 337 of file boundBoxI.H.

References Foam::add().

Here is the call graph for this function:

◆ add() [4/9]

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

Extend to include two additional points.

Definition at line 344 of file boundBoxI.H.

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

Here is the call graph for this function:

◆ add() [5/9]

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

Extend to include the points.

Definition at line 351 of file boundBoxI.H.

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

Here is the call graph for this function:

◆ add() [6/9]

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

Extend to include the points from the temporary point field.

Definition at line 360 of file boundBoxI.H.

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

Here is the call graph for this function:

◆ add() [7/9]

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

Extend to include the points.

References points().

Here is the call graph for this function:

◆ add() [8/9]

template<unsigned N>
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.

References points().

Here is the call graph for this function:

◆ add() [9/9]

template<class IntContainer>
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

References delta, and points().

Here is the call graph for this function:

◆ grow() [1/2]

◆ grow() [2/2]

void grow ( const vector & delta)
inline

Expand box by adjusting min/max by specified amounts.

Definition at line 374 of file boundBoxI.H.

References delta.

◆ inflate() [1/3]

◆ inflate() [2/3]

void inflate ( Random & rndGen,
const scalar factor )

Expand box slightly by expanding all dimensions with factor*span*(random 0-1) and guarantees factor*mag(span) minimum width in any direction.

References rndGen.

◆ inflate() [3/3]

void inflate ( Random & r,
const scalar factor,
const scalar delta )

As per two parameter version but with additional grow() by given amount in each dimension.

References delta.

◆ reduce() [1/2]

void reduce ( )

Inplace parallel reduction of min/max values, using UPstream::worldComm.

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

Here is the caller graph for this function:

◆ reduce() [2/2]

void reduce ( int communicator)

Inplace parallel reduction of min/max values, using the specified communicator.

◆ returnReduce() [1/2]

boundBox returnReduce ( const boundBox & bb)
static

Perform a reduction on a copy and return the result, using UPstream::worldComm.

References boundBox().

Here is the call graph for this function:

◆ returnReduce() [2/2]

boundBox returnReduce ( const boundBox & bb,
int communicator )
static

Perform a reduction on a copy and return the result, using the specified communicator.

References boundBox().

Here is the call graph for this function:

◆ intersects() [1/2]

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

Referenced by cuttingPlane::checkOverlap(), and treeBoundBox::subOverlaps().

Here is the caller graph for this function:

◆ intersects() [2/2]

bool intersects ( const triPointRef & tri) const

Does triangle intersect this bounding box or is contained within this bounding box.

Note
results may be unreliable when it coincides almost exactly with a box face

References boundBox(), and centre().

Here is the call graph for this function:

◆ overlaps() [1/2]

bool overlaps ( const boundBox & bb) const
inline

◆ overlaps() [2/2]

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

Overlaps boundingSphere (centre + sqr(radius))?

Definition at line 445 of file boundBoxI.H.

References box_sphere_overlaps(), and centre().

Here is the call graph for this function:

◆ contains() [1/5]

bool contains ( const point & pt) const
inline

◆ contains() [2/5]

bool contains ( const boundBox & bb) const
inline

Fully contains other boundingBox?

Definition at line 466 of file boundBoxI.H.

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

Here is the call graph for this function:

◆ containsInside()

bool containsInside ( const point & pt) const
inline

Contains point? (inside only).

Definition at line 472 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).

References points().

Here is the call graph for this function:

◆ contains() [4/5]

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).

References points().

Here is the call graph for this function:

◆ contains() [5/5]

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

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

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

References points().

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).

References points().

Referenced by searchableRotatedBox::overlaps(), 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]

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).

References points().

Here is the call graph for this function:

◆ containsAny() [3/3]

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

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

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

References points().

Here is the call graph for this function:

◆ nearest()

point nearest ( const point & p) const

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

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

References p.

Referenced by treeBoundBox::calcExtremities().

Here is the caller graph for this function:

◆ operator&=()

void operator&= ( const boundBox & bb)

Restrict min/max to union with other box.

References boundBox().

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.

Can be used in a reduction operation.

Definition at line 671 of file boundBox.H.

References add(), and boundBox().

Here is the call graph for this function:

◆ intersect()

bool intersect ( const boundBox & bb)
inline

Deprecated(2022-10) - use 'operator&=' to avoid confusion with other intersects() methods.

Deprecated
(2022-10) - use 'operator&='

Definition at line 689 of file boundBox.H.

References boundBox(), good(), and intersect().

Referenced by intersect().

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

◆ midpoint()

point midpoint ( ) const
inline

Identical to centre().

Definition at line 694 of file boundBox.H.

References centre().

Here is the call graph for this function:

◆ hexCorner() [2/10]

template<>
Foam::point hexCorner ( ) const
inline

Definition at line 32 of file boundBoxI.H.

◆ hexCorner() [3/10]

template<>
Foam::point hexCorner ( ) const
inline

Definition at line 39 of file boundBoxI.H.

◆ hexCorner() [4/10]

template<>
Foam::point hexCorner ( ) const
inline

Definition at line 46 of file boundBoxI.H.

◆ hexCorner() [5/10]

template<>
Foam::point hexCorner ( ) const
inline

Definition at line 53 of file boundBoxI.H.

◆ hexCorner() [6/10]

template<>
Foam::point hexCorner ( ) const
inline

Definition at line 60 of file boundBoxI.H.

◆ hexCorner() [7/10]

template<>
Foam::point hexCorner ( ) const
inline

Definition at line 67 of file boundBoxI.H.

◆ hexCorner() [8/10]

template<>
Foam::point hexCorner ( ) const
inline

Definition at line 74 of file boundBoxI.H.

◆ hexCorner() [9/10]

template<>
Foam::point hexCorner ( ) const
inline

Definition at line 81 of file boundBoxI.H.

◆ hexCorner() [10/10]

template<Foam::direction CornerNumber>
Foam::point hexCorner ( ) const
inline

Definition at line 92 of file boundBoxI.H.

◆ operator>>

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

References boundBox().

◆ operator<<

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

References boundBox(), and os().

Member Data Documentation

◆ greatBox

const boundBox greatBox
static

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

Definition at line 126 of file boundBox.H.

◆ invertedBox

const boundBox invertedBox
static

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

Definition at line 131 of file boundBox.H.

Referenced by boundBox(), null(), treeBoundBox::null(), and reset().

◆ faceNormals

const FixedList<vector, 6> faceNormals
static

The unit normal per face.

Definition at line 136 of file boundBox.H.

Referenced by searchableBox::getNormal().


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