treeBoundBoxI.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 \*---------------------------------------------------------------------------*/
25 
26 #include "treeBoundBox.H"
27 #include "randomGenerator.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
32 :
33  boundBox()
34 {}
35 
36 
38 :
39  boundBox(min, max)
40 {}
41 
42 
44 :
45  boundBox(bb)
46 {}
47 
48 
50 :
51  boundBox(dict)
52 {}
53 
54 
56 :
57  boundBox(is)
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
62 
63 inline Foam::scalar Foam::treeBoundBox::typDim() const
64 {
65  return avgDim();
66 }
67 
68 
70 {
71  return point
72  (
73  (octant & octantBit::rightHalf) ? max().x() : min().x(),
74  (octant & octantBit::topHalf) ? max().y() : min().y(),
75  (octant & octantBit::frontHalf) ? max().z() : min().z()
76  );
77 }
78 
79 
80 // Returns octant in which point resides. Reverse of subBbox.
82 {
83  return subOctant(midpoint(), pt);
84 }
85 
86 
87 // Returns octant in which point resides. Reverse of subBbox.
88 // Precalculated midpoint
90 (
91  const point& mid,
92  const point& pt
93 )
94 {
95  direction octant = 0;
96 
97  if (pt.x() > mid.x())
98  {
99  octant |= octantBit::rightHalf;
100  }
101 
102  if (pt.y() > mid.y())
103  {
104  octant |= octantBit::topHalf;
105  }
106 
107  if (pt.z() > mid.z())
108  {
109  octant |= octantBit::frontHalf;
110  }
111 
112  return octant;
113 }
114 
115 
116 // Returns octant in which point resides. Reverse of subBbox.
117 // Flags point exactly on edge.
119 (
120  const point& pt,
121  bool& onEdge
122 ) const
123 {
124  return subOctant(midpoint(), pt, onEdge);
125 }
126 
127 
128 // Returns octant in which point resides. Reverse of subBbox.
129 // Precalculated midpoint
131 (
132  const point& mid,
133  const point& pt,
134  bool& onEdge
135 )
136 {
137  direction octant = 0;
138  onEdge = false;
139 
140  if (pt.x() > mid.x())
141  {
142  octant |= octantBit::rightHalf;
143  }
144  else if (pt.x() == mid.x())
145  {
146  onEdge = true;
147  }
148 
149  if (pt.y() > mid.y())
150  {
151  octant |= octantBit::topHalf;
152  }
153  else if (pt.y() == mid.y())
154  {
155  onEdge = true;
156  }
157 
158  if (pt.z() > mid.z())
159  {
160  octant |= octantBit::frontHalf;
161  }
162  else if (pt.z() == mid.z())
163  {
164  onEdge = true;
165  }
166 
167  return octant;
168 }
169 
170 
171 // Returns octant in which intersection resides.
172 // Precalculated midpoint. If the point is on the dividing line between
173 // the octants the direction vector determines which octant to use
174 // (i.e. in which octant the point would be if it were moved along dir)
176 (
177  const point& mid,
178  const vector& dir,
179  const point& pt,
180  bool& onEdge
181 )
182 {
183  direction octant = 0;
184  onEdge = false;
185 
186  if (pt.x() > mid.x())
187  {
188  octant |= octantBit::rightHalf;
189  }
190  else if (pt.x() == mid.x())
191  {
192  onEdge = true;
193  if (dir.x() > 0)
194  {
195  octant |= octantBit::rightHalf;
196  }
197  }
198 
199  if (pt.y() > mid.y())
200  {
201  octant |= octantBit::topHalf;
202  }
203  else if (pt.y() == mid.y())
204  {
205  onEdge = true;
206  if (dir.y() > 0)
207  {
208  octant |= octantBit::topHalf;
209  }
210  }
211 
212  if (pt.z() > mid.z())
213  {
214  octant |= octantBit::frontHalf;
215  }
216  else if (pt.z() == mid.z())
217  {
218  onEdge = true;
219  if (dir.z() > 0)
220  {
221  octant |= octantBit::frontHalf;
222  }
223  }
224 
225  return octant;
226 }
227 
228 
229 // Returns reference to octantOrder which defines the
230 // order to do the search.
232 (
233  const point& pt,
234  FixedList<direction,8>& octantOrder
235 ) const
236 {
237  vector dist = midpoint() - pt;
238 
239  direction octant = 0;
240 
241  if (dist.x() < 0)
242  {
243  octant |= octantBit::rightHalf;
244  dist.x() *= -1;
245  }
246 
247  if (dist.y() < 0)
248  {
249  octant |= octantBit::topHalf;
250  dist.y() *= -1;
251  }
252 
253  if (dist.z() < 0)
254  {
255  octant |= octantBit::frontHalf;
256  dist.z() *= -1;
257  }
258 
259  direction min = 0;
260  direction mid = 0;
261  direction max = 0;
262 
263  if (dist.x() < dist.y())
264  {
265  if (dist.y() < dist.z())
266  {
267  min = octantBit::rightHalf;
268  mid = octantBit::topHalf;
269  max = octantBit::frontHalf;
270  }
271  else if (dist.z() < dist.x())
272  {
273  min = octantBit::frontHalf;
274  mid = octantBit::rightHalf;
275  max = octantBit::topHalf;
276  }
277  else
278  {
279  min = octantBit::rightHalf;
280  mid = octantBit::frontHalf;
281  max = octantBit::topHalf;
282  }
283  }
284  else
285  {
286  if (dist.z() < dist.y())
287  {
288  min = octantBit::frontHalf;
289  mid = octantBit::topHalf;
290  max = octantBit::rightHalf;
291  }
292  else if (dist.x() < dist.z())
293  {
294  min = octantBit::topHalf;
295  mid = octantBit::rightHalf;
296  max = octantBit::frontHalf;
297  }
298  else
299  {
300  min = octantBit::topHalf;
301  mid = octantBit::frontHalf;
302  max = octantBit::rightHalf;
303  }
304  }
305 
306  // Primary subOctant
307  octantOrder[0] = octant;
308  // subOctants joined to the primary by faces.
309  octantOrder[1] = octant ^ min;
310  octantOrder[2] = octant ^ mid;
311  octantOrder[3] = octant ^ max;
312  // subOctants joined to the primary by edges.
313  octantOrder[4] = octantOrder[1] ^ mid;
314  octantOrder[5] = octantOrder[1] ^ max;
315  octantOrder[6] = octantOrder[2] ^ max;
316  // subOctants joined to the primary by corner.
317  octantOrder[7] = octantOrder[4] ^ max;
318 }
319 
320 
322 {
323  // Numbers that don't approximate rational fractions with which to make the
324  // box asymmetric. These are between one and two.
325  static const vector a((sqrt(5.0) + 1)/2, sqrt(2.0), (sqrt(13.0) - 1)/2);
326  static const vector b(a.y(), a.z(), a.x());
327 
328  treeBoundBox bb(*this);
329 
330  const scalar delta = s*Foam::mag(bb.span());
331  bb.min() -= Foam::max(delta*a, vector::uniform(rootVSmall));
332  bb.max() += Foam::max(delta*b, vector::uniform(rootVSmall));
333 
334  return bb;
335 }
336 
337 
338 // ************************************************************************* //
scalar y
scalar delta
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
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:168
const Cmpt & z() const
Definition: VectorI.H:87
const Cmpt & y() const
Definition: VectorI.H:81
const Cmpt & x() const
Definition: VectorI.H:75
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:59
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
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:90
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:90
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.
point corner(const direction) const
Corner point given octant.
Definition: treeBoundBoxI.H:69
direction subOctant(const point &pt) const
Returns octant number given point and the calculated midpoint.
Definition: treeBoundBoxI.H:81
scalar typDim() const
Typical dimension length,height,width.
Definition: treeBoundBoxI.H:63
treeBoundBox()
Construct null setting points to zero.
Definition: treeBoundBoxI.H:31
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))
volScalarField & b
Definition: createFields.H:25
vector point
Point is a vector.
Definition: point.H:41
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
uint8_t direction
Definition: direction.H:45
dictionary dict