treeBoundBox Class Reference

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

Inheritance diagram for treeBoundBox:
Collaboration diagram for treeBoundBox:

Classes

struct  edgeId
 Edges codes. More...
 
struct  faceBit
 Bits used for face coding. More...
 
struct  faceId
 Face codes. More...
 
struct  octantBit
 Bits used for octant/point coding. More...
 

Public Member Functions

 treeBoundBox ()
 Construct null setting points to zero. More...
 
 treeBoundBox (const boundBox &bb)
 Construct from a boundBox. More...
 
 treeBoundBox (const point &min, const point &max)
 Construct from components. More...
 
 treeBoundBox (const UList< point > &)
 Construct as the bounding box of the given pointField. More...
 
 treeBoundBox (const UList< point > &, const labelUList &indices)
 Construct as subset of points. More...
 
template<unsigned Size>
 treeBoundBox (const UList< point > &, const FixedList< label, Size > &indices)
 Construct as subset of points. More...
 
 treeBoundBox (Istream &)
 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) const
 Corner point given octant. More...
 
treeBoundBox subBbox (const direction) const
 Sub box given by octant number. 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 &) const
 Code position of point on bounding box faces. More...
 
direction posBits (const point &) 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 &) const
 Returns distance point to furthest away corner. More...
 
label distanceCmp (const point &, const treeBoundBox &other) const
 Compare distance to point with other bounding box. More...
 
treeBoundBox extend (const scalar s) const
 Return asymetrically extended bounding box, with guaranteed. More...
 
void writeOBJ (const fileName &fName) const
 
- Public Member Functions inherited from boundBox
 boundBox ()
 Construct null, setting points to zero. More...
 
 boundBox (const point &min, const point &max)
 Construct from components. More...
 
 boundBox (const UList< point > &, const bool doReduce=true)
 Construct as the bounding box of the given points. More...
 
 boundBox (const tmp< pointField > &, const bool doReduce=true)
 Construct as the bounding box of the given temporary pointField. More...
 
 boundBox (const UList< point > &, const labelUList &indices, const bool doReduce=true)
 Construct bounding box as subset of the pointField. More...
 
template<unsigned Size>
 boundBox (const UList< point > &, const FixedList< label, Size > &indices, const bool doReduce=true)
 Construct bounding box as subset of the pointField. More...
 
 boundBox (Istream &)
 Construct from Istream. More...
 
const pointmin () const
 Minimum point defining the bounding box. More...
 
const pointmax () const
 Maximum point defining the bounding box. More...
 
pointmin ()
 Minimum point defining the bounding box, non-const access. More...
 
pointmax ()
 Maximum point defining the bounding box, non-const access. More...
 
point midpoint () const
 The midpoint of the bounding box. More...
 
vector span () const
 The bounding box span (from minimum to maximum) More...
 
scalar mag () const
 The magnitude of the bounding box span. More...
 
scalar volume () const
 The volume of the bound box. More...
 
scalar minDim () const
 Smallest length/height/width dimension. More...
 
scalar maxDim () const
 Largest length/height/width dimension. More...
 
scalar avgDim () const
 Average length/height/width dimension. More...
 
tmp< pointFieldpoints () const
 Return corner points in an order corresponding to a 'hex' cell. More...
 
void inflate (const scalar s)
 Inflate box by factor*mag(span) in all dimensions. More...
 
bool overlaps (const boundBox &) const
 Overlaps/touches boundingBox? More...
 
bool overlaps (const point &, const scalar radiusSqr) const
 Overlaps boundingSphere (centre + sqr(radius))? More...
 
bool contains (const point &) const
 Contains point? (inside or on edge) More...
 
bool contains (const boundBox &) const
 Fully contains other boundingBox? More...
 
bool containsInside (const point &) const
 Contains point? (inside only) More...
 
bool contains (const UList< point > &) const
 Contains all of the points? (inside or on edge) More...
 
bool contains (const UList< point > &, const labelUList &indices) const
 Contains all of the points? (inside or on edge) More...
 
template<unsigned Size>
bool contains (const UList< point > &, const FixedList< label, Size > &indices) const
 Contains all of the points? (inside or on edge) More...
 
bool containsAny (const UList< point > &) const
 Contains any of the points? (inside or on edge) More...
 
bool containsAny (const UList< point > &, const labelUList &indices) const
 Contains any of the points? (inside or on edge) More...
 
template<unsigned Size>
bool containsAny (const UList< point > &, const FixedList< label, Size > &indices) const
 Contains any of the points? (inside or on edge) More...
 
point nearest (const point &) const
 Return the nearest point on the boundBox to the supplied point. More...
 

Static Public Member Functions

static 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 Member Functions inherited from boundBox
static faceList faces ()
 Return faces with correct point order. More...
 

Static Public Attributes

static const scalar great
 The great value used for greatBox and invertedBox. More...
 
static const treeBoundBox greatBox
 As per boundBox::greatBox, but with great instead of vGreat. More...
 
static const treeBoundBox invertedBox
 As per boundBox::invertedBox, but with great instead of vGreat. More...
 
static const faceList faces
 Face to point addressing. More...
 
static const edgeList edges
 Edge to point addressing. More...
 
static const FixedList< vector, 6 > faceNormals
 Per face the unit normal. More...
 
- Static Public Attributes inherited from boundBox
static const scalar great
 The great value used for greatBox and invertedBox. More...
 
static const boundBox greatBox
 A very large boundBox: min/max == -/+ vGreat. More...
 
static const boundBox invertedBox
 A very large inverted boundBox: min/max == +/- vGreat. More...
 

Friends

bool operator== (const treeBoundBox &, const treeBoundBox &)
 
bool operator!= (const treeBoundBox &, const treeBoundBox &)
 
Istreamoperator>> (Istream &is, treeBoundBox &)
 
Ostreamoperator<< (Ostream &os, const treeBoundBox &)
 

Detailed Description

Standard boundBox + 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.

Source files

Definition at line 87 of file treeBoundBox.H.

Constructor & Destructor Documentation

◆ treeBoundBox() [1/7]

treeBoundBox ( )
inline

Construct null setting points to zero.

Definition at line 31 of file treeBoundBoxI.H.

Referenced by treeBoundBox::treeBoundBox().

Here is the caller graph for this function:

◆ treeBoundBox() [2/7]

treeBoundBox ( const boundBox bb)
inlineexplicit

Construct from a boundBox.

Definition at line 43 of file treeBoundBoxI.H.

◆ treeBoundBox() [3/7]

treeBoundBox ( const point min,
const point max 
)
inline

Construct from components.

Definition at line 37 of file treeBoundBoxI.H.

◆ treeBoundBox() [4/7]

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 130 of file treeBoundBox.C.

References UList< T >::empty(), Foam::endl(), treeBoundBox::treeBoundBox(), and WarningInFunction.

Here is the call graph for this function:

◆ treeBoundBox() [5/7]

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 146 of file treeBoundBox.C.

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

Here is the call graph for this function:

◆ treeBoundBox() [6/7]

treeBoundBox ( const UList< point > &  points,
const FixedList< label, Size > &  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 UList< T >::empty(), Foam::endl(), and WarningInFunction.

Here is the call graph for this function:

◆ treeBoundBox() [7/7]

treeBoundBox ( Istream is)
inline

Construct from Istream.

Definition at line 49 of file treeBoundBoxI.H.

Member Function Documentation

◆ typDim()

Foam::scalar typDim ( ) const
inline

Typical dimension length,height,width.

Definition at line 57 of file treeBoundBoxI.H.

References boundBox::avgDim().

Here is the call graph for this function:

◆ points()

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

Vertex coordinates. In octant coding.

Definition at line 166 of file treeBoundBox.C.

References treeBoundBox::corner(), forAll, and tmp< T >::ref().

Referenced by searchableBox::boundingSpheres(), searchableBox::coordinates(), indexedOctree< Foam::treeDataFace >::findLine(), dynamicIndexedOctree< Foam::dynamicTreeDataPoint >::findNearest(), searchableBox::points(), and treeBoundBox::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 given octant.

Definition at line 63 of file treeBoundBoxI.H.

References treeBoundBox::octantBit::frontHalf, boundBox::max(), boundBox::min(), treeBoundBox::octantBit::rightHalf, treeBoundBox::octantBit::topHalf, x, and y.

Referenced by treeBoundBox::points().

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

◆ subBbox() [1/2]

◆ subBbox() [2/2]

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

◆ subOctant() [1/5]

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

Returns octant number given point and the calculated midpoint.

Definition at line 75 of file treeBoundBoxI.H.

References boundBox::midpoint().

Referenced by treeBoundBox::subOctant().

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

◆ subOctant() [2/5]

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

Returns octant number given point and midpoint.

Definition at line 84 of file treeBoundBoxI.H.

References treeBoundBox::octantBit::frontHalf, treeBoundBox::octantBit::rightHalf, treeBoundBox::subOctant(), treeBoundBox::octantBit::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 113 of file treeBoundBoxI.H.

References boundBox::midpoint(), and treeBoundBox::subOctant().

Here is the call graph for this function:

◆ 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 125 of file treeBoundBoxI.H.

References treeBoundBox::octantBit::frontHalf, treeBoundBox::octantBit::rightHalf, treeBoundBox::subOctant(), treeBoundBox::octantBit::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 170 of file treeBoundBoxI.H.

References treeBoundBox::octantBit::frontHalf, treeBoundBox::octantBit::rightHalf, treeBoundBox::searchOrder(), treeBoundBox::octantBit::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 226 of file treeBoundBoxI.H.

References treeBoundBox::octantBit::frontHalf, boundBox::max(), boundBox::midpoint(), boundBox::min(), treeBoundBox::octantBit::rightHalf, treeBoundBox::octantBit::topHalf, Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().

Referenced by treeBoundBox::subOctant().

Here is the call graph for this function:
Here is the caller 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 236 of file treeBoundBox.C.

References treeBoundBox::faceBit::back, treeBoundBox::faceBit::bottom, treeBoundBox::faceBits(), treeBoundBox::faceBit::front, treeBoundBox::faceBit::left, Foam::mag(), boundBox::max(), boundBox::min(), treeBoundBox::posBits(), treeBoundBox::faceBit::right, s(), treeBoundBox::faceBit::top, Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().

Referenced by indexedOctree< Foam::treeDataFace >::findLine(), dynamicIndexedOctree< Foam::dynamicTreeDataPoint >::findNearest(), indexedOctree< Foam::treeDataFace >::findNearest(), triangleFuncs::intersectBb(), treeBoundBox::intersects(), NamedEnum< compressibleField, 8 >::names(), treeDataEdge::overlaps(), and treeBoundBox::subBbox().

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 383 of file treeBoundBox.C.

References treeBoundBox::intersects().

Here is the call graph for this function:

◆ contains()

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 394 of file treeBoundBox.C.

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

Referenced by dynamicTreeDataPoint::findNearest(), dynamicIndexedOctree< Foam::dynamicTreeDataPoint >::findNearest(), indexedOctree< Foam::treeDataFace >::findNearest(), conformalVoronoiMesh::maxSurfaceProtrusion(), NamedEnum< compressibleField, 8 >::names(), dynamicTreeDataPoint::overlaps(), and treeDataPoint::overlaps().

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

◆ faceBits()

◆ posBits()

◆ 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 500 of file treeBoundBox.C.

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

Referenced by treeBoundBox::distanceCmp(), treeBoundBox::maxDist(), and treeBoundBox::posBits().

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 547 of file treeBoundBox.C.

References treeBoundBox::calcExtremities(), treeBoundBox::distanceCmp(), and 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 557 of file treeBoundBox.C.

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

Referenced by treeBoundBox::maxDist().

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

◆ extend()

Foam::treeBoundBox extend ( const scalar  s) const
inline

◆ writeOBJ()

void writeOBJ ( const fileName fName) const

Friends And Related Function Documentation

◆ operator==

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

◆ operator!=

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

◆ operator>>

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

◆ operator<<

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

Member Data Documentation

◆ great

const Foam::scalar great
static

The great value used for greatBox and invertedBox.

Definition at line 102 of file treeBoundBox.H.

◆ greatBox

const Foam::treeBoundBox greatBox
static

As per boundBox::greatBox, but with great instead of vGreat.

Definition at line 105 of file treeBoundBox.H.

◆ invertedBox

const Foam::treeBoundBox invertedBox
static

As per boundBox::invertedBox, but with great instead of vGreat.

Definition at line 108 of file treeBoundBox.H.

Referenced by patchToPatch::srcBox(), and patchToPatch::tgtBox().

◆ faces

const Foam::faceList faces
static

Face to point addressing.

Definition at line 175 of file treeBoundBox.H.

◆ edges

const Foam::edgeList edges
static

Edge to point addressing.

Definition at line 178 of file treeBoundBox.H.

Referenced by treeBoundBox::writeOBJ().

◆ faceNormals

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

Per face the unit normal.

Definition at line 181 of file treeBoundBox.H.

Referenced by searchableBox::getNormal().


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