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 | Copyright (C) 2011-2016 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 
93 private:
94 
95  //- To initialise edges.
96  static edgeList calcEdges(const label[12][2]);
97 
98  //- To initialise faceNormals.
99  static FixedList<vector, 6> calcFaceNormals();
100 
101 public:
102 
103  // Static data members
104 
105  //- The great value used for greatBox and invertedBox
106  static const scalar great;
107 
108  //- As per boundBox::greatBox, but with GREAT instead of VGREAT
109  static const treeBoundBox greatBox;
110 
111  //- As per boundBox::invertedBox, but with GREAT instead of VGREAT
112  static const treeBoundBox invertedBox;
113 
114  //- Bits used for octant/point coding.
115  // Every octant/corner point is the combination of three faces.
116  enum octantBit
117  {
118  RIGHTHALF = 0x1 << 0,
119  TOPHALF = 0x1 << 1,
120  FRONTHALF = 0x1 << 2
121  };
122 
123  //- Face codes
124  enum faceId
125  {
126  LEFT = 0,
127  RIGHT = 1,
128  BOTTOM = 2,
129  TOP = 3,
130  BACK = 4,
131  FRONT = 5
132  };
133 
134  //- Bits used for face coding
135  enum faceBit
136  {
137  NOFACE = 0,
138  LEFTBIT = 0x1 << LEFT, //1
139  RIGHTBIT = 0x1 << RIGHT, //2
140  BOTTOMBIT = 0x1 << BOTTOM, //4
141  TOPBIT = 0x1 << TOP, //8
142  BACKBIT = 0x1 << BACK, //16
143  FRONTBIT = 0x1 << FRONT, //32
144  };
145 
146  //- Edges codes.
147  // E01 = edge between 0 and 1.
148  enum edgeId
149  {
150  E01 = 0,
151  E13 = 1,
152  E23 = 2,
153  E02 = 3,
155  E45 = 4,
156  E57 = 5,
157  E67 = 6,
158  E46 = 7,
160  E04 = 8,
161  E15 = 9,
162  E37 = 10,
163  E26 = 11
164  };
165 
166  //- Face to point addressing
167  static const faceList faces;
168 
169  //- Edge to point addressing
170  static const edgeList edges;
171 
172  //- Per face the unit normal
173  static const FixedList<vector, 6> faceNormals;
174 
175 
176  // Constructors
177 
178  //- Construct null setting points to zero
179  inline treeBoundBox();
180 
181  //- Construct from a boundBox
182  explicit inline treeBoundBox(const boundBox& bb);
183 
184  //- Construct from components
185  inline treeBoundBox(const point& min, const point& max);
186 
187  //- Construct as the bounding box of the given pointField.
188  // Local processor domain only (no reduce as in boundBox)
189  explicit treeBoundBox(const UList<point>&);
190 
191  //- Construct as subset of points
192  // Local processor domain only (no reduce as in boundBox)
193  treeBoundBox(const UList<point>&, const labelUList& indices);
194 
195  //- Construct as subset of points
196  // The indices could be from edge/triFace etc.
197  // Local processor domain only (no reduce as in boundBox)
198  template<unsigned Size>
200  (
201  const UList<point>&,
202  const FixedList<label, Size>& indices
203  );
204 
205 
206  //- Construct from Istream
207  inline treeBoundBox(Istream&);
208 
209 
210  // Member functions
211 
212  // Access
213 
214  //- Typical dimension length,height,width
215  inline scalar typDim() const;
216 
217  //- Vertex coordinates. In octant coding.
218  tmp<pointField> points() const;
219 
220 
221  // Check
222 
223  //- Corner point given octant
224  inline point corner(const direction) const;
225 
226  //- Sub box given by octant number. Midpoint calculated.
227  treeBoundBox subBbox(const direction) const;
228 
229  //- Sub box given by octant number. Midpoint provided.
230  treeBoundBox subBbox(const point& mid, const direction) const;
231 
232  //- Returns octant number given point and the calculated midpoint.
233  inline direction subOctant
234  (
235  const point& pt
236  ) const;
237 
238  //- Returns octant number given point and midpoint.
239  static inline direction subOctant
240  (
241  const point& mid,
242  const point& pt
243  );
244 
245  //- Returns octant number given point and the calculated midpoint.
246  // onEdge set if the point is on edge of subOctant
247  inline direction subOctant
248  (
249  const point& pt,
250  bool& onEdge
251  ) const;
252 
253  //- Returns octant number given point and midpoint.
254  // onEdge set if the point is on edge of subOctant
255  static inline direction subOctant
256  (
257  const point& mid,
258  const point& pt,
259  bool& onEdge
260  );
261 
262  //- Returns octant number given intersection and midpoint.
263  // onEdge set if the point is on edge of subOctant
264  // If onEdge, the direction vector determines which octant to use
265  // (acc. to which octant the point would be if it were moved
266  // along dir)
267  static inline direction subOctant
268  (
269  const point& mid,
270  const vector& dir,
271  const point& pt,
272  bool& onEdge
273  );
274 
275  //- Calculates optimal order to look for nearest to point.
276  // First will be the octant containing the point,
277  // second the octant with boundary nearest to the point etc.
278  inline void searchOrder
279  (
280  const point& pt,
281  FixedList<direction, 8>& octantOrder
282  ) const;
283 
284  //- Overlaps other bounding box?
285  using boundBox::overlaps;
286 
287  //- Intersects segment; set point to intersection position and face,
288  // return true if intersection found.
289  // (pt argument used during calculation even if not intersecting).
290  // Calculates intersections from outside supplied vector
291  // (overallStart, overallVec). This is so when
292  // e.g. tracking through lots of consecutive boxes
293  // (typical octree) we're not accumulating truncation errors. Set
294  // to start, (end-start) if not used.
295  bool intersects
296  (
297  const point& overallStart,
298  const vector& overallVec,
299  const point& start,
300  const point& end,
301  point& pt,
302  direction& ptBits
303  ) const;
304 
305  //- Like above but does not return faces point is on
306  bool intersects
307  (
308  const point& start,
309  const point& end,
310  point& pt
311  ) const;
312 
313  //- Contains point or other bounding box?
314  using boundBox::contains;
315 
316  //- Contains point (inside or on edge) and moving in direction
317  // dir would cause it to go inside.
318  bool contains(const vector& dir, const point&) const;
319 
320  //- Code position of point on bounding box faces
321  direction faceBits(const point&) const;
322 
323  //- Position of point relative to bounding box
324  direction posBits(const point&) const;
325 
326  //- Calculate nearest and furthest (to point) vertex coords of
327  // bounding box
328  void calcExtremities
329  (
330  const point& pt,
331  point& nearest,
332  point& furthest
333  ) const;
334 
335  //- Returns distance point to furthest away corner.
336  scalar maxDist(const point&) const;
337 
338  //- Compare distance to point with other bounding box
339  // return:
340  // -1 : all vertices of my bounding box are nearer than any of
341  // other
342  // +1 : all vertices of my bounding box are further away than
343  // any of other
344  // 0 : none of the above.
345  label distanceCmp(const point&, const treeBoundBox& other) const;
346 
347  //- Return slightly wider bounding box
348  // Extends all dimensions with s*span*Random::scalar01()
349  // and guarantees in any direction s*mag(span) minimum width
350  inline treeBoundBox extend(Random&, const scalar s) const;
351 
352  // Friend Operators
353 
354  friend bool operator==(const treeBoundBox&, const treeBoundBox&);
355  friend bool operator!=(const treeBoundBox&, const treeBoundBox&);
356 
357  // IOstream operator
358 
359  friend Istream& operator>>(Istream& is, treeBoundBox&);
360  friend Ostream& operator<<(Ostream& os, const treeBoundBox&);
361 };
362 
363 
364 //- Data associated with treeBoundBox type are contiguous
365 template<>
366 inline bool contiguous<treeBoundBox>() {return contiguous<boundBox>();}
367 
368 
369 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
370 
371 } // End namespace Foam
372 
373 #include "treeBoundBoxI.H"
374 
375 #ifdef NoRepository
376  #include "treeBoundBoxTemplates.C"
377 #endif
378 
379 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
380 
381 #endif
382 
383 // ************************************************************************* //
static const FixedList< vector, 6 > faceNormals
Per face the unit normal.
Definition: treeBoundBox.H:172
edgeId
Edges codes.
Definition: treeBoundBox.H:147
scalar maxDist(const point &) const
Returns distance point to furthest away corner.
Definition: treeBoundBox.C:546
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
uint8_t direction
Definition: direction.H:46
void calcExtremities(const point &pt, point &nearest, point &furthest) const
Calculate nearest and furthest (to point) vertex coords of.
Definition: treeBoundBox.C:499
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:53
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
static const treeBoundBox greatBox
As per boundBox::greatBox, but with GREAT instead of VGREAT.
Definition: treeBoundBox.H:108
bool contiguous< treeBoundBox >()
Data associated with treeBoundBox type are contiguous.
Definition: treeBoundBox.H:365
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
static faceList faces()
Return faces with correct point order.
Definition: boundBox.C:167
treeBoundBox extend(Random &, const scalar s) const
Return slightly wider bounding box.
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
direction faceBits(const point &) const
Code position of point on bounding box faces.
Definition: treeBoundBox.C:431
friend Istream & operator>>(Istream &is, treeBoundBox &)
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:235
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))
direction posBits(const point &) const
Position of point relative to bounding box.
Definition: treeBoundBox.C:464
friend bool operator==(const treeBoundBox &, const treeBoundBox &)
scalar typDim() const
Typical dimension length,height,width.
Definition: treeBoundBoxI.H:57
void searchOrder(const point &pt, FixedList< direction, 8 > &octantOrder) const
Calculates optimal order to look for nearest to point.
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 &)
direction subOctant(const point &pt) const
Returns octant number given point and the calculated midpoint.
Definition: treeBoundBoxI.H:75
static const edgeList edges
Edge to point addressing.
Definition: treeBoundBox.H:169
faceBit
Bits used for face coding.
Definition: treeBoundBox.H:134
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
Simple random number generator.
Definition: Random.H:49
faceId
Face codes.
Definition: treeBoundBox.H:123
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
label distanceCmp(const point &, const treeBoundBox &other) const
Compare distance to point with other bounding box.
Definition: treeBoundBox.C:556
octantBit
Bits used for octant/point coding.
Definition: treeBoundBox.H:115
bool contains(const vector &dir, const point &) const
Contains point (inside or on edge) and moving in direction.
Definition: treeBoundBox.C:393
point nearest(const point &) const
Return the nearest point on the boundBox to the supplied point.
Definition: boundBox.C:303
point corner(const direction) const
Corner point given octant.
Definition: treeBoundBoxI.H:63
Direction is an 8-bit unsigned integer type used to represent the Cartesian directions etc...
tmp< pointField > points() const
Vertex coordinates. In octant coding.
Definition: treeBoundBox.C:165
Ostream & operator<<(Ostream &, const ensightPart &)
static const treeBoundBox invertedBox
As per boundBox::invertedBox, but with GREAT instead of VGREAT.
Definition: treeBoundBox.H:111
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
A class for managing temporary objects.
Definition: PtrList.H:54
treeBoundBox subBbox(const direction) const
Sub box given by octant number. Midpoint calculated.
Definition: treeBoundBox.C:179
bool operator!=(const particle &, const particle &)
Definition: particle.C:145
bool contains(const point &) const
Contains point? (inside or on edge)
Definition: boundBoxI.H:170
bool overlaps(const boundBox &) const
Overlaps/touches boundingBox?
Definition: boundBoxI.H:120
static const scalar great
The great value used for greatBox and invertedBox.
Definition: treeBoundBox.H:105
Namespace for OpenFOAM.