treeBoundBox.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::treeBoundBox
26 
27 Description
28  Standard boundBox + extra functionality for use in octree.
29 
30  Numbering of corner points is according to octant numbering.
31 
32  On the back plane (z=0):
33 
34  \verbatim
35  Y
36  ^
37  |
38  +--------+
39  |2 3|
40  | |
41  | |
42  | |
43  |0 1|
44  +--------+->X
45  \endverbatim
46 
47  For the front plane add 4 to the point labels.
48 
49 
50 SourceFiles
51  treeBoundBoxI.H
52  treeBoundBox.C
53  treeBoundBoxTemplates.C
54 
55 \*---------------------------------------------------------------------------*/
56 
57 #ifndef treeBoundBox_H
58 #define treeBoundBox_H
59 
60 #include "boundBox.H"
61 #include "direction.H"
62 #include "pointField.H"
63 #include "faceList.H"
64 
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66 
67 namespace Foam
68 {
69 
70 class randomGenerator;
71 
72 
73 // Forward declaration of friend functions and operators
74 
75 class treeBoundBox;
76 
77 bool operator==(const treeBoundBox&, const treeBoundBox&);
78 bool operator!=(const treeBoundBox&, const treeBoundBox&);
79 
80 Istream& operator>>(Istream& is, treeBoundBox&);
81 Ostream& operator<<(Ostream& os, const treeBoundBox&);
82 
83 
84 /*---------------------------------------------------------------------------*\
85  Class treeBoundBox Declaration
86 \*---------------------------------------------------------------------------*/
87 
88 class treeBoundBox
89 :
90  public boundBox
91 {
92  //- To initialise edges.
93  static edgeList calcEdges(const label[12][2]);
94 
95  //- To initialise faceNormals.
96  static FixedList<vector, 6> calcFaceNormals();
97 
98 public:
99 
100  // Static Data Members
101 
102  //- The great value used for greatBox and invertedBox
103  static const scalar great;
104 
105  //- As per boundBox::greatBox, but with great instead of vGreat
106  static const treeBoundBox greatBox;
107 
108  //- As per boundBox::invertedBox, but with great instead of vGreat
109  static const treeBoundBox invertedBox;
110 
111  //- Bits used for octant/point coding.
112  // Every octant/corner point is the combination of three faces.
113  struct octantBit
114  {
115  enum
116  {
117  rightHalf = 0x1 << 0,
118  topHalf = 0x1 << 1,
119  frontHalf = 0x1 << 2
120  };
121  };
122 
123  //- Face codes
124  struct faceId
125  {
126  enum
127  {
128  left = 0,
129  right = 1,
130  bottom = 2,
131  top = 3,
132  back = 4,
133  front = 5
134  };
135  };
136 
137  //- Bits used for face coding
138  struct faceBit
139  {
140  enum
141  {
142  noFace = 0,
143  left = 0x1 << faceId::left, //1
144  right = 0x1 << faceId::right, //2
145  bottom = 0x1 << faceId::bottom, //4
146  top = 0x1 << faceId::top, //8
147  back = 0x1 << faceId::back, //16
148  front = 0x1 << faceId::front, //32
149  };
150  };
151 
152  //- Edges codes.
153  // e01 = edge between 0 and 1.
154  struct edgeId
155  {
156  enum
157  {
158  e01 = 0,
159  e13 = 1,
160  e23 = 2,
161  e02 = 3,
162 
163  e45 = 4,
164  e57 = 5,
165  e67 = 6,
166  e46 = 7,
167 
168  e04 = 8,
169  e15 = 9,
170  e37 = 10,
171  e26 = 11
172  };
173  };
174 
175  //- Face to point addressing
176  static const faceList faces;
177 
178  //- Edge to point addressing
179  static const edgeList edges;
180 
181  //- Per face the unit normal
182  static const FixedList<vector, 6> faceNormals;
183 
184 
185  // Constructors
186 
187  //- Construct null setting points to zero
188  inline treeBoundBox();
189 
190  //- Construct from a boundBox
191  explicit inline treeBoundBox(const boundBox& bb);
192 
193  //- Construct from components
194  inline treeBoundBox(const point& min, const point& max);
195 
196  //- Construct as the bounding box of the given pointField.
197  // Local processor domain only (no reduce as in boundBox)
198  explicit treeBoundBox(const UList<point>&);
199 
200  //- Construct as subset of points
201  // Local processor domain only (no reduce as in boundBox)
202  treeBoundBox(const UList<point>&, const labelUList& indices);
203 
204  //- Construct as subset of points
205  // The indices could be from edge/triFace etc.
206  // Local processor domain only (no reduce as in boundBox)
207  template<unsigned Size>
209  (
210  const UList<point>&,
211  const FixedList<label, Size>& indices
212  );
213 
214  //- Construct from dictionary
215  inline treeBoundBox(const dictionary&);
216 
217  //- Construct from Istream
218  inline treeBoundBox(Istream&);
219 
220 
221  // Member Functions
222 
223  // Access
224 
225  //- Typical dimension length,height,width
226  inline scalar typDim() const;
227 
228  //- Vertex coordinates. In octant coding.
229  tmp<pointField> points() const;
230 
231 
232  // Check
233 
234  //- Corner point given octant
235  inline point corner(const direction) const;
236 
237  //- Sub box given by octant number. Midpoint calculated.
238  treeBoundBox subBbox(const direction) const;
239 
240  //- Sub box given by octant number. Midpoint provided.
241  treeBoundBox subBbox(const point& mid, const direction) const;
242 
243  //- Returns octant number given point and the calculated midpoint.
244  inline direction subOctant
245  (
246  const point& pt
247  ) const;
248 
249  //- Returns octant number given point and midpoint.
250  static inline direction subOctant
251  (
252  const point& mid,
253  const point& pt
254  );
255 
256  //- Returns octant number given point and the calculated midpoint.
257  // onEdge set if the point is on edge of subOctant
258  inline direction subOctant
259  (
260  const point& pt,
261  bool& onEdge
262  ) const;
263 
264  //- Returns octant number given point and midpoint.
265  // onEdge set if the point is on edge of subOctant
266  static inline direction subOctant
267  (
268  const point& mid,
269  const point& pt,
270  bool& onEdge
271  );
272 
273  //- Returns octant number given intersection and midpoint.
274  // onEdge set if the point is on edge of subOctant
275  // If onEdge, the direction vector determines which octant to use
276  // (acc. to which octant the point would be if it were moved
277  // along dir)
278  static inline direction subOctant
279  (
280  const point& mid,
281  const vector& dir,
282  const point& pt,
283  bool& onEdge
284  );
285 
286  //- Calculates optimal order to look for nearest to point.
287  // First will be the octant containing the point,
288  // second the octant with boundary nearest to the point etc.
289  inline void searchOrder
290  (
291  const point& pt,
292  FixedList<direction, 8>& octantOrder
293  ) const;
294 
295  //- Overlaps other bounding box?
296  using boundBox::overlaps;
297 
298  //- Intersects segment; set point to intersection position and face,
299  // return true if intersection found.
300  // (pt argument used during calculation even if not intersecting).
301  // Calculates intersections from outside supplied vector
302  // (overallStart, overallVec). This is so when
303  // e.g. tracking through lots of consecutive boxes
304  // (typical octree) we're not accumulating truncation errors. Set
305  // to start, (end-start) if not used.
306  bool intersects
307  (
308  const point& overallStart,
309  const vector& overallVec,
310  const point& start,
311  const point& end,
312  point& pt,
313  direction& ptBits
314  ) const;
315 
316  //- Like above but does not return faces point is on
317  bool intersects
318  (
319  const point& start,
320  const point& end,
321  point& pt
322  ) const;
323 
324  //- Contains point or other bounding box?
325  using boundBox::contains;
326 
327  //- Contains point (inside or on edge) and moving in direction
328  // dir would cause it to go inside.
329  bool contains(const vector& dir, const point&) const;
330 
331  //- Code position of point on bounding box faces
332  direction faceBits(const point&) const;
333 
334  //- Position of point relative to bounding box
335  direction posBits(const point&) const;
336 
337  //- Calculate nearest and furthest (to point) vertex coords of
338  // bounding box
339  void calcExtremities
340  (
341  const point& pt,
342  point& nearest,
343  point& furthest
344  ) const;
345 
346  //- Returns distance point to furthest away corner.
347  scalar maxDist(const point&) const;
348 
349  //- Compare distance to point with other bounding box
350  // return:
351  // -1 : all vertices of my bounding box are nearer than any of
352  // other
353  // +1 : all vertices of my bounding box are further away than
354  // any of other
355  // 0 : none of the above.
356  label distanceCmp(const point&, const treeBoundBox& other) const;
357 
358  //- Return asymmetrically extended bounding box, with guaranteed
359  // minimum width of s*mag(span) in any direction
360  inline treeBoundBox extend(const scalar s) const;
361 
362  // Write
363 
364  void writeOBJ(const fileName& fName) const;
365 
366 
367  // Friend Operators
368 
369  friend bool operator==(const treeBoundBox&, const treeBoundBox&);
370  friend bool operator!=(const treeBoundBox&, const treeBoundBox&);
371 
372 
373  // IOstream operator
374 
376  friend Ostream& operator<<(Ostream& os, const treeBoundBox&);
377 };
378 
379 
380 //- Data associated with treeBoundBox type are contiguous
381 template<>
382 inline bool contiguous<treeBoundBox>() {return contiguous<boundBox>();}
383 
384 
385 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
386 
387 } // End namespace Foam
388 
389 #include "treeBoundBoxI.H"
390 
391 #ifdef NoRepository
392  #include "treeBoundBoxTemplates.C"
393 #endif
394 
395 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
396 
397 #endif
398 
399 // ************************************************************************* //
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:78
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
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:74
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:59
point nearest(const point &) const
Return the nearest point on the boundBox to the supplied point.
Definition: boundBox.C:292
const point & min() const
Minimum point defining the bounding box.
Definition: boundBoxI.H:60
const point & max() const
Maximum point defining the bounding box.
Definition: boundBoxI.H:66
bool overlaps(const boundBox &) const
Overlaps/touches boundingBox?
Definition: boundBoxI.H:126
bool contains(const point &) const
Contains point? (inside or on edge)
Definition: boundBoxI.H:176
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
A class for handling file names.
Definition: fileName.H:82
A class for managing temporary objects.
Definition: tmp.H:55
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:90
void writeOBJ(const fileName &fName) const
Definition: treeBoundBox.C:612
void calcExtremities(const point &pt, point &nearest, point &furthest) const
Calculate nearest and furthest (to point) vertex coords of.
Definition: treeBoundBox.C:493
friend Istream & operator>>(Istream &is, treeBoundBox &)
friend Ostream & operator<<(Ostream &os, const treeBoundBox &)
scalar maxDist(const point &) const
Returns distance point to furthest away corner.
Definition: treeBoundBox.C:540
label distanceCmp(const point &, const treeBoundBox &other) const
Compare distance to point with other bounding box.
Definition: treeBoundBox.C:550
direction faceBits(const point &) const
Code position of point on bounding box faces.
Definition: treeBoundBox.C:425
void searchOrder(const point &pt, FixedList< direction, 8 > &octantOrder) const
Calculates optimal order to look for nearest to point.
treeBoundBox extend(const scalar s) const
Return asymmetrically extended bounding box, with guaranteed.
friend bool operator!=(const treeBoundBox &, const treeBoundBox &)
treeBoundBox subBbox(const direction) const
Sub box given by octant number. Midpoint calculated.
Definition: treeBoundBox.C:173
static const scalar great
The great value used for greatBox and invertedBox.
Definition: treeBoundBox.H:102
static const edgeList edges
Edge to point addressing.
Definition: treeBoundBox.H:178
point corner(const direction) const
Corner point given octant.
Definition: treeBoundBoxI.H:69
bool contains(const point &) const
Contains point or other bounding box?
Definition: boundBoxI.H:176
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:229
static const treeBoundBox invertedBox
As per boundBox::invertedBox, but with great instead of vGreat.
Definition: treeBoundBox.H:108
static const treeBoundBox greatBox
As per boundBox::greatBox, but with great instead of vGreat.
Definition: treeBoundBox.H:105
direction subOctant(const point &pt) const
Returns octant number given point and the calculated midpoint.
Definition: treeBoundBoxI.H:81
static const FixedList< vector, 6 > faceNormals
Per face the unit normal.
Definition: treeBoundBox.H:181
scalar typDim() const
Typical dimension length,height,width.
Definition: treeBoundBoxI.H:63
static const faceList faces
Face to point addressing.
Definition: treeBoundBox.H:175
friend bool operator==(const treeBoundBox &, const treeBoundBox &)
direction posBits(const point &) const
Position of point relative to bounding box.
Definition: treeBoundBox.C:458
treeBoundBox()
Construct null setting points to zero.
Definition: treeBoundBoxI.H:31
tmp< pointField > points() const
Vertex coordinates. In octant coding.
Definition: treeBoundBox.C:159
Direction is an 8-bit unsigned integer type used to represent the Cartesian directions etc.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Namespace for OpenFOAM.
bool operator!=(const particle &, const particle &)
Definition: particle.C:1210
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
bool contiguous< boundBox >()
Data associated with boundBox type are contiguous.
Definition: boundBox.H:254
bool contiguous< treeBoundBox >()
Data associated with treeBoundBox type are contiguous.
Definition: treeBoundBox.H:381
Istream & operator>>(Istream &, pistonPointEdgeData &)
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
uint8_t direction
Definition: direction.H:45
Bits used for face coding.
Definition: treeBoundBox.H:138
Bits used for octant/point coding.
Definition: treeBoundBox.H:113