Loading...
Searching...
No Matches
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-2023 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 * * * * * * * * * * * * * * //
37{}
38
43{}
44
49{}
50
55{}
56
57
60 boundBox(is)
61{}
62
63
64// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65
66inline Foam::point Foam::treeBoundBox::corner(const direction octant) const
67{
68 return point
69 (
70 (octant & RIGHTHALF) ? max().x() : min().x(),
71 (octant & TOPHALF) ? max().y() : min().y(),
72 (octant & FRONTHALF) ? max().z() : min().z()
73 );
74}
75
76
77// Returns octant in which point resides. Reverse of subBbox.
80 return subOctant(centre(), pt);
81}
82
83
84// Returns octant in which point resides. Reverse of subBbox.
85// Precalculated midpoint
87(
88 const point& mid,
89 const point& pt
90)
91{
92 direction octant = 0;
93
94 if (pt.x() > mid.x())
95 {
97 }
98
99 if (pt.y() > mid.y())
100 {
101 octant |= treeBoundBox::TOPHALF;
102 }
103
104 if (pt.z() > mid.z())
105 {
106 octant |= treeBoundBox::FRONTHALF;
107 }
109 return octant;
110}
111
112
113// Returns octant in which point resides. Reverse of subBbox.
114// Flags point exactly on edge.
116(
117 const point& pt,
118 bool& onEdge
119) const
121 return subOctant(centre(), pt, onEdge);
122}
123
124
125// Returns octant in which point resides. Reverse of subBbox.
126// Precalculated midpoint
128(
129 const point& mid,
130 const point& pt,
131 bool& onEdge
132)
133{
134 direction octant = 0;
135 onEdge = false;
136
137 if (pt.x() > mid.x())
138 {
139 octant |= treeBoundBox::RIGHTHALF;
140 }
141 else if (pt.x() == mid.x())
142 {
143 onEdge = true;
144 }
145
146 if (pt.y() > mid.y())
147 {
148 octant |= treeBoundBox::TOPHALF;
149 }
150 else if (pt.y() == mid.y())
151 {
152 onEdge = true;
153 }
154
155 if (pt.z() > mid.z())
156 {
157 octant |= treeBoundBox::FRONTHALF;
158 }
159 else if (pt.z() == mid.z())
160 {
161 onEdge = true;
162 }
163
164 return octant;
166
167
168// Returns octant in which intersection resides.
169// Precalculated midpoint. If the point is on the dividing line between
170// the octants the direction vector determines which octant to use
171// (i.e. in which octant the point would be if it were moved along dir)
173(
174 const point& mid,
175 const vector& dir,
176 const point& pt,
177 bool& onEdge
178)
179{
180 direction octant = 0;
181 onEdge = false;
182
183 if (pt.x() > mid.x())
184 {
185 octant |= treeBoundBox::RIGHTHALF;
186 }
187 else if (pt.x() == mid.x())
188 {
189 onEdge = true;
190 if (dir.x() > 0)
191 {
192 octant |= treeBoundBox::RIGHTHALF;
193 }
194 }
195
196 if (pt.y() > mid.y())
197 {
198 octant |= treeBoundBox::TOPHALF;
199 }
200 else if (pt.y() == mid.y())
201 {
202 onEdge = true;
203 if (dir.y() > 0)
204 {
205 octant |= treeBoundBox::TOPHALF;
206 }
207 }
208
209 if (pt.z() > mid.z())
210 {
211 octant |= treeBoundBox::FRONTHALF;
212 }
213 else if (pt.z() == mid.z())
214 {
215 onEdge = true;
216 if (dir.z() > 0)
217 {
218 octant |= treeBoundBox::FRONTHALF;
220 }
221
222 return octant;
223}
224
225
227(
228 const point& pt,
229 FixedList<direction, 8>& octantOrder
230) const
231{
232 vector dist = centre() - pt;
233
234 direction octant = 0;
235
236 if (dist.x() < 0)
237 {
238 octant |= treeBoundBox::RIGHTHALF;
239 dist.x() *= -1;
240 }
241
242 if (dist.y() < 0)
243 {
244 octant |= treeBoundBox::TOPHALF;
245 dist.y() *= -1;
246 }
247
248 if (dist.z() < 0)
249 {
250 octant |= treeBoundBox::FRONTHALF;
251 dist.z() *= -1;
252 }
253
254 direction min = 0;
255 direction mid = 0;
256 direction max = 0;
257
258 if (dist.x() < dist.y())
259 {
260 if (dist.y() < dist.z())
261 {
265 }
266 else if (dist.z() < dist.x())
267 {
271 }
272 else
273 {
277 }
278 }
279 else
280 {
281 if (dist.z() < dist.y())
282 {
286 }
287 else if (dist.x() < dist.z())
288 {
292 }
293 else
294 {
298 }
299 }
300
301 // Primary subOctant
302 octantOrder[0] = octant;
303 // subOctants joined to the primary by faces.
304 octantOrder[1] = octant ^ min;
305 octantOrder[2] = octant ^ mid;
306 octantOrder[3] = octant ^ max;
307 // subOctants joined to the primary by edges.
308 octantOrder[4] = octantOrder[1] ^ mid;
309 octantOrder[5] = octantOrder[1] ^ max;
310 octantOrder[6] = octantOrder[2] ^ max;
311 // subOctants joined to the primary by corner.
312 octantOrder[7] = octantOrder[4] ^ max;
313}
314
315
319 FixedList<direction, 8> octantOrder;
320 searchOrder(pt, octantOrder);
321 return octantOrder;
322}
323
324
326(
327 const direction octant,
328 const point& centre,
329 const scalar radiusSqr
330) const
331{
332 // Accelerated version of
333 // subBbox(octant).overlaps(centre, radiusSqr)
334
335 // NB: ordering of corners is irrelevant
336 return box_sphere_overlaps
337 (
338 this->centre(),
339 corner(octant),
340 centre,
341 radiusSqr
342 );
343}
344
345
347(
348 const linePointRef& ln,
350) const
351{
352 return intersects(ln.start(), ln.end(), pt);
353}
354
355
357(
358 Random& rndGen,
359 const scalar factor
360) const
361{
362 treeBoundBox bb(*this);
364 bb.inflate(rndGen, factor);
365
366 return bb;
367}
368
369
371(
372 Random& rndGen,
373 const scalar factor,
374 const scalar delta
375) const
376{
377 treeBoundBox bb(*this);
378
379 bb.inflate(rndGen, factor, delta);
381 return bb;
382}
383
384
385// * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * * //
387inline bool Foam::operator==(const treeBoundBox& a, const treeBoundBox& b)
388{
389 return static_cast<const boundBox&>(a) == static_cast<const boundBox&>(b);
390}
391
392
393inline bool Foam::operator!=(const treeBoundBox& a, const treeBoundBox& b)
394{
395 return !(a == b);
396}
397
398// ************************************************************************* //
scalar y
scalar delta
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
An ordered pair of two objects of type <T> with first() and second() elements.
Definition Pair.H:66
Random number generator.
Definition Random.H:56
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector.H:135
const Cmpt & z() const noexcept
Access to the vector z component.
Definition Vector.H:145
const Cmpt & y() const noexcept
Access to the vector y component.
Definition Vector.H:140
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
const point & max() const noexcept
Maximum describing the bounding box.
Definition boundBoxI.H:168
static bool box_sphere_overlaps(const point &corner0, const point &corner1, const point &centre, const scalar radiusSqr)
Test for overlap of box and boundingSphere (centre + sqr(radius)).
Definition boundBoxI.H:403
const point & min() const noexcept
Minimum describing the bounding box.
Definition boundBoxI.H:162
boundBox()
Default construct: an inverted bounding box.
Definition boundBoxI.H:101
void inflate(const scalar factor)
Expand box by factor*mag(span) in all dimensions.
Definition boundBoxI.H:381
point centre() const
The centre (midpoint) of the bounding box.
Definition boundBoxI.H:186
Standard boundBox with extra functionality for use in octree.
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.
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,.
point corner(const direction octant) const
Corner point of given octant.
treeBoundBox()=default
Default construct: an inverted bounding box.
bool subOverlaps(const direction octant, const boundBox &bb) const
Does sub-octant overlap/touch boundingBox?
direction subOctant(const point &pt) const
Returns octant number given point and the calculated midpoint.
@ TOPHALF
2: positive y-direction
@ RIGHTHALF
1: positive x-direction
@ FRONTHALF
4: positive z-direction
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
Definition pTraits.H:51
volScalarField & p
Namespace for OpenFOAM.
bool operator!=(const eddy &a, const eddy &b)
Definition eddy.H:297
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
line< point, const point & > linePointRef
A line using referred points.
Definition line.H:66
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
uint8_t direction
Definition direction.H:49
vector point
Point is a vector.
Definition point.H:37
Vector< scalar > vector
Definition vector.H:57
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition POSIX.C:1239
volScalarField & b
Random rndGen