A bounding box defined in terms of min/max extrema points. More...
#include <boundBox.H>


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. | |
| boundBox & | operator= (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 point & | min () const noexcept |
| Minimum describing the bounding box. | |
| const point & | max () const noexcept |
| Maximum describing the bounding box. | |
| point & | min () noexcept |
| Minimum describing the bounding box, non-const access. | |
| point & | max () 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< pointField > | hexCorners () const |
| Corner points in an order corresponding to a 'hex' cell. | |
| tmp< pointField > | points () const |
| Corner points in an order corresponding to a 'hex' cell. | |
| tmp< pointField > | faceCentres () 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 ¢re, 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 boundBox & | null () noexcept |
| The null boundBox is the same as an inverted box. | |
| static const Foam::faceList & | hexFaces () |
| 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 ¢re, const scalar radiusSqr) |
| Test for overlap of box and boundingSphere (centre + sqr(radius)). | |
Friends | |
| Istream & | operator>> (Istream &is, boundBox &bb) |
| Ostream & | operator<< (Ostream &os, const boundBox &bb) |
A bounding box defined in terms of min/max extrema points.
Definition at line 70 of file boundBox.H.
| enum directionBit : direction |
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.
|
inline |
Default construct: an inverted bounding box.
Definition at line 101 of file boundBoxI.H.
References invertedBox, max(), and min().
Referenced by add(), boundBox(), boundBox(), contains(), intersect(), intersects(), null(), operator&=(), operator+=(), operator<<, operator=(), treeBoundBox::operator=(), operator>>, overlaps(), returnReduce(), returnReduce(), treeBoundBox::subOverlaps(), treeBoundBox::treeBoundBox(), treeBoundBox::treeBoundBox(), treeBoundBox::treeBoundBox(), treeBoundBox::treeBoundBox(), treeBoundBox::treeBoundBox(), and treeBoundBox::treeBoundBox().


|
default |
| boundBox | ( | const boundBox & | bb, |
| bool | doReduce ) |
Copy construct with specified global reduction.
References boundBox(), max(), min(), and p.

|
inlineexplicit |
Construct a bounding box containing a single initial point.
Definition at line 115 of file boundBoxI.H.
References p.
|
inlineexplicit |
Construct a 0/1 unit bounding box.
Definition at line 108 of file boundBoxI.H.
Construct from bound box min/max points.
Definition at line 122 of file boundBoxI.H.

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.
Construct as the bounding box of the given points.
Does parallel communication (doReduce = true)
References points().

|
explicit |
Construct as the bounding box of the given temporary pointField.
Does parallel communication (doReduce = true)
| 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().

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

|
inlineexplicit |
Construct from Istream.
Definition at line 136 of file boundBoxI.H.
References Foam::operator>>().

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


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


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


|
static |
The boundBox faces as a hexCell, using hexCorner points. Same as hexCell::modelFaces().
|
inlinestaticconstexprnoexcept |
Number of points for boundBox and HEX.
Definition at line 158 of file boundBox.H.
References Foam::noexcept.
|
inlinestaticconstexprnoexcept |
Number of edges for boundBox and HEX.
Definition at line 163 of file boundBox.H.
References Foam::noexcept.
|
inlinestaticconstexprnoexcept |
Number of faces for boundBox and HEX.
Definition at line 168 of file boundBox.H.
References Foam::noexcept.
|
inline |
Bounding box is inverted, contains no points.
Definition at line 144 of file boundBoxI.H.
Referenced by fieldExtents::calcFieldExtents(), cuttingSurfaceBase::cellSelection(), and good().

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


|
inline |
Bounding box is non-inverted - same as good().
Definition at line 283 of file boundBox.H.
References good().

|
inlinenoexcept |
Minimum describing the bounding box.
Definition at line 162 of file boundBoxI.H.
References Foam::noexcept.
Referenced by PDRblock::blockMeshDict(), boundBox(), boundBox(), boundBox(), eddy::bounds(), fieldExtents::calcFieldExtents(), voxelMeshSearch::centre(), advancingFrontAMI::checkPatches(), contains(), AABBTree< Type >::createBoxes(), faceCentre(), inverseDistance::fill(), voxelMeshSearch::fill(), voxelMeshSearch::fill(), boundaryMesh::getNearest(), voxelMeshSearch::index(), inverseDistance::index3(), voxelMeshSearch::index3(), uniformBin::initialise(), waveMakerPointPatchVectorField::initialiseGeometry(), waveModel::initialiseGeometry(), projectVertex::operator point(), Foam::operator<(), Foam::operator==(), overlaps(), inverseDistance::overlaps(), searchableSphere::overlaps(), voxelMeshSearch::overlaps(), inverseDistance::position(), powerLawLopesdaCostaZone::powerLawLopesdaCostaZone(), Foam::readBoxDim(), reset(), PDRblock::reset(), searchableDisk::searchableDisk(), searchableSphere::searchableSphere(), searchableSurfaceCollection::searchableSurfaceCollection(), treeBoundBox::searchOrder(), treeBoundBox::treeBoundBox(), inverseDistance::update(), trackingInverseDistance::update(), topoSet::writeDebug(), and voxelMeshSearch::writeGrid().

|
inlinenoexcept |
Maximum describing the bounding box.
Definition at line 168 of file boundBoxI.H.
References Foam::noexcept.
Referenced by PDRblock::blockMeshDict(), boundBox(), boundBox(), boundBox(), eddy::bounds(), fieldExtents::calcFieldExtents(), advancingFrontAMI::checkPatches(), contains(), faceCentre(), inverseDistance::fill(), voxelMeshSearch::fill(), voxelMeshSearch::fill(), boundaryMesh::getNearest(), uniformBin::initialise(), waveMakerPointPatchVectorField::initialiseGeometry(), waveModel::initialiseGeometry(), Kmesh::Kmesh(), projectVertex::operator point(), Foam::operator<(), Foam::operator==(), orientedSurface::orientedSurface(), overlaps(), inverseDistance::overlaps(), searchableSphere::overlaps(), voxelMeshSearch::overlaps(), powerLawLopesdaCostaZone::powerLawLopesdaCostaZone(), Foam::readBoxDim(), reset(), PDRblock::reset(), searchableDisk::searchableDisk(), searchableSphere::searchableSphere(), searchableSurfaceCollection::searchableSurfaceCollection(), treeBoundBox::searchOrder(), treeBoundBox::treeBoundBox(), and topoSet::writeDebug().

|
inlinenoexcept |
Minimum describing the bounding box, non-const access.
Definition at line 174 of file boundBoxI.H.
References Foam::noexcept.
|
inlinenoexcept |
Maximum describing the bounding box, non-const access.
Definition at line 180 of file boundBoxI.H.
References Foam::noexcept.
|
inline |
The centre (midpoint) of the bounding box.
Definition at line 186 of file boundBoxI.H.
Referenced by PDRblock::blockMeshDict(), box_sphere_overlaps(), searchableBox::findNearest(), searchableBox::findNearest(), searchableBox::findNearestOnEdge(), advancingFrontAMI::findTargetFace(), intersects(), midpoint(), overlaps(), treeBoundBox::searchOrder(), treeBoundBox::subOctant(), treeBoundBox::subOctant(), treeBoundBox::subOverlaps(), and treeBoundBox::subOverlaps().

|
inline |
The bounding box span (from minimum to maximum).
Definition at line 192 of file boundBoxI.H.
Referenced by avgDim(), PDRblock::blockMeshDict(), booleanSurface::booleanSurface(), voxelMeshSearch::centre(), advancingFrontAMI::checkPatches(), AABBTree< Type >::createBoxes(), voxelMeshSearch::index(), inverseDistance::index3(), voxelMeshSearch::index3(), uniformBin::initialise(), Kmesh::Kmesh(), maxDim(), minDim(), orientedSurface::orientedSurface(), porosityModel::porosityModel(), inverseDistance::position(), energySpectrum::read(), inverseDistance::update(), trackingInverseDistance::update(), sampledMeshedSurface::update(), points0MotionSolver::updateMesh(), volume(), and voxelMeshSearch::writeGrid().

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


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

|
inline |
The volume of the bound box.
Definition at line 210 of file boundBoxI.H.
References Foam::cmptProduct(), and span().

|
inline |
Smallest length/height/width dimension.
Definition at line 216 of file boundBoxI.H.
References Foam::cmptMin(), and span().

|
inline |
Largest length/height/width dimension.
Definition at line 222 of file boundBoxI.H.
References Foam::cmptMax(), and span().

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


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

|
inline |
Count the number of positive, non-zero dimensions.
Definition at line 274 of file boundBoxI.H.
References VectorSpace< Form, Cmpt, Ncmpts >::nComponents.
Return corner point [0..7] corresponding to a 'hex' cell.
References hexCorner().
Referenced by hexCorner().


| tmp< pointField > hexCorners | ( | ) | const |
Corner points in an order corresponding to a 'hex' cell.
References hexCorners().
Referenced by hexCorners(), points(), and ensightCells::writeBox().


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


| tmp< pointField > faceCentres | ( | ) | const |
Face midpoints.
|
inline |
Reset to an inverted box.
Definition at line 295 of file boundBoxI.H.
References invertedBox.
Referenced by PatchTools::calcBounds(), fieldExtents::calcFieldExtents(), clear(), searchableExtrudedCircle::searchableExtrudedCircle(), searchableSurfaceCollection::searchableSurfaceCollection(), and sampledMeshedSurface::update().

|
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.
|
inline |
Reset min/max to be identical to the specified point.
Definition at line 309 of file boundBoxI.H.
Reset min/max to specified values.
Definition at line 316 of file boundBoxI.H.

|
inline |
Same as reset() - reset to an inverted box.
Definition at line 419 of file boundBox.H.
References reset().
Referenced by cuttingSurfaceBase::cellSelection().


|
inline |
Extend to include the second box.
Definition at line 323 of file boundBoxI.H.
References boundBox(), Foam::max(), and Foam::min().
Referenced by searchableSurfacesQueries::bounds(), treeDataCell::bounds(), treeDataCell::bounds(), treeDataEdge::bounds(), treeDataEdge::bounds(), treeDataFace::bounds(), treeDataPoint::bounds(), Foam::boundsImpl(), Foam::boundsImpl(), PatchTools::calcBounds(), fieldExtents::calcFieldExtents(), cellBox::calcSrcBox(), cellBox::calcTgtBox(), cuttingSurfaceBase::cellSelection(), AABBTree< Type >::createBoxes(), patchFieldProbe::findElements(), uniformBin::initialise(), operator+=(), porosityModel::porosityModel(), and triSurface::writeStats().


|
inline |
Extend to include the point.
Definition at line 330 of file boundBoxI.H.
References Foam::max(), and Foam::min().

Extend to include two additional points.
Definition at line 337 of file boundBoxI.H.
References Foam::add().

Extend to include two additional points.
Definition at line 344 of file boundBoxI.H.
References Foam::add(), and points().

Extend to include the points.
Definition at line 351 of file boundBoxI.H.
References Foam::add(), p, and points().

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

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

|
inline |
Expand box by adjusting min/max by specified amount in each dimension.
Definition at line 367 of file boundBoxI.H.
References delta.
Referenced by faceAreaWeightAMI2D::calculate(), Foam::createTree(), dynamicTreeDataPoint::findNearest(), inflate(), inverseDistance::markBoundaries(), trackingInverseDistance::markBoundaries(), inverseDistance::markDonors(), trackingInverseDistance::markDonors(), inverseDistance::markPatchesAsHoles(), trackingInverseDistance::markPatchesAsHoles(), treeDataEdge::findNearestOp::operator()(), treeDataPoint::findNearestOp::operator()(), searchableExtrudedCircle::searchableExtrudedCircle(), and sampledMeshedSurface::update().

|
inline |
Expand box by adjusting min/max by specified amounts.
Definition at line 374 of file boundBoxI.H.
References delta.
|
inline |
Expand box by factor*mag(span) in all dimensions.
Definition at line 381 of file boundBoxI.H.
References grow(), and Foam::mag().
Referenced by AABBTree< Type >::AABBTree(), polyMesh::cellTree(), advancingFrontAMI::checkPatches(), AABBTree< Type >::createBoxes(), AMIInterpolation::createTree(), edgeSlipDisplacementPointPatchVectorField::edgeTree(), triSurfaceMesh::edgeTree(), treeBoundBox::extend(), treeBoundBox::extend(), patchFieldProbe::findElements(), meshToMeshMethod::maskCells(), surfaceFeatures::nearestFeatEdge(), surfaceFeatures::nearestSurfEdge(), primitiveMesh::pointInCellBB(), pointAttractionDisplacementPointPatchVectorField::pointTree(), refinementFeatures::regionEdgeTrees(), Foam::simpleGeometricFilter(), triSurfaceSearch::tree(), triSurfaceRegionSearch::treeByRegion(), sampledMeshedSurface::update(), and propellerInfo::updateSampleDiskCells().


| 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.
| void inflate | ( | Random & | r, |
| const scalar | factor, | ||
| const scalar | delta ) |
| 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().

| void reduce | ( | int | communicator | ) |
Inplace parallel reduction of min/max values, using the specified communicator.
Perform a reduction on a copy and return the result, using UPstream::worldComm.
References boundBox().

Perform a reduction on a copy and return the result, using the specified communicator.
References boundBox().

| bool intersects | ( | const plane & | pln | ) | const |
Does plane intersect this bounding box.
There is an intersection if the plane segments the corner points
Referenced by cuttingPlane::checkOverlap(), and treeBoundBox::subOverlaps().

| bool intersects | ( | const triPointRef & | tri | ) | const |
Does triangle intersect this bounding box or is contained within this bounding box.
References boundBox(), and centre().

|
inline |
Overlaps/touches boundingBox?
Definition at line 439 of file boundBoxI.H.
References boundBox(), box_box_overlaps(), max(), and min().
Referenced by cuttingSurfaceBase::checkOverlap(), inverseDistance::markBoundaries(), trackingInverseDistance::markBoundaries(), meshToMeshMethod::maskCells(), faceAreaWeightAMI2D::overlappingTgtFaces(), searchableBox::overlaps(), searchablePlate::overlaps(), searchableRotatedBox::overlaps(), searchableSphere::overlaps(), and treeBoundBox::subHalf().


|
inline |
Overlaps boundingSphere (centre + sqr(radius))?
Definition at line 445 of file boundBoxI.H.
References box_sphere_overlaps(), and centre().

|
inline |
Contains point? (inside or on edge).
Definition at line 455 of file boundBoxI.H.
References Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().
Referenced by isoSurfaceBase::blockCells(), cuttingSurfaceBase::cellSelection(), advancingFrontAMI::checkPatches(), contains(), treeBoundBox::intersects(), trackingInverseDistance::markBoundaries(), primitiveMesh::pointInCellBB(), powerLawLopesdaCostaZone::powerLawLopesdaCostaZone(), and Foam::simpleGeometricFilter().


|
inline |
Fully contains other boundingBox?
Definition at line 466 of file boundBoxI.H.
References boundBox(), contains(), max(), and min().

|
inline |
Contains point? (inside only).
Definition at line 472 of file boundBoxI.H.
References Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().

Contains all points? (inside or on edge).
References points().

| bool contains | ( | const UList< point > & | points, |
| const FixedList< label, N > & | indices ) const |
Contains all of the (subsetted) points? (inside or on edge).
References points().

| bool contains | ( | const UList< point > & | points, |
| const IntContainer & | indices ) const |
Contains all of the (subsetted) points? (inside or on edge).
| IntContainer | A container with an iterator that dereferences to an label |
References points().

Contains any of the points? (inside or on edge).
References points().
Referenced by searchableRotatedBox::overlaps(), treeDataFace::overlaps(), treeDataPrimitivePatch< PatchType >::overlaps(), and triSurfaceTools::triangulate().


| bool containsAny | ( | const UList< point > & | points, |
| const FixedList< label, N > & | indices ) const |
Contains any of the (subsetted) points? (inside or on edge).
References points().

| bool containsAny | ( | const UList< point > & | points, |
| const IntContainer & | indices ) const |
Contains any of the (subsetted) points? (inside or on edge).
| IntContainer | A container with an iterator that dereferences to an label |
References points().

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

| void operator&= | ( | const boundBox & | bb | ) |
Restrict min/max to union with other box.
References boundBox().

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

|
inline |
Deprecated(2022-10) - use 'operator&=' to avoid confusion with other intersects() methods.
Definition at line 689 of file boundBox.H.
References boundBox(), good(), and intersect().
Referenced by intersect().


|
inline |
Identical to centre().
Definition at line 694 of file boundBox.H.
References centre().

|
inline |
Definition at line 32 of file boundBoxI.H.
|
inline |
Definition at line 39 of file boundBoxI.H.
|
inline |
Definition at line 46 of file boundBoxI.H.
|
inline |
Definition at line 53 of file boundBoxI.H.
|
inline |
Definition at line 60 of file boundBoxI.H.
|
inline |
Definition at line 67 of file boundBoxI.H.
|
inline |
Definition at line 74 of file boundBoxI.H.
|
inline |
Definition at line 81 of file boundBoxI.H.
|
inline |
Definition at line 92 of file boundBoxI.H.
References boundBox().
References boundBox(), and os().
|
static |
A large boundBox: min/max == -/+ ROOTVGREAT.
Definition at line 126 of file boundBox.H.
|
static |
A large inverted boundBox: min/max == +/- ROOTVGREAT.
Definition at line 131 of file boundBox.H.
Referenced by boundBox(), null(), treeBoundBox::null(), and reset().
The unit normal per face.
Definition at line 136 of file boundBox.H.
Referenced by searchableBox::getNormal().