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-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
27 Class
28  Foam::treeBoundBox
29 
30 Description
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 
56 Note
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 
61 SeeAlso
62  Foam::boundBox
63 
64 SourceFiles
65  treeBoundBoxI.H
66  treeBoundBox.C
67  treeBoundBoxTemplates.C
68 
69 \*---------------------------------------------------------------------------*/
70 
71 #ifndef Foam_treeBoundBox_H
72 #define Foam_treeBoundBox_H
73 
74 #include "boundBox.H"
75 
76 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
77 
78 namespace Foam
79 {
80 
81 // Forward Declarations
82 class treeBoundBox;
83 
86 
87 // List types
89 
90 
91 /*---------------------------------------------------------------------------*\
92  Class treeBoundBox Declaration
93 \*---------------------------------------------------------------------------*/
94 
95 class treeBoundBox
96 :
97  public boundBox
98 {
99 public:
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).
107  {
108  RIGHTHALF = 1,
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
117  {
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;
199  }
200 
201 
202  // Constructors
203 
204  //- Construct a bounding box containing a single initial point
205  inline explicit treeBoundBox(const point& p);
206 
207  //- Construct from bound box min/max points
208  inline treeBoundBox(const point& min, const point& max);
209 
212 
213  //- Construct from bound box min/max points
214  inline explicit treeBoundBox(const Pair<point>& bb);
215 
216  //- Construct as the bounding box of the given pointField.
217  // Local processor domain only (no reduce as in boundBox)
218  explicit treeBoundBox(const UList<point>& points);
220  //- Construct as subset of points
221  // Local processor domain only (no reduce as in boundBox)
222  treeBoundBox(const UList<point>& points, const labelUList& indices);
223 
224  //- Construct as subset of points
225  // The indices could be from edge/triFace etc.
226  // Local processor domain only (no reduce as in boundBox)
227  template<unsigned N>
229  (
230  const UList<point>& points,
231  const FixedList<label, N>& indices
232  );
233 
234  //- Construct from Istream
235  inline explicit treeBoundBox(Istream& is);
236 
237 
238  // Member Functions
239 
240  // Access
241 
242  //- Vertex coordinates. In octant coding.
243  tmp<pointField> points() const;
244 
245 
246  // Octant handling
247 
248  //- Corner point of given octant
249  inline point corner(const direction octant) const;
250 
251  //- Sub-box of given octant. Midpoint calculated.
252  treeBoundBox subBbox(const direction octant) const;
253 
254  //- Sub-box given by octant number. Midpoint provided.
255  treeBoundBox subBbox(const point& mid, const direction) const;
256 
257  //- Sub-box half for given face
258  treeBoundBox subHalf(const direction whichFace) const;
259 
260  //- Sub-box half for given face with prescribed mid point value.
261  //- Eg, subHalf(scalar, LEFT)
263  (
264  const scalar mid,
265  const direction whichFace
266  ) const;
267 
268  //- Returns octant number given point and the calculated midpoint.
269  inline direction subOctant
270  (
271  const point& pt
272  ) const;
273 
274  //- Returns octant number given point and midpoint.
275  static inline direction subOctant
276  (
277  const point& mid,
278  const point& pt
279  );
280 
281  //- Returns octant number given point and the calculated midpoint.
282  // onEdge set if the point is on edge of subOctant
283  inline direction subOctant
284  (
285  const point& pt,
286  bool& onEdge
287  ) const;
288 
289  //- Returns octant number given point and midpoint.
290  // onEdge set if the point is on edge of subOctant
291  static inline direction subOctant
292  (
293  const point& mid,
294  const point& pt,
295  bool& onEdge
296  );
297 
298  //- Returns octant number given intersection and midpoint.
299  // onEdge set if the point is on edge of subOctant
300  // If onEdge, the direction vector determines which octant to use
301  // (acc. to which octant the point would be if it were moved
302  // along dir)
303  static inline direction subOctant
304  (
305  const point& mid,
306  const vector& dir,
307  const point& pt,
308  bool& onEdge
309  );
310 
311  //- Calculates optimal order to look for nearest to point.
312  // First will be the octant containing the point,
313  // second the octant with boundary nearest to the point etc.
314  inline void searchOrder
315  (
316  const point& pt,
317  FixedList<direction, 8>& octantOrder
318  ) const;
319 
320  //- Return optimal search order to look for nearest to point.
321  inline FixedList<direction, 8> searchOrder(const point& pt) const;
322 
323 
324  // Check
325 
326  //- Overlaps with other bounding box, sphere etc?
327  using boundBox::overlaps;
328 
329  //- Does sub-octant overlap/touch boundingBox?
330  bool subOverlaps
331  (
332  const direction octant,
333  const boundBox& bb
334  ) const;
335 
336  //- Does sub-octant overlap boundingSphere (centre + sqr(radius))?
337  inline bool subOverlaps
338  (
339  const direction octant,
340  const point& centre,
341  const scalar radiusSqr
342  ) const;
343 
344  //- intersects other bounding box, sphere etc?
345  using boundBox::intersects;
346 
347  //- Intersects segment; set point to intersection position and face,
348  // return true if intersection found.
349  // (pt argument used during calculation even if not intersecting).
350  // Calculates intersections from outside supplied vector
351  // (overallStart, overallVec). This is so when
352  // e.g. tracking through lots of consecutive boxes
353  // (typical octree) we're not accumulating truncation errors. Set
354  // to start, (end-start) if not used.
355  bool intersects
356  (
357  const point& overallStart,
358  const vector& overallVec,
359  const point& start,
360  const point& end,
361  point& pt,
362  direction& ptBits
363  ) const;
364 
365  //- Like above but does not return faces point is on
366  bool intersects
367  (
368  const point& start,
369  const point& end,
370  point& pt
371  ) const;
372 
373  //- Like above but does not return faces point is on
374  inline bool intersects(const linePointRef& ln, point& pt) const;
375 
376  //- Contains point or other bounding box?
377  using boundBox::contains;
378 
379  //- Contains point (inside or on edge) and moving in direction
380  // dir would cause it to go inside.
381  bool contains(const vector& dir, const point&) const;
382 
383  //- Code position of point on bounding box faces
384  direction faceBits(const point& pt) const;
385 
386  //- Position of point relative to bounding box
387  direction posBits(const point& pt) const;
388 
389  //- Calculate nearest and furthest (to point) vertex coords of
390  // bounding box
391  void calcExtremities
392  (
393  const point& pt,
394  point& nearest,
395  point& furthest
396  ) const;
397 
398  //- Returns distance point to furthest away corner.
399  scalar maxDist(const point& pt) const;
400 
401  //- Compare distance to point with other bounding box
402  // return:
403  // -1 : all vertices of my bounding box are nearer than any of
404  // other
405  // +1 : all vertices of my bounding box are further away than
406  // any of other
407  // 0 : none of the above.
408  label distanceCmp(const point& pt, const treeBoundBox& other) const;
409 
410  //- Return slightly wider bounding box
411  // Extends all dimensions with s*span*Random::sample01<scalar>()
412  // and guarantees in any direction s*mag(span) minimum width
413  inline treeBoundBox extend(Random& rndGen, const scalar s) const;
414 
415  //- As per two parameter version but with additional
416  //- grow() by given amount in each dimension
417  inline treeBoundBox extend
418  (
419  Random& rndGen,
420  const scalar s,
421  const scalar delta
422  ) const;
423 
424 
425  // IOstream Operators
426 
427  friend Istream& operator>>(Istream& is, treeBoundBox& bb);
428  friend Ostream& operator<<(Ostream& os, const treeBoundBox& bb);
429 
430 
431  // Housekeeping
432 
433  //- Typical dimension length,height,width. Identical to avgDim()
434  scalar typDim() const { return avgDim(); }
435 };
436 
437 
438 // * * * * * * * * * * * * * * * * * Traits * * * * * * * * * * * * * * * * //
439 
440 //- Contiguous data for treeBoundBox
441 template<> struct is_contiguous<treeBoundBox> : is_contiguous<boundBox> {};
442 
443 //- Contiguous scalar data for treeBoundBox
444 template<> struct is_contiguous_scalar<treeBoundBox>
445 :
446  is_contiguous_scalar<boundBox>
447 {};
448 
449 
450 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
451 
452 inline bool operator==(const treeBoundBox& a, const treeBoundBox& b);
453 inline bool operator!=(const treeBoundBox& a, const treeBoundBox& b);
454 
455 
456 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
457 
458 } // End namespace Foam
459 
460 #include "treeBoundBoxI.H"
461 
462 #ifdef NoRepository
463  #include "treeBoundBoxTemplates.C"
464 #endif
465 
466 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
467 
468 #endif
469 
470 // ************************************************************************* //
edgeId
Edges codes.
Definition: treeBoundBox.H:150
scalar delta
direction faceBits(const point &pt) const
Code position of point on bounding box faces.
Definition: treeBoundBox.C:451
direction subOctant(const point &pt) const
Returns octant number given point and the calculated midpoint.
Definition: treeBoundBoxI.H:65
treeBoundBox subHalf(const direction whichFace) const
Sub-box half for given face.
Definition: treeBoundBox.C:185
uint8_t direction
Definition: direction.H:46
A line primitive.
Definition: line.H:52
static const treeBoundBox & null() noexcept
The null treeBoundBox is the same as an inverted box.
Definition: treeBoundBox.H:187
bool subOverlaps(const direction octant, const boundBox &bb) const
Does sub-octant overlap/touch boundingBox?
Definition: treeBoundBox.C:211
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:101
bool contains(const point &pt) const
Contains point? (inside or on edge)
Definition: boundBoxI.H:441
treeBoundBox & operator=(const treeBoundBox &)=default
Copy assignment.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Random rndGen
Definition: createFields.H:23
octantBit
Bits used for octant/point directional encoding.
Definition: treeBoundBox.H:104
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
void searchOrder(const point &pt, FixedList< direction, 8 > &octantOrder) const
Calculates optimal order to look for nearest to point.
List< treeBoundBox > treeBoundBoxList
A List of treeBoundBox.
Definition: treeBoundBox.H:83
label distanceCmp(const point &pt, const treeBoundBox &other) const
Compare distance to point with other bounding box.
Definition: treeBoundBox.C:558
const point & min() const noexcept
Minimum describing the bounding box.
Definition: boundBoxI.H:155
bool intersects(const plane &pln) const
Does plane intersect this bounding box.
Definition: boundBox.C:193
static const faceList faces
Face to point addressing, using octant corner points.
Definition: treeBoundBox.H:174
1: positive x-direction
Definition: treeBoundBox.H:106
scalar avgDim() const
Average length/height/width dimension.
Definition: boundBoxI.H:221
point corner(const direction octant) const
Corner point of given octant.
Definition: treeBoundBoxI.H:53
direction posBits(const point &pt) const
Position of point relative to bounding box.
Definition: treeBoundBox.C:486
const point & max() const noexcept
Maximum describing the bounding box.
Definition: boundBoxI.H:161
friend Ostream & operator<<(Ostream &os, const treeBoundBox &bb)
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,.
Definition: treeBoundBox.C:255
point nearest(const point &p) const
Return the nearest point on the boundBox to the supplied point.
Definition: boundBox.C:404
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Istream & operator>>(Istream &, directionInfo &)
static const edgeList edges
Edge to point addressing, using octant corner points.
Definition: treeBoundBox.H:179
treeBoundBox subBbox(const direction octant) const
Sub-box of given octant. Midpoint calculated.
Definition: treeBoundBox.C:204
treeBoundBox()=default
Default construct: an inverted bounding box.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:99
Random number generator.
Definition: Random.H:55
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: POSIX.C:1190
void calcExtremities(const point &pt, point &nearest, point &furthest) const
Calculate nearest and furthest (to point) vertex coords of.
Definition: treeBoundBox.C:522
bool contains(const vector &dir, const point &) const
Contains point (inside or on edge) and moving in direction.
Definition: treeBoundBox.C:413
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:55
const direction noexcept
Definition: Scalar.H:258
scalar maxDist(const point &pt) const
Returns distance point to furthest away corner.
Definition: treeBoundBox.C:548
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:193
OBJstream os(runTime.globalPath()/outputName)
4: positive z-direction
Definition: treeBoundBox.H:108
faceBit
Bits used for face encoding.
Definition: treeBoundBox.H:134
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:76
static const boundBox invertedBox
A large inverted boundBox: min/max == +/- ROOTVGREAT.
Definition: boundBox.H:124
vector point
Point is a vector.
Definition: point.H:37
friend Istream & operator>>(Istream &is, treeBoundBox &bb)
tmp< pointField > points() const
Vertex coordinates. In octant coding.
Definition: treeBoundBox.C:87
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:90
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:299
faceId
Face codes. Identical order and meaning as per hex cellmodel and boundBox, but have a different point...
Definition: treeBoundBox.H:117
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
point centre() const
The centre (midpoint) of the bounding box.
Definition: boundBoxI.H:179
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))
treeBoundBox(const boundBox &bb)
Copy construct from a boundBox.
Definition: treeBoundBox.H:214
bool overlaps(const boundBox &bb) const
Overlaps/touches boundingBox?
Definition: boundBoxI.H:425
scalar typDim() const
Typical dimension length,height,width. Identical to avgDim()
Definition: treeBoundBox.H:543
treeBoundBox extend(Random &rndGen, const scalar s) const
Return slightly wider bounding box.
2: positive y-direction
Definition: treeBoundBox.H:107
Namespace for OpenFOAM.