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-2019 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 Random;
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,
163  e45 = 4,
164  e57 = 5,
165  e67 = 6,
166  e46 = 7,
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 
215  //- Construct from Istream
216  inline treeBoundBox(Istream&);
217 
218 
219  // Member Functions
220 
221  // Access
222 
223  //- Typical dimension length,height,width
224  inline scalar typDim() const;
225 
226  //- Vertex coordinates. In octant coding.
227  tmp<pointField> points() const;
228 
229 
230  // Check
231 
232  //- Corner point given octant
233  inline point corner(const direction) const;
234 
235  //- Sub box given by octant number. Midpoint calculated.
236  treeBoundBox subBbox(const direction) const;
237 
238  //- Sub box given by octant number. Midpoint provided.
239  treeBoundBox subBbox(const point& mid, const direction) const;
240 
241  //- Returns octant number given point and the calculated midpoint.
242  inline direction subOctant
243  (
244  const point& pt
245  ) const;
246 
247  //- Returns octant number given point and midpoint.
248  static inline direction subOctant
249  (
250  const point& mid,
251  const point& pt
252  );
253 
254  //- Returns octant number given point and the calculated midpoint.
255  // onEdge set if the point is on edge of subOctant
256  inline direction subOctant
257  (
258  const point& pt,
259  bool& onEdge
260  ) const;
261 
262  //- Returns octant number given point and midpoint.
263  // onEdge set if the point is on edge of subOctant
264  static inline direction subOctant
265  (
266  const point& mid,
267  const point& pt,
268  bool& onEdge
269  );
270 
271  //- Returns octant number given intersection and midpoint.
272  // onEdge set if the point is on edge of subOctant
273  // If onEdge, the direction vector determines which octant to use
274  // (acc. to which octant the point would be if it were moved
275  // along dir)
276  static inline direction subOctant
277  (
278  const point& mid,
279  const vector& dir,
280  const point& pt,
281  bool& onEdge
282  );
283 
284  //- Calculates optimal order to look for nearest to point.
285  // First will be the octant containing the point,
286  // second the octant with boundary nearest to the point etc.
287  inline void searchOrder
288  (
289  const point& pt,
290  FixedList<direction, 8>& octantOrder
291  ) const;
292 
293  //- Overlaps other bounding box?
294  using boundBox::overlaps;
295 
296  //- Intersects segment; set point to intersection position and face,
297  // return true if intersection found.
298  // (pt argument used during calculation even if not intersecting).
299  // Calculates intersections from outside supplied vector
300  // (overallStart, overallVec). This is so when
301  // e.g. tracking through lots of consecutive boxes
302  // (typical octree) we're not accumulating truncation errors. Set
303  // to start, (end-start) if not used.
304  bool intersects
305  (
306  const point& overallStart,
307  const vector& overallVec,
308  const point& start,
309  const point& end,
310  point& pt,
311  direction& ptBits
312  ) const;
313 
314  //- Like above but does not return faces point is on
315  bool intersects
316  (
317  const point& start,
318  const point& end,
319  point& pt
320  ) const;
321 
322  //- Contains point or other bounding box?
323  using boundBox::contains;
324 
325  //- Contains point (inside or on edge) and moving in direction
326  // dir would cause it to go inside.
327  bool contains(const vector& dir, const point&) const;
328 
329  //- Code position of point on bounding box faces
330  direction faceBits(const point&) const;
331 
332  //- Position of point relative to bounding box
333  direction posBits(const point&) const;
334 
335  //- Calculate nearest and furthest (to point) vertex coords of
336  // bounding box
337  void calcExtremities
338  (
339  const point& pt,
340  point& nearest,
341  point& furthest
342  ) const;
343 
344  //- Returns distance point to furthest away corner.
345  scalar maxDist(const point&) const;
346 
347  //- Compare distance to point with other bounding box
348  // return:
349  // -1 : all vertices of my bounding box are nearer than any of
350  // other
351  // +1 : all vertices of my bounding box are further away than
352  // any of other
353  // 0 : none of the above.
354  label distanceCmp(const point&, const treeBoundBox& other) const;
355 
356  //- Return asymetrically extended bounding box, with guaranteed
357  // minimum width of s*mag(span) in any direction
358  inline treeBoundBox extend(const scalar s) const;
359 
360  // Write
361 
362  void writeOBJ(const fileName& fName) const;
363 
364 
365  // Friend Operators
366 
367  friend bool operator==(const treeBoundBox&, const treeBoundBox&);
368  friend bool operator!=(const treeBoundBox&, const treeBoundBox&);
369 
370 
371  // IOstream operator
372 
373  friend Istream& operator>>(Istream& is, treeBoundBox&);
374  friend Ostream& operator<<(Ostream& os, const treeBoundBox&);
375 };
376 
377 
378 //- Data associated with treeBoundBox type are contiguous
379 template<>
380 inline bool contiguous<treeBoundBox>() {return contiguous<boundBox>();}
381 
382 
383 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
384 
385 } // End namespace Foam
386 
387 #include "treeBoundBoxI.H"
388 
389 #ifdef NoRepository
390  #include "treeBoundBoxTemplates.C"
391 #endif
392 
393 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
394 
395 #endif
396 
397 // ************************************************************************* //
static const FixedList< vector, 6 > faceNormals
Per face the unit normal.
Definition: treeBoundBox.H:181
bool overlaps(const boundBox &) const
Overlaps/touches boundingBox?
Definition: boundBoxI.H:120
scalar typDim() const
Typical dimension length,height,width.
Definition: treeBoundBoxI.H:57
direction subOctant(const point &pt) const
Returns octant number given point and the calculated midpoint.
Definition: treeBoundBoxI.H:75
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
A class for handling file names.
Definition: fileName.H:79
point nearest(const point &) const
Return the nearest point on the boundBox to the supplied point.
Definition: boundBox.C:303
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:54
Bits used for octant/point coding.
Definition: treeBoundBox.H:112
point corner(const direction) const
Corner point given octant.
Definition: treeBoundBoxI.H:63
static const treeBoundBox greatBox
As per boundBox::greatBox, but with great instead of vGreat.
Definition: treeBoundBox.H:105
treeBoundBox subBbox(const direction) const
Sub box given by octant number. Midpoint calculated.
Definition: treeBoundBox.C:180
bool contiguous< treeBoundBox >()
Data associated with treeBoundBox type are contiguous.
Definition: treeBoundBox.H:379
uint8_t direction
Definition: direction.H:45
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
direction faceBits(const point &) const
Code position of point on bounding box faces.
Definition: treeBoundBox.C:432
static faceList faces()
Return faces with correct point order.
Definition: boundBox.C:167
bool contains(const point &) const
Contains point? (inside or on edge)
Definition: boundBoxI.H:170
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
void searchOrder(const point &pt, FixedList< direction, 8 > &octantOrder) const
Calculates optimal order to look for nearest to point.
direction posBits(const point &) const
Position of point relative to bounding box.
Definition: treeBoundBox.C:465
friend Istream & operator>>(Istream &is, treeBoundBox &)
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.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
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:236
friend bool operator==(const treeBoundBox &, const treeBoundBox &)
treeBoundBox()
Construct null setting points to zero.
Definition: treeBoundBoxI.H:31
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Istream & operator>>(Istream &, directionInfo &)
static const edgeList edges
Edge to point addressing.
Definition: treeBoundBox.H:178
friend bool operator!=(const treeBoundBox &, const treeBoundBox &)
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:60
void calcExtremities(const point &pt, point &nearest, point &furthest) const
Calculate nearest and furthest (to point) vertex coords of.
Definition: treeBoundBox.C:500
bool contains(const vector &dir, const point &) const
Contains point (inside or on edge) and moving in direction.
Definition: treeBoundBox.C:394
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
scalar maxDist(const point &) const
Returns distance point to furthest away corner.
Definition: treeBoundBox.C:547
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
Direction is an 8-bit unsigned integer type used to represent the Cartesian directions etc...
label distanceCmp(const point &, const treeBoundBox &other) const
Compare distance to point with other bounding box.
Definition: treeBoundBox.C:557
Ostream & operator<<(Ostream &, const ensightPart &)
static const treeBoundBox invertedBox
As per boundBox::invertedBox, but with great instead of vGreat.
Definition: treeBoundBox.H:108
tmp< pointField > points() const
Vertex coordinates. In octant coding.
Definition: treeBoundBox.C:166
friend Ostream & operator<<(Ostream &os, const treeBoundBox &)
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:87
bool contiguous< boundBox >()
Data associated with boundBox type are contiguous.
Definition: boundBox.H:251
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
A class for managing temporary objects.
Definition: PtrList.H:53
bool operator!=(const particle &, const particle &)
Definition: particle.C:1204
static const scalar great
The great value used for greatBox and invertedBox.
Definition: treeBoundBox.H:102
Bits used for face coding.
Definition: treeBoundBox.H:137
treeBoundBox extend(const scalar s) const
Return asymetrically extended bounding box, with guaranteed.
Namespace for OpenFOAM.
void writeOBJ(const fileName &fName) const
Definition: treeBoundBox.C:619