Standard boundBox + extra functionality for use in octree. More...
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< pointField > | points () 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 point & | min () const |
Minimum describing the bounding box. More... | |
const point & | max () const |
Maximum describing the bounding box. More... | |
point & | min () |
Minimum describing the bounding box, non-const access. More... | |
point & | max () |
Maximum describing 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< pointField > | points () 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 &) |
Istream & | operator>> (Istream &is, treeBoundBox &) |
Ostream & | operator<< (Ostream &os, const treeBoundBox &) |
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.
Definition at line 87 of file treeBoundBox.H.
|
inline |
Construct null setting points to zero.
Definition at line 31 of file treeBoundBoxI.H.
Referenced by treeBoundBox::treeBoundBox().
|
inlineexplicit |
Construct from a boundBox.
Definition at line 43 of file treeBoundBoxI.H.
|
inline |
Construct from components.
Definition at line 37 of file treeBoundBoxI.H.
|
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.
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.
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.
|
inline |
Construct from Istream.
Definition at line 49 of file treeBoundBoxI.H.
|
inline |
Typical dimension length,height,width.
Definition at line 57 of file treeBoundBoxI.H.
References boundBox::avgDim().
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().
|
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().
Foam::treeBoundBox subBbox | ( | const direction | octant | ) | const |
Sub box given by octant number. Midpoint calculated.
Definition at line 180 of file treeBoundBox.C.
References boundBox::midpoint().
Referenced by indexedOctree< Foam::treeDataFace >::findLine(), dynamicIndexedOctree< Foam::dynamicTreeDataPoint >::findNearest(), indexedOctree< Foam::treeDataFace >::findNearest(), dynamicIndexedOctree< Foam::dynamicTreeDataPoint >::insertIndex(), dynamicIndexedOctree< Foam::dynamicTreeDataPoint >::overlaps(), indexedOctree< Foam::treeDataFace >::overlaps(), dynamicIndexedOctree< Foam::dynamicTreeDataPoint >::print(), indexedOctree< Foam::treeDataFace >::print(), and dynamicIndexedOctree< Foam::dynamicTreeDataPoint >::removeIndex().
Foam::treeBoundBox subBbox | ( | const point & | mid, |
const direction | octant | ||
) | const |
Sub box given by octant number. Midpoint provided.
Definition at line 187 of file treeBoundBox.C.
References Foam::abort(), Foam::FatalError, FatalErrorInFunction, treeBoundBox::octantBit::frontHalf, treeBoundBox::intersects(), boundBox::max(), boundBox::min(), treeBoundBox::octantBit::rightHalf, treeBoundBox::octantBit::topHalf, Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().
|
inline |
Returns octant number given point and the calculated midpoint.
Definition at line 75 of file treeBoundBoxI.H.
References boundBox::midpoint().
Referenced by dynamicIndexedOctree< Foam::dynamicTreeDataPoint >::overlaps(), and treeBoundBox::subOctant().
|
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().
|
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().
|
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().
|
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().
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().
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().
Like above but does not return faces point is on.
Definition at line 383 of file treeBoundBox.C.
References treeBoundBox::intersects().
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().
Foam::direction faceBits | ( | const point & | pt | ) | const |
Code position of point on bounding box faces.
Definition at line 432 of file treeBoundBox.C.
References treeBoundBox::faceBit::back, treeBoundBox::faceBit::bottom, treeBoundBox::faceBit::front, treeBoundBox::faceBit::left, boundBox::max(), boundBox::min(), treeBoundBox::faceBit::right, 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(), and treeBoundBox::intersects().
Foam::direction posBits | ( | const point & | pt | ) | const |
Position of point relative to bounding box.
Definition at line 465 of file treeBoundBox.C.
References treeBoundBox::faceBit::back, treeBoundBox::faceBit::bottom, treeBoundBox::calcExtremities(), treeBoundBox::faceBit::front, treeBoundBox::faceBit::left, boundBox::max(), boundBox::min(), treeBoundBox::faceBit::right, treeBoundBox::faceBit::top, Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().
Referenced by treeDataPrimitivePatch< PatchType >::findIntersection(), indexedOctree< Foam::treeDataFace >::findLine(), dynamicIndexedOctree< Foam::dynamicTreeDataPoint >::findNearest(), indexedOctree< Foam::treeDataFace >::findNearest(), treeBoundBox::intersects(), treeDataCell::findIntersectOp::operator()(), and treeDataFace::findIntersectOp::operator()().
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().
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().
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().
|
inline |
Return asymetrically extended bounding box, with guaranteed.
minimum width of s*mag(span) in any direction
Definition at line 315 of file treeBoundBoxI.H.
References Foam::constant::physicoChemical::b, delta, Foam::mag(), Foam::max(), Foam::sqrt(), VectorSpace< Vector< scalar >, scalar, 3 >::uniform(), Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().
Referenced by meshSearch::boundaryTree(), meshSearch::cellTree(), triSurfaceMesh::edgeTree(), patchProbes::findElements(), NamedEnum< compressibleField, 8 >::names(), listPlusEqOp< T >::operator()(), refinementFeatures::regionEdgeTrees(), triSurfaceSearch::tree(), and triSurfaceRegionSearch::treeByRegion().
void writeOBJ | ( | const fileName & | fName | ) | const |
Definition at line 619 of file treeBoundBox.C.
References Foam::e, treeBoundBox::edges, Foam::endl(), forAll, Foam::Info, OFstream::name(), Foam::nl, treeBoundBox::points(), Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().
|
friend |
|
friend |
|
friend |
|
friend |
|
static |
The great value used for greatBox and invertedBox.
Definition at line 102 of file treeBoundBox.H.
|
static |
As per boundBox::greatBox, but with great instead of vGreat.
Definition at line 105 of file treeBoundBox.H.
|
static |
As per boundBox::invertedBox, but with great instead of vGreat.
Definition at line 108 of file treeBoundBox.H.
|
static |
Face to point addressing.
Definition at line 175 of file treeBoundBox.H.
|
static |
Edge to point addressing.
Definition at line 178 of file treeBoundBox.H.
Referenced by treeBoundBox::writeOBJ().
|
static |
Per face the unit normal.
Definition at line 181 of file treeBoundBox.H.
Referenced by searchableBox::getNormal().