treeBoundBoxI.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 OpenFOAM Foundation
9 Copyright (C) 2017-2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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\*---------------------------------------------------------------------------*/
28
29#include "treeBoundBox.H"
30#include "Random.H"
31
32// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33
35:
36 boundBox()
37{}
38
39
41:
42 boundBox(bb)
43{}
44
45
47:
48 boundBox(pt)
49{}
50
51
53:
55{}
56
57
59:
60 boundBox(is)
61{}
62
63
64// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65
66inline Foam::scalar Foam::treeBoundBox::typDim() const
67{
68 return avgDim();
69}
70
71
73{
74 return point
75 (
76 (octant & RIGHTHALF) ? max().x() : min().x(),
77 (octant & TOPHALF) ? max().y() : min().y(),
78 (octant & FRONTHALF) ? max().z() : min().z()
79 );
80}
81
82
83// Returns octant in which point resides. Reverse of subBbox.
85{
86 return subOctant(centre(), pt);
87}
88
89
90// Returns octant in which point resides. Reverse of subBbox.
91// Precalculated midpoint
93(
94 const point& mid,
95 const point& pt
96)
97{
98 direction octant = 0;
99
100 if (pt.x() > mid.x())
101 {
102 octant |= treeBoundBox::RIGHTHALF;
103 }
104
105 if (pt.y() > mid.y())
106 {
107 octant |= treeBoundBox::TOPHALF;
108 }
109
110 if (pt.z() > mid.z())
111 {
112 octant |= treeBoundBox::FRONTHALF;
113 }
114
115 return octant;
116}
117
118
119// Returns octant in which point resides. Reverse of subBbox.
120// Flags point exactly on edge.
122(
123 const point& pt,
124 bool& onEdge
125) const
126{
127 return subOctant(centre(), pt, onEdge);
128}
129
130
131// Returns octant in which point resides. Reverse of subBbox.
132// Precalculated midpoint
134(
135 const point& mid,
136 const point& pt,
137 bool& onEdge
138)
139{
140 direction octant = 0;
141 onEdge = false;
142
143 if (pt.x() > mid.x())
144 {
145 octant |= treeBoundBox::RIGHTHALF;
146 }
147 else if (pt.x() == mid.x())
148 {
149 onEdge = true;
150 }
151
152 if (pt.y() > mid.y())
153 {
154 octant |= treeBoundBox::TOPHALF;
155 }
156 else if (pt.y() == mid.y())
157 {
158 onEdge = true;
159 }
160
161 if (pt.z() > mid.z())
162 {
163 octant |= treeBoundBox::FRONTHALF;
164 }
165 else if (pt.z() == mid.z())
166 {
167 onEdge = true;
168 }
169
170 return octant;
171}
172
173
174// Returns octant in which intersection resides.
175// Precalculated midpoint. If the point is on the dividing line between
176// the octants the direction vector determines which octant to use
177// (i.e. in which octant the point would be if it were moved along dir)
179(
180 const point& mid,
181 const vector& dir,
182 const point& pt,
183 bool& onEdge
184)
185{
186 direction octant = 0;
187 onEdge = false;
188
189 if (pt.x() > mid.x())
190 {
191 octant |= treeBoundBox::RIGHTHALF;
192 }
193 else if (pt.x() == mid.x())
194 {
195 onEdge = true;
196 if (dir.x() > 0)
197 {
198 octant |= treeBoundBox::RIGHTHALF;
199 }
200 }
201
202 if (pt.y() > mid.y())
203 {
204 octant |= treeBoundBox::TOPHALF;
205 }
206 else if (pt.y() == mid.y())
207 {
208 onEdge = true;
209 if (dir.y() > 0)
210 {
211 octant |= treeBoundBox::TOPHALF;
212 }
213 }
214
215 if (pt.z() > mid.z())
216 {
217 octant |= treeBoundBox::FRONTHALF;
218 }
219 else if (pt.z() == mid.z())
220 {
221 onEdge = true;
222 if (dir.z() > 0)
223 {
224 octant |= treeBoundBox::FRONTHALF;
225 }
226 }
227
228 return octant;
229}
230
231
232// Returns reference to octantOrder which defines the
233// order to do the search.
235(
236 const point& pt,
237 FixedList<direction,8>& octantOrder
238) const
239{
240 vector dist = centre() - pt;
241
242 direction octant = 0;
243
244 if (dist.x() < 0)
245 {
246 octant |= treeBoundBox::RIGHTHALF;
247 dist.x() *= -1;
248 }
249
250 if (dist.y() < 0)
251 {
252 octant |= treeBoundBox::TOPHALF;
253 dist.y() *= -1;
254 }
255
256 if (dist.z() < 0)
257 {
258 octant |= treeBoundBox::FRONTHALF;
259 dist.z() *= -1;
260 }
261
262 direction min = 0;
263 direction mid = 0;
264 direction max = 0;
265
266 if (dist.x() < dist.y())
267 {
268 if (dist.y() < dist.z())
269 {
273 }
274 else if (dist.z() < dist.x())
275 {
279 }
280 else
281 {
285 }
286 }
287 else
288 {
289 if (dist.z() < dist.y())
290 {
294 }
295 else if (dist.x() < dist.z())
296 {
300 }
301 else
302 {
306 }
307 }
308
309 // Primary subOctant
310 octantOrder[0] = octant;
311 // subOctants joined to the primary by faces.
312 octantOrder[1] = octant ^ min;
313 octantOrder[2] = octant ^ mid;
314 octantOrder[3] = octant ^ max;
315 // subOctants joined to the primary by edges.
316 octantOrder[4] = octantOrder[1] ^ mid;
317 octantOrder[5] = octantOrder[1] ^ max;
318 octantOrder[6] = octantOrder[2] ^ max;
319 // subOctants joined to the primary by corner.
320 octantOrder[7] = octantOrder[4] ^ max;
321}
322
323
325(
326 Random& rndGen,
327 const scalar s
328) const
329{
330 treeBoundBox bb(*this);
331
332 vector newSpan = bb.span();
333
334 // Make 3D
335 scalar minSpan = s * Foam::mag(newSpan);
336
337 for (direction dir = 0; dir < vector::nComponents; ++dir)
338 {
339 newSpan[dir] = Foam::max(newSpan[dir], minSpan);
340 }
341
342 bb.min() -= cmptMultiply(s*rndGen.sample01<vector>(), newSpan);
343 bb.max() += cmptMultiply(s*rndGen.sample01<vector>(), newSpan);
344
345 return bb;
346}
347
348
349// * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * * //
350
351inline bool Foam::operator==(const treeBoundBox& a, const treeBoundBox& b)
352{
353 return static_cast<const boundBox&>(a) == static_cast<const boundBox&>(b);
354}
355
356
357inline bool Foam::operator!=(const treeBoundBox& a, const treeBoundBox& b)
358{
359 return !(a == b);
360}
361
362// ************************************************************************* //
scalar y
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: FixedList.H:81
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
Random number generator.
Definition: Random.H:60
Type sample01()
Return a sample whose components lie in the range [0,1].
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:64
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:127
static constexpr direction nComponents
Number of components in bool is 1.
Definition: bool.H:98
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:89
treeBoundBox extend(Random &rndGen, const scalar s) const
Return slightly wider bounding box.
void searchOrder(const point &pt, FixedList< direction, 8 > &octantOrder) const
Calculates optimal order to look for nearest to point.
point corner(const direction octant) const
Corner point of given octant.
Definition: treeBoundBoxI.H:72
direction subOctant(const point &pt) const
Returns octant number given point and the calculated midpoint.
Definition: treeBoundBoxI.H:84
@ TOPHALF
2: positive y-direction
Definition: treeBoundBox.H:99
@ RIGHTHALF
1: positive x-direction
Definition: treeBoundBox.H:98
@ FRONTHALF
4: positive z-direction
Definition: treeBoundBox.H:100
scalar typDim() const
Typical dimension length,height,width.
Definition: treeBoundBoxI.H:66
treeBoundBox()
Construct without any points - an inverted bounding box.
Definition: treeBoundBoxI.H:34
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))
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:239
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
vector point
Point is a vector.
Definition: point.H:43
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
uint8_t direction
Definition: direction.H:56
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
volScalarField & b
Definition: createFields.H:27
Random rndGen
Definition: createFields.H:23