treeBoundBox Class Reference

Standard boundBox with extra functionality for use in octree. More...

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

Public Types

enum  octantBit : direction { RIGHTHALF = 0x1, TOPHALF = 0x2, FRONTHALF = 0x4 }
 Bits used for octant/point encoding. More...
 
enum  faceId {
  LEFT = 0, RIGHT = 1, BOTTOM = 2, TOP = 3,
  BACK = 4, FRONT = 5
}
 
enum  faceBit : direction {
  NOFACE = 0, LEFTBIT = 0x1 << LEFT, RIGHTBIT = 0x1 << RIGHT, BOTTOMBIT = 0x1 << BOTTOM,
  TOPBIT = 0x1 << TOP, BACKBIT = 0x1 << BACK, FRONTBIT = 0x1 << FRONT
}
 Bits used for face encoding. More...
 
enum  edgeId {
  E01 = 0, E13 = 1, E23 = 2, E02 = 3,
  E45 = 4, E57 = 5, E67 = 6, E46 = 7,
  E04 = 8, E15 = 9, E37 = 10, E26 = 11
}
 Edges codes. More...
 
- Public Types inherited from boundBox
enum  directionBit : direction { XDIR = 0x1, YDIR = 0x2, ZDIR = 0x4 }
 Bits used for (x/y/z) direction encoding. More...
 

Public Member Functions

 treeBoundBox ()
 Construct without any points - an inverted bounding box. More...
 
 treeBoundBox (const boundBox &bb)
 Construct from a boundBox. More...
 
 treeBoundBox (const point &pt)
 Construct a bounding box containing a single initial point. More...
 
 treeBoundBox (const point &min, const point &max)
 Construct from components. More...
 
 treeBoundBox (const UList< point > &points)
 Construct as the bounding box of the given pointField. More...
 
 treeBoundBox (const UList< point > &points, const labelUList &indices)
 Construct as subset of points. More...
 
template<unsigned N>
 treeBoundBox (const UList< point > &points, const FixedList< label, N > &indices)
 Construct as subset of points. More...
 
 treeBoundBox (Istream &is)
 Construct from Istream. More...
 
scalar typDim () const
 Typical dimension length,height,width. More...
 
tmp< pointFieldpoints () const
 Vertex coordinates. In octant coding. More...
 
point corner (const direction octant) const
 Corner point of given octant. More...
 
treeBoundBox subBbox (const direction octant) const
 Sub-box of given octant. Midpoint calculated. More...
 
treeBoundBox subBbox (const point &mid, const direction) const
 Sub-box given by octant number. Midpoint provided. More...
 
direction subOctant (const point &pt) const
 Returns octant number given point and the calculated midpoint. More...
 
direction subOctant (const point &pt, bool &onEdge) const
 Returns octant number given point and the calculated midpoint. More...
 
void searchOrder (const point &pt, FixedList< direction, 8 > &octantOrder) const
 Calculates optimal order to look for nearest to point. More...
 
bool intersects (const point &overallStart, const vector &overallVec, const point &start, const point &end, point &pt, direction &ptBits) const
 Intersects segment; set point to intersection position and face,. More...
 
bool intersects (const point &start, const point &end, point &pt) const
 Like above but does not return faces point is on. More...
 
bool contains (const vector &dir, const point &) const
 Contains point (inside or on edge) and moving in direction. More...
 
direction faceBits (const point &pt) const
 Code position of point on bounding box faces. More...
 
direction posBits (const point &pt) const
 Position of point relative to bounding box. More...
 
void calcExtremities (const point &pt, point &nearest, point &furthest) const
 Calculate nearest and furthest (to point) vertex coords of. More...
 
scalar maxDist (const point &pt) const
 Returns distance point to furthest away corner. More...
 
label distanceCmp (const point &pt, const treeBoundBox &other) const
 Compare distance to point with other bounding box. More...
 
treeBoundBox extend (Random &rndGen, const scalar s) const
 Return slightly wider bounding box. More...
 
bool overlaps (const boundBox &bb) const
 Overlaps other bounding box? More...
 
bool overlaps (const point &centre, const scalar radiusSqr) const
 Overlaps other bounding box? More...
 
bool contains (const point &pt) const
 Contains point or other bounding box? More...
 
bool contains (const boundBox &bb) const
 Contains point or other bounding box? More...
 
bool contains (const UList< point > &points) const
 Contains point or other bounding box? More...
 
template<unsigned N>
bool contains (const UList< point > &points, const FixedList< label, N > &indices) const
 Contains point or other bounding box? More...
 
template<class IntContainer >
bool contains (const UList< point > &points, const IntContainer &indices) const
 Contains point or other bounding box? More...
 
- Public Member Functions inherited from boundBox
 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 Member Functions

static direction subOctant (const point &mid, const point &pt)
 Returns octant number given point and midpoint. More...
 
static direction subOctant (const point &mid, const point &pt, bool &onEdge)
 Returns octant number given point and midpoint. More...
 
static direction subOctant (const point &mid, const vector &dir, const point &pt, bool &onEdge)
 Returns octant number given intersection and midpoint. More...
 

Static Public Attributes

static const faceList faces
 Face to point addressing. More...
 
static const edgeList edges
 Edge to point addressing. More...
 
- Static Public Attributes inherited from boundBox
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, treeBoundBox &bb)
 
Ostreamoperator<< (Ostream &os, const treeBoundBox &bb)
 

Detailed Description

Standard boundBox with extra functionality for use in octree.

Numbering of corner points is according to octant numbering.

On the back plane (z=0):

    Y
    ^
    |
    +--------+
    |2      3|
    |        |
    |        |
    |        |
    |0      1|
    +--------+->X

For the front plane add 4 to the point labels.

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.
Source files

Definition at line 86 of file treeBoundBox.H.

Member Enumeration Documentation

◆ octantBit

Bits used for octant/point encoding.

Every octant/corner point is the combination of three directions.

Enumerator
RIGHTHALF 

1: positive x-direction

TOPHALF 

2: positive y-direction

FRONTHALF 

4: positive z-direction

Definition at line 96 of file treeBoundBox.H.

◆ faceId

enum faceId

Face codes. Identical order and meaning as per hex cellmodel and boundBox, but have a different point ordering.

Enumerator
LEFT 

0: x-min, left

RIGHT 

1: x-max, right

BOTTOM 

2: y-min, bottom

TOP 

3: y-max, top

BACK 

4: z-min, back

FRONT 

5: z-max, front

Definition at line 105 of file treeBoundBox.H.

◆ faceBit

Bits used for face encoding.

Enumerator
NOFACE 
LEFTBIT 

1: x-min, left

RIGHTBIT 

2: x-max, right

BOTTOMBIT 

4: y-min, bottom

TOPBIT 

8: y-max, top

BACKBIT 

16: z-min, back

FRONTBIT 

32: z-max, front

Definition at line 116 of file treeBoundBox.H.

◆ edgeId

enum edgeId

Edges codes.

E01 = edge between 0 and 1.

Enumerator
E01 
E13 
E23 
E02 
E45 
E57 
E67 
E46 
E04 
E15 
E37 
E26 

Definition at line 129 of file treeBoundBox.H.

Constructor & Destructor Documentation

◆ treeBoundBox() [1/8]

treeBoundBox ( )
inline

Construct without any points - an inverted bounding box.

Definition at line 34 of file treeBoundBoxI.H.

◆ treeBoundBox() [2/8]

treeBoundBox ( const boundBox bb)
inlineexplicit

Construct from a boundBox.

Definition at line 40 of file treeBoundBoxI.H.

◆ treeBoundBox() [3/8]

treeBoundBox ( const point pt)
inlineexplicit

Construct a bounding box containing a single initial point.

Definition at line 46 of file treeBoundBoxI.H.

◆ treeBoundBox() [4/8]

treeBoundBox ( const point min,
const point max 
)
inline

Construct from components.

Definition at line 52 of file treeBoundBoxI.H.

◆ treeBoundBox() [5/8]

treeBoundBox ( const UList< point > &  points)
explicit

Construct as the bounding box of the given pointField.

Local processor domain only (no reduce as in boundBox)

Definition at line 62 of file treeBoundBox.C.

References Foam::nl, treeBoundBox::points(), and WarningInFunction.

Here is the call graph for this function:

◆ treeBoundBox() [6/8]

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

Construct as subset of points.

Local processor domain only (no reduce as in boundBox)

Definition at line 75 of file treeBoundBox.C.

References UList< T >::empty(), Foam::nl, points, and WarningInFunction.

Here is the call graph for this function:

◆ treeBoundBox() [7/8]

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

Construct as subset of points.

The indices could be from edge/triFace etc. Local processor domain only (no reduce as in boundBox)

Definition at line 35 of file treeBoundBoxTemplates.C.

References Foam::nl, points, and WarningInFunction.

◆ treeBoundBox() [8/8]

treeBoundBox ( Istream is)
inline

Construct from Istream.

Definition at line 58 of file treeBoundBoxI.H.

Member Function Documentation

◆ typDim()

Foam::scalar typDim ( ) const
inline

Typical dimension length,height,width.

Definition at line 66 of file treeBoundBoxI.H.

◆ points()

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

Vertex coordinates. In octant coding.

Definition at line 92 of file treeBoundBox.C.

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

Referenced by searchableBox::boundingSpheres(), searchableBox::coordinates(), triangleFuncs::intersectBb(), searchableBox::points(), treeBoundBox::treeBoundBox(), OBJstream::write(), box::writeBoxes(), AABBTree< Type >::writeOBJ(), and Foam::meshTools::writeOBJ().

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

◆ corner()

Foam::point corner ( const direction  octant) const
inline

Corner point of given octant.

Definition at line 72 of file treeBoundBoxI.H.

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

Here is the call graph for this function:

◆ subBbox() [1/2]

Foam::treeBoundBox subBbox ( const direction  octant) const

Sub-box of given octant. Midpoint calculated.

Definition at line 106 of file treeBoundBox.C.

Referenced by dynamicIndexedOctree< Foam::dynamicTreeDataPoint >::insertIndex(), dynamicIndexedOctree< Foam::dynamicTreeDataPoint >::print(), indexedOctree< Foam::treeDataPrimitivePatch< PatchType > >::print(), box::refineBox(), and dynamicIndexedOctree< Foam::dynamicTreeDataPoint >::removeIndex().

Here is the caller graph for this function:

◆ subBbox() [2/2]

Foam::treeBoundBox subBbox ( const point mid,
const direction  octant 
) const

Sub-box given by octant number. Midpoint provided.

Definition at line 113 of file treeBoundBox.C.

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

Here is the call graph for this function:

◆ subOctant() [1/5]

Foam::direction subOctant ( const point pt) const
inline

Returns octant number given point and the calculated midpoint.

Definition at line 84 of file treeBoundBoxI.H.

◆ subOctant() [2/5]

Foam::direction subOctant ( const point mid,
const point pt 
)
inlinestatic

Returns octant number given point and midpoint.

Definition at line 93 of file treeBoundBoxI.H.

References treeBoundBox::FRONTHALF, treeBoundBox::RIGHTHALF, treeBoundBox::TOPHALF, Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().

Here is the call graph for this function:

◆ subOctant() [3/5]

Foam::direction subOctant ( const point pt,
bool onEdge 
) const
inline

Returns octant number given point and the calculated midpoint.

onEdge set if the point is on edge of subOctant

Definition at line 122 of file treeBoundBoxI.H.

◆ subOctant() [4/5]

Foam::direction subOctant ( const point mid,
const point pt,
bool onEdge 
)
inlinestatic

Returns octant number given point and midpoint.

onEdge set if the point is on edge of subOctant

Definition at line 134 of file treeBoundBoxI.H.

References treeBoundBox::FRONTHALF, treeBoundBox::RIGHTHALF, treeBoundBox::TOPHALF, Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().

Here is the call graph for this function:

◆ subOctant() [5/5]

Foam::direction subOctant ( const point mid,
const vector dir,
const point pt,
bool onEdge 
)
inlinestatic

Returns octant number given intersection and midpoint.

onEdge set if the point is on edge of subOctant If onEdge, the direction vector determines which octant to use (acc. to which octant the point would be if it were moved along dir)

Definition at line 179 of file treeBoundBoxI.H.

References treeBoundBox::FRONTHALF, treeBoundBox::RIGHTHALF, treeBoundBox::TOPHALF, Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().

Here is the call graph for this function:

◆ searchOrder()

void searchOrder ( const point pt,
FixedList< direction, 8 > &  octantOrder 
) const
inline

Calculates optimal order to look for nearest to point.

First will be the octant containing the point, second the octant with boundary nearest to the point etc.

Definition at line 235 of file treeBoundBoxI.H.

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

Here is the call graph for this function:

◆ intersects() [1/2]

bool intersects ( const point overallStart,
const vector overallVec,
const point start,
const point end,
point pt,
direction ptBits 
) const

Intersects segment; set point to intersection position and face,.

return true if intersection found. (pt argument used during calculation even if not intersecting). Calculates intersections from outside supplied vector (overallStart, overallVec). This is so when e.g. tracking through lots of consecutive boxes (typical octree) we're not accumulating truncation errors. Set to start, (end-start) if not used.

Definition at line 162 of file treeBoundBox.C.

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

Referenced by indexedOctree< Foam::treeDataPrimitivePatch< PatchType > >::findLine(), triangleFuncs::intersectBb(), searchableRotatedBox::overlaps(), treeDataEdge::overlaps(), and streamLineBase::trimToBox().

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

◆ intersects() [2/2]

bool intersects ( const point start,
const point end,
point pt 
) const

Like above but does not return faces point is on.

Definition at line 309 of file treeBoundBox.C.

References stdFoam::end().

Here is the call graph for this function:

◆ contains() [1/6]

bool contains ( const vector dir,
const point pt 
) const

Contains point (inside or on edge) and moving in direction.

dir would cause it to go inside.

Definition at line 320 of file treeBoundBox.C.

References Foam::max(), Foam::min(), and VectorSpace< Vector< Cmpt >, Cmpt, 3 >::nComponents.

Referenced by waveMethod::calculate(), surfaceFeatures::deleteBox(), dynamicTreeDataPoint::findNearest(), dynamicTreeDataPoint::overlaps(), treeDataPoint::overlaps(), AABBTree< Type >::pointInside(), streamLineBase::trimToBox(), and propellerInfo::updateSampleDiskCells().

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

◆ faceBits()

Foam::direction faceBits ( const point pt) const

Code position of point on bounding box faces.

Definition at line 358 of file treeBoundBox.C.

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

Referenced by indexedOctree< Foam::treeDataPrimitivePatch< PatchType > >::findLine().

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

◆ posBits()

Foam::direction posBits ( const point pt) const

Position of point relative to bounding box.

Definition at line 393 of file treeBoundBox.C.

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

Referenced by treeDataPrimitivePatch< PatchType >::findIntersection(), indexedOctree< Foam::treeDataPrimitivePatch< PatchType > >::findLine(), treeDataCell::findIntersectOp::operator()(), and treeDataFace::findIntersectOp::operator()().

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

◆ calcExtremities()

void calcExtremities ( const point pt,
point nearest,
point furthest 
) const

Calculate nearest and furthest (to point) vertex coords of.

bounding box

Definition at line 429 of file treeBoundBox.C.

References Foam::mag(), Foam::max(), Foam::min(), and VectorSpace< Vector< Cmpt >, Cmpt, 3 >::nComponents.

Referenced by treeBoundBox::distanceCmp().

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

◆ maxDist()

Foam::scalar maxDist ( const point pt) const

Returns distance point to furthest away corner.

Definition at line 455 of file treeBoundBox.C.

References Foam::mag().

Here is the call graph for this function:

◆ distanceCmp()

Foam::label distanceCmp ( const point pt,
const treeBoundBox other 
) const

Compare distance to point with other bounding box.

return: -1 : all vertices of my bounding box are nearer than any of other +1 : all vertices of my bounding box are further away than any of other 0 : none of the above.

Definition at line 465 of file treeBoundBox.C.

References treeBoundBox::calcExtremities(), Foam::sqr(), Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().

Here is the call graph for this function:

◆ extend()

Foam::treeBoundBox extend ( Random rndGen,
const scalar  s 
) const
inline

Return slightly wider bounding box.

Extends all dimensions with s*span*Random::sample01<scalar>() and guarantees in any direction s*mag(span) minimum width

Definition at line 325 of file treeBoundBoxI.H.

References Foam::cmptMultiply(), Foam::mag(), boundBox::max(), Foam::max(), boundBox::min(), VectorSpace< Vector< scalar >, scalar, 3 >::nComponents, rndGen, s, Random::sample01(), and boundBox::span().

Referenced by polyMesh::cellTree(), meshRefinement::distribute(), triSurfaceMesh::edgeTree(), patchProbes::findElements(), refinementFeatures::regionEdgeTrees(), triSurfaceSearch::tree(), and triSurfaceRegionSearch::treeByRegion().

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

◆ overlaps() [1/2]

bool overlaps
inline

Overlaps other bounding box?

Definition at line 221 of file boundBoxI.H.

Referenced by trackingInverseDistance::markDonors(), inverseDistance::markDonors(), trackingInverseDistance::markPatchesAsHoles(), inverseDistance::markPatchesAsHoles(), AABBTree< Type >::overlaps(), treeDataFace::overlaps(), treeDataCell::overlaps(), and treeDataPrimitivePatch< PatchType >::overlaps().

Here is the caller graph for this function:

◆ overlaps() [2/2]

bool overlaps
inline

Overlaps other bounding box?

Definition at line 233 of file boundBoxI.H.

◆ contains() [2/6]

bool contains
inline

Contains point or other bounding box?

Definition at line 271 of file boundBoxI.H.

◆ contains() [3/6]

bool contains
inline

Contains point or other bounding box?

Definition at line 282 of file boundBoxI.H.

◆ contains() [4/6]

bool contains

Contains point or other bounding box?

Definition at line 230 of file boundBox.C.

◆ contains() [5/6]

bool contains ( unsigned  N)
inline

Contains point or other bounding box?

Definition at line 120 of file boundBoxTemplates.C.

◆ contains() [6/6]

bool contains ( class IntContainer  )
inline

Contains point or other bounding box?

Definition at line 149 of file boundBoxTemplates.C.

Friends And Related Function Documentation

◆ operator>>

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

◆ operator<<

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

Member Data Documentation

◆ faces

◆ edges


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