Loading...
Searching...
No Matches
treeBoundBox.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2017-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::treeBoundBox
29
30Description
31 Standard boundBox with extra functionality for use in octree.
32
33 Numbering of corner points is according to octant numbering as shown below.
34 For vertex numbering in the sequence 0 to 7:
35 - faces 0 and 1 are x-min and x-max, respectively;
36 - faces 2 and 3 are y-min and y-max, respectively;
37 - faces 4 and 5 are z-min and z-max, respectively.
38 .
39
40 \verbatim
41 Octant vertex ordering | Hex cell ordering
42 (treeBoundBox) | (boundBox, blockMesh, hexCell)
43 |
44 6 ---- 7 | 7 ---- 6
45 f5 |\ :\ f3 | f5 |\ :\ f3
46 | | 4 ---- 5 \ | | | 4 ---- 5 \
47 | 2.|....3 | \ | | 3.|....2 | \
48 | \| \| f2 | | \| \| f2
49 f4 0 ---- 1 | f4 0 ---- 1
50 Y Z |
51 \ | f0 ------ f1 | f0 ------ f1
52 \| |
53 o--- X |
54 \endverbatim
55
56Note
57 When a bounding box is created without any points, it creates an inverted
58 bounding box. Points can be added later and the bounding box will grow to
59 include them.
60
61SeeAlso
62 Foam::boundBox
63
64SourceFiles
65 treeBoundBoxI.H
66 treeBoundBox.cxx
67 treeBoundBox.txx
68
69\*---------------------------------------------------------------------------*/
70
71#ifndef Foam_treeBoundBox_H
72#define Foam_treeBoundBox_H
73
74#include "boundBox.H"
75
76// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
77
78namespace Foam
81// Forward Declarations
82class treeBoundBox;
86
87// List types
89
91/*---------------------------------------------------------------------------*\
92 Class treeBoundBox Declaration
93\*---------------------------------------------------------------------------*/
94
95class treeBoundBox
96:
97 public boundBox
98{
99public:
100
101 // Enumerations
102
103 //- Bits used for octant/point directional encoding.
104 // Every octant/corner point is the combination of three directions.
105 // Defined as (1 << vector::components).
109 TOPHALF = 2,
110 FRONTHALF = 4
111 };
112
113 //- Face codes. Identical order and meaning as per hex cellmodel
114 //- and boundBox, but have a different point ordering.
115 // Defined as (2 * vector::components + positive).
116 enum faceId : direction
118 LEFT = 0,
119 RIGHT = 1,
120 BOTTOM = 2,
121 TOP = 3,
122 BACK = 4,
123 FRONT = 5
124 };
125
126 //- Bits used for face encoding.
127 // Defined as (1 << (2 * vector::components + positive)).
128 // For example, the positive Z-face:
129 // (1 << (2 * vector::Z + 1))
130 enum faceBit : direction
131 {
132 NOFACE = 0,
133 LEFTBIT = (1 << LEFT),
134 RIGHTBIT = (1 << RIGHT),
135 BOTTOMBIT = (1 << BOTTOM),
136 TOPBIT = (1 << TOP),
137 BACKBIT = (1 << BACK),
138 FRONTBIT = (1 << FRONT)
139 };
141 //- Edges codes.
142 // E01 = edge between 0 and 1.
143 enum edgeId
144 {
145 E01 = 0,
146 E13 = 1,
147 E23 = 2,
148 E02 = 3,
149
150 E45 = 4,
151 E57 = 5,
152 E67 = 6,
153 E46 = 7,
155 E04 = 8,
156 E15 = 9,
157 E37 = 10,
158 E26 = 11
159 };
161
162 // Static Data Members
164 //- Face to point addressing, using octant corner points.
165 static const faceList faces;
166
167 //- Edge to point addressing, using octant corner points.
168 static const edgeList edges;
169
170
171 // Static Methods
172
173 //- The null treeBoundBox is the same as an inverted box
174 static const treeBoundBox& null() noexcept
175 {
176 return
177 reinterpret_cast<const treeBoundBox&>(boundBox::invertedBox);
178 }
180
181 // Standard (Generated) Methods
182
183 //- Default construct: an inverted bounding box
184 treeBoundBox() = default;
185
186 //- Copy construct
187 treeBoundBox(const treeBoundBox&) = default;
188
189 //- Copy assignment
190 treeBoundBox& operator=(const treeBoundBox&) = default;
191
192 //- Copy construct from a boundBox
193 explicit treeBoundBox(const boundBox& bb) : boundBox(bb) {}
194
195 //- Copy assignment from a boundBox
196 void operator=(const boundBox& bb)
197 {
198 static_cast<boundBox&>(*this) = bb;
200
201
202 // Constructors
203
204 //- Construct a 0/1 unit bounding box
205 inline explicit treeBoundBox(Foam::zero_one);
206
207 //- Construct a bounding box containing a single initial point
208 inline explicit treeBoundBox(const point& p);
210 //- Construct from bound box min/max points
211 inline treeBoundBox(const point& min, const point& max);
212
215
216 //- Construct from bound box min/max points
217 inline explicit treeBoundBox(const Pair<point>& bb);
218
219 //- Construct as the bounding box of the given pointField.
220 // Local processor domain only (no reduce as in boundBox)
221 explicit treeBoundBox(const UList<point>& points);
222
223 //- Construct as subset of points
224 // Local processor domain only (no reduce as in boundBox)
225 treeBoundBox(const UList<point>& points, const labelUList& indices);
226
227 //- Construct as subset of points
228 // The indices could be from edge/triFace etc.
229 // Local processor domain only (no reduce as in boundBox)
230 template<unsigned N>
232 (
233 const UList<point>& points,
234 const FixedList<label, N>& indices
235 );
236
237 //- Construct from Istream
238 inline explicit treeBoundBox(Istream& is);
239
240
241 // Member Functions
242
243 // Access
244
245 //- Vertex coordinates. In octant coding.
246 tmp<pointField> points() const;
247
248
249 // Octant handling
250
251 //- Corner point of given octant
252 inline point corner(const direction octant) const;
253
254 //- Sub-box of given octant. Midpoint calculated.
255 treeBoundBox subBbox(const direction octant) const;
256
257 //- Sub-box given by octant number. Midpoint provided.
258 treeBoundBox subBbox(const point& mid, const direction) const;
259
260 //- Sub-box half for given face
261 treeBoundBox subHalf(const direction whichFace) const;
263 //- Sub-box half for given face with prescribed mid point value.
264 //- Eg, subHalf(scalar, LEFT)
266 (
267 const scalar mid,
268 const direction whichFace
269 ) const;
270
271 //- Returns octant number given point and the calculated midpoint.
272 inline direction subOctant
273 (
274 const point& pt
275 ) const;
276
277 //- Returns octant number given point and midpoint.
278 static inline direction subOctant
279 (
280 const point& mid,
281 const point& pt
282 );
283
284 //- Returns octant number given point and the calculated midpoint.
285 // onEdge set if the point is on edge of subOctant
286 inline direction subOctant
287 (
288 const point& pt,
289 bool& onEdge
290 ) const;
291
292 //- Returns octant number given point and midpoint.
293 // onEdge set if the point is on edge of subOctant
294 static inline direction subOctant
295 (
296 const point& mid,
297 const point& pt,
298 bool& onEdge
299 );
300
301 //- Returns octant number given intersection and midpoint.
302 // onEdge set if the point is on edge of subOctant
303 // If onEdge, the direction vector determines which octant to use
304 // (acc. to which octant the point would be if it were moved
305 // along dir)
306 static inline direction subOctant
307 (
308 const point& mid,
309 const vector& dir,
310 const point& pt,
311 bool& onEdge
312 );
314 //- Calculates optimal order to look for nearest to point.
315 // First will be the octant containing the point,
316 // second the octant with boundary nearest to the point etc.
317 inline void searchOrder
318 (
319 const point& pt,
320 FixedList<direction, 8>& octantOrder
321 ) const;
322
323 //- Return optimal search order to look for nearest to point.
324 inline FixedList<direction, 8> searchOrder(const point& pt) const;
325
326
327 // Check
328
329 //- Overlaps with other bounding box, sphere etc?
330 using boundBox::overlaps;
331
332 //- Does sub-octant overlap/touch boundingBox?
333 bool subOverlaps
334 (
335 const direction octant,
336 const boundBox& bb
337 ) const;
338
339 //- Does sub-octant overlap boundingSphere (centre + sqr(radius))?
340 inline bool subOverlaps
341 (
342 const direction octant,
343 const point& centre,
344 const scalar radiusSqr
345 ) const;
346
347 //- intersects other bounding box, sphere etc?
349
350 //- Intersects segment; set point to intersection position and face,
351 // return true if intersection found.
352 // (pt argument used during calculation even if not intersecting).
353 // Calculates intersections from outside supplied vector
354 // (overallStart, overallVec). This is so when
355 // e.g. tracking through lots of consecutive boxes
356 // (typical octree) we're not accumulating truncation errors. Set
357 // to start, (end-start) if not used.
358 bool intersects
359 (
360 const point& overallStart,
361 const vector& overallVec,
362 const point& start,
363 const point& end,
364 point& pt,
365 direction& ptBits
366 ) const;
367
368 //- Like above but does not return faces point is on
369 bool intersects
370 (
371 const point& start,
372 const point& end,
373 point& pt
374 ) const;
375
376 //- Like above but does not return faces point is on
377 inline bool intersects(const linePointRef& ln, point& pt) const;
378
379 //- Contains point or other bounding box?
380 using boundBox::contains;
381
382 //- Contains point (inside or on edge) and moving in direction
383 // dir would cause it to go inside.
384 bool contains(const vector& dir, const point&) const;
385
386 //- Code position of point on bounding box faces
387 direction faceBits(const point& pt) const;
388
389 //- Position of point relative to bounding box
390 direction posBits(const point& pt) const;
391
392 //- Calculate nearest and furthest (to point) vertex coords of
393 // bounding box
394 void calcExtremities
395 (
396 const point& pt,
397 point& nearest,
398 point& furthest
399 ) const;
400
401 //- Returns distance point to furthest away corner.
402 scalar maxDist(const point& pt) const;
403
404 //- Compare distance to point with other bounding box
405 // return:
406 // -1 : all vertices of my bounding box are nearer than any of
407 // other
408 // +1 : all vertices of my bounding box are further away than
409 // any of other
410 // 0 : none of the above.
411 label distanceCmp(const point& pt, const treeBoundBox& other) const;
412
413 //- Return slightly wider bounding box
414 // Extends all dimensions with s*span*Random::sample01<scalar>()
415 // and guarantees in any direction s*mag(span) minimum width
416 inline treeBoundBox extend(Random& rndGen, const scalar s) const;
417
418 //- As per two parameter version but with additional
419 //- grow() by given amount in each dimension
420 inline treeBoundBox extend
421 (
422 Random& rndGen,
423 const scalar s,
424 const scalar delta
425 ) const;
426
427
428 // IOstream Operators
429
430 friend Istream& operator>>(Istream& is, treeBoundBox& bb);
431 friend Ostream& operator<<(Ostream& os, const treeBoundBox& bb);
432
433
434 // Housekeeping
435
436 //- Typical dimension length,height,width. Identical to avgDim()
437 scalar typDim() const { return avgDim(); }
438};
439
440
441// * * * * * * * * * * * * * * * * * Traits * * * * * * * * * * * * * * * * //
442
443//- Contiguous data for treeBoundBox
444template<> struct is_contiguous<treeBoundBox> : is_contiguous<boundBox> {};
445
446//- Contiguous scalar data for treeBoundBox
447template<> struct is_contiguous_scalar<treeBoundBox>
448:
449 is_contiguous_scalar<boundBox>
450{};
451
452
453// * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
455inline bool operator==(const treeBoundBox& a, const treeBoundBox& b);
456inline bool operator!=(const treeBoundBox& a, const treeBoundBox& b);
457
458
459// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
460
461} // End namespace Foam
462
463#include "treeBoundBoxI.H"
464
465#ifdef NoRepository
466 #include "treeBoundBox.txx"
467#endif
468
469// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
470
471#endif
472
473// ************************************************************************* //
scalar delta
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
An ordered pair of two objects of type <T> with first() and second() elements.
Definition Pair.H:66
Random number generator.
Definition Random.H:56
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
const point & max() const noexcept
Maximum describing the bounding box.
Definition boundBoxI.H:168
bool overlaps(const boundBox &bb) const
Overlaps/touches boundingBox?
Definition boundBoxI.H:439
static const boundBox invertedBox
A large inverted boundBox: min/max == +/- ROOTVGREAT.
Definition boundBox.H:131
const point & min() const noexcept
Minimum describing the bounding box.
Definition boundBoxI.H:162
bool contains(const point &pt) const
Contains point? (inside or on edge).
Definition boundBoxI.H:455
point nearest(const point &p) const
Return the nearest point on the boundBox to the supplied point.
scalar avgDim() const
Average length/height/width dimension.
Definition boundBoxI.H:228
boundBox()
Default construct: an inverted bounding box.
Definition boundBoxI.H:101
bool intersects(const plane &pln) const
Does plane intersect this bounding box.
point centre() const
The centre (midpoint) of the bounding box.
Definition boundBoxI.H:186
A class for managing temporary objects.
Definition tmp.H:75
Standard boundBox with extra functionality for use in octree.
void calcExtremities(const point &pt, point &nearest, point &furthest) const
Calculate nearest and furthest (to point) vertex coords of.
treeBoundBox(const boundBox &bb)
Copy construct from a boundBox.
scalar typDim() const
Typical dimension length,height,width. Identical to avgDim().
faceBit
Bits used for face encoding.
@ FRONTBIT
32: z-max face
@ TOPBIT
8: y-max face
@ LEFTBIT
1: x-min face
@ BACKBIT
16: z-min face
@ RIGHTBIT
2: x-max face
@ BOTTOMBIT
4: y-min face
treeBoundBox extend(Random &rndGen, const scalar s) const
Return slightly wider bounding box.
static const edgeList edges
Edge to point addressing, using octant corner points.
static const faceList faces
Face to point addressing, using octant corner points.
void searchOrder(const point &pt, FixedList< direction, 8 > &octantOrder) const
Calculates optimal order to look for nearest to point.
tmp< pointField > points() const
Vertex coordinates. In octant coding.
treeBoundBox(const UList< point > &points, const labelUList &indices)
Construct as subset of points.
friend Istream & operator>>(Istream &is, treeBoundBox &bb)
treeBoundBox(const UList< point > &points, const FixedList< label, N > &indices)
Construct as subset of points.
label distanceCmp(const point &pt, const treeBoundBox &other) const
Compare distance to point with other bounding box.
friend Ostream & operator<<(Ostream &os, const treeBoundBox &bb)
scalar maxDist(const point &pt) const
Returns distance point to furthest away corner.
treeBoundBox subBbox(const direction octant) const
Sub-box of given octant. Midpoint calculated.
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,.
point corner(const direction octant) const
Corner point of given octant.
treeBoundBox()=default
Default construct: an inverted bounding box.
void operator=(const boundBox &bb)
Copy assignment from a boundBox.
bool subOverlaps(const direction octant, const boundBox &bb) const
Does sub-octant overlap/touch boundingBox?
bool contains(const vector &dir, const point &) const
Contains point (inside or on edge) and moving in direction.
static const treeBoundBox & null() noexcept
The null treeBoundBox is the same as an inverted box.
faceId
Face codes. Identical order and meaning as per hex cellmodel and boundBox, but have a different point...
@ TOP
3: y-max face
@ FRONT
5: z-max face
@ BOTTOM
2: y-min face
@ BACK
4: z-min face
@ LEFT
0: x-min face
@ RIGHT
1: x-max face
direction subOctant(const point &pt) const
Returns octant number given point and the calculated midpoint.
treeBoundBox subHalf(const direction whichFace) const
Sub-box half for given face.
octantBit
Bits used for octant/point directional encoding.
@ TOPHALF
2: positive y-direction
@ RIGHTHALF
1: positive x-direction
@ FRONTHALF
4: positive z-direction
direction faceBits(const point &pt) const
Code position of point on bounding box faces.
treeBoundBox(const treeBoundBox &)=default
Copy construct.
treeBoundBox subBbox(const point &mid, const direction) const
Sub-box given by octant number. Midpoint provided.
bool intersects(const point &start, const point &end, point &pt) const
Like above but does not return faces point is on.
direction posBits(const point &pt) const
Position of point relative to bounding box.
treeBoundBox(const UList< point > &points)
Construct as the bounding box of the given pointField.
treeBoundBox subHalf(const scalar mid, const direction whichFace) const
Sub-box half for given face with prescribed mid point value. Eg, subHalf(scalar, LEFT).
treeBoundBox & operator=(const treeBoundBox &)=default
Copy assignment.
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
Definition pTraits.H:51
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
label faceId(-1)
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
bool operator!=(const eddy &a, const eddy &b)
Definition eddy.H:297
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< treeBoundBox > treeBoundBoxList
A List of treeBoundBox.
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces).
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Istream & operator>>(Istream &, directionInfo &)
line< point, const point & > linePointRef
A line using referred points.
Definition line.H:66
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
uint8_t direction
Definition direction.H:49
vector point
Point is a vector.
Definition point.H:37
const direction noexcept
Definition scalarImpl.H:265
Vector< scalar > vector
Definition vector.H:57
UList< label > labelUList
A UList of labels.
Definition UList.H:75
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition POSIX.C:1239
volScalarField & b
A template class to specify if a data type is composed solely of Foam::scalar elements.
Definition contiguous.H:87
A template class to specify that a data type can be considered as being contiguous in memory.
Definition contiguous.H:70
Random rndGen