treeBoundBox.C
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-2016 OpenFOAM Foundation
9 Copyright (C) 2017-2020 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
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32
34({
35 face({0, 4, 6, 2}), // 0: x-min, left
36 face({1, 3, 7, 5}), // 1: x-max, right
37 face({0, 1, 5, 4}), // 2: y-min, bottom
38 face({2, 6, 7, 3}), // 3: y-max, top
39 face({0, 2, 3, 1}), // 4: z-min, back
40 face({4, 5, 7, 6}) // 5: z-max, front
41});
42
44({
45 {0, 1}, // 0
46 {1, 3},
47 {2, 3}, // 2
48 {0, 2},
49 {4, 5}, // 4
50 {5, 7},
51 {6, 7}, // 6
52 {4, 6},
53 {0, 4}, // 8
54 {1, 5},
55 {3, 7}, // 10
56 {2, 6}
57});
58
59
60// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
61
63:
64 boundBox(points, false)
65{
66 if (points.empty())
67 {
69 << "No bounding box for zero-sized pointField" << nl;
70 }
71}
72
73
75(
76 const UList<point>& points,
77 const labelUList& indices
78)
79:
80 boundBox(points, indices, false)
81{
82 if (points.empty() || indices.empty())
83 {
85 << "No bounding box for zero-sized pointField" << nl;
86 }
87}
88
89
90// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91
93{
94 auto tpts = tmp<pointField>::New(8);
95 auto& pts = tpts.ref();
96
97 forAll(pts, octant)
98 {
99 pts[octant] = corner(octant);
100 }
101
102 return tpts;
103}
104
105
107{
108 return subBbox(centre(), octant);
109}
110
111
113(
114 const point& mid,
115 const direction octant
116) const
117{
118 if (octant > 7)
119 {
121 << "octant should be [0..7]"
122 << abort(FatalError);
123 }
124
125 // start with a copy of this bounding box and adjust limits accordingly
126 treeBoundBox subBb(*this);
127 point& bbMin = subBb.min();
128 point& bbMax = subBb.max();
129
130 if (octant & treeBoundBox::RIGHTHALF)
131 {
132 bbMin.x() = mid.x(); // mid -> max
133 }
134 else
135 {
136 bbMax.x() = mid.x(); // min -> mid
137 }
138
139 if (octant & treeBoundBox::TOPHALF)
140 {
141 bbMin.y() = mid.y(); // mid -> max
142 }
143 else
144 {
145 bbMax.y() = mid.y(); // min -> mid
146 }
147
148 if (octant & treeBoundBox::FRONTHALF)
149 {
150 bbMin.z() = mid.z(); // mid -> max
151 }
152 else
153 {
154 bbMax.z() = mid.z(); // min -> mid
155 }
156
157 return subBb;
158}
159
160
162(
163 const point& overallStart,
164 const vector& overallVec,
165 const point& start,
166 const point& end,
167 point& pt,
168 direction& ptOnFaces
169) const
170{
171 // Sutherlands algorithm:
172 // loop
173 // - start = intersection of line with one of the planes bounding
174 // the bounding box
175 // - stop if start inside bb (return true)
176 // - stop if start and end in same 'half' (e.g. both above bb)
177 // (return false)
178 //
179 // Uses posBits to efficiently determine 'half' in which start and end
180 // point are.
181 //
182 // Note:
183 // - sets coordinate to exact position: e.g. pt.x() = min().x()
184 // since plane intersect routine might have truncation error.
185 // This makes sure that posBits tests 'inside'
186
187 const direction endBits = posBits(end);
188 pt = start;
189
190 // Allow maximum of 3 clips.
191 for (label i = 0; i < 4; ++i)
192 {
193 direction ptBits = posBits(pt);
194
195 if (ptBits == 0)
196 {
197 // pt inside bb
198 ptOnFaces = faceBits(pt);
199 return true;
200 }
201
202 if ((ptBits & endBits) != 0)
203 {
204 // pt and end in same block outside of bb
205 ptOnFaces = faceBits(pt);
206 return false;
207 }
208
209 if (ptBits & LEFTBIT)
210 {
211 // Intersect with plane V=min, n=-1,0,0
212 if (Foam::mag(overallVec.x()) > VSMALL)
213 {
214 scalar s = (min().x() - overallStart.x())/overallVec.x();
215 pt.x() = min().x();
216 pt.y() = overallStart.y() + overallVec.y()*s;
217 pt.z() = overallStart.z() + overallVec.z()*s;
218 }
219 else
220 {
221 // Vector not in x-direction. But still intersecting bb planes.
222 // So must be close - just snap to bb.
223 pt.x() = min().x();
224 }
225 }
226 else if (ptBits & RIGHTBIT)
227 {
228 // Intersect with plane V=max, n=1,0,0
229 if (Foam::mag(overallVec.x()) > VSMALL)
230 {
231 scalar s = (max().x() - overallStart.x())/overallVec.x();
232 pt.x() = max().x();
233 pt.y() = overallStart.y() + overallVec.y()*s;
234 pt.z() = overallStart.z() + overallVec.z()*s;
235 }
236 else
237 {
238 pt.x() = max().x();
239 }
240 }
241 else if (ptBits & BOTTOMBIT)
242 {
243 // Intersect with plane V=min, n=0,-1,0
244 if (Foam::mag(overallVec.y()) > VSMALL)
245 {
246 scalar s = (min().y() - overallStart.y())/overallVec.y();
247 pt.x() = overallStart.x() + overallVec.x()*s;
248 pt.y() = min().y();
249 pt.z() = overallStart.z() + overallVec.z()*s;
250 }
251 else
252 {
253 pt.x() = min().y();
254 }
255 }
256 else if (ptBits & TOPBIT)
257 {
258 // Intersect with plane V=max, n=0,1,0
259 if (Foam::mag(overallVec.y()) > VSMALL)
260 {
261 scalar s = (max().y() - overallStart.y())/overallVec.y();
262 pt.x() = overallStart.x() + overallVec.x()*s;
263 pt.y() = max().y();
264 pt.z() = overallStart.z() + overallVec.z()*s;
265 }
266 else
267 {
268 pt.y() = max().y();
269 }
270 }
271 else if (ptBits & BACKBIT)
272 {
273 // Intersect with plane V=min, n=0,0,-1
274 if (Foam::mag(overallVec.z()) > VSMALL)
275 {
276 scalar s = (min().z() - overallStart.z())/overallVec.z();
277 pt.x() = overallStart.x() + overallVec.x()*s;
278 pt.y() = overallStart.y() + overallVec.y()*s;
279 pt.z() = min().z();
280 }
281 else
282 {
283 pt.z() = min().z();
284 }
285 }
286 else if (ptBits & FRONTBIT)
287 {
288 // Intersect with plane V=max, n=0,0,1
289 if (Foam::mag(overallVec.z()) > VSMALL)
290 {
291 scalar s = (max().z() - overallStart.z())/overallVec.z();
292 pt.x() = overallStart.x() + overallVec.x()*s;
293 pt.y() = overallStart.y() + overallVec.y()*s;
294 pt.z() = max().z();
295 }
296 else
297 {
298 pt.z() = max().z();
299 }
300 }
301 }
302
303 // Can end up here if the end point is on the edge of the boundBox
304 return true;
305}
306
307
309(
310 const point& start,
311 const point& end,
312 point& pt
313) const
314{
315 direction ptBits;
316 return intersects(start, end-start, start, end, pt, ptBits);
317}
318
319
320bool Foam::treeBoundBox::contains(const vector& dir, const point& pt) const
321{
322 // Compare all components against min and max of bb
323
324 for (direction cmpt=0; cmpt < point::nComponents; ++cmpt)
325 {
326 if (pt[cmpt] < min()[cmpt])
327 {
328 return false;
329 }
330 else if (pt[cmpt] == min()[cmpt])
331 {
332 // On edge. Outside if direction points outwards.
333 if (dir[cmpt] < 0)
334 {
335 return false;
336 }
337 }
338
339 if (pt[cmpt] > max()[cmpt])
340 {
341 return false;
342 }
343 else if (pt[cmpt] == max()[cmpt])
344 {
345 // On edge. Outside if direction points outwards.
346 if (dir[cmpt] > 0)
347 {
348 return false;
349 }
350 }
351 }
352
353 // All components inside bb
354 return true;
355}
356
357
359{
360 direction octant = 0;
361
362 if (pt.x() == min().x())
363 {
364 octant |= LEFTBIT;
365 }
366 else if (pt.x() == max().x())
367 {
368 octant |= RIGHTBIT;
369 }
370
371 if (pt.y() == min().y())
372 {
373 octant |= BOTTOMBIT;
374 }
375 else if (pt.y() == max().y())
376 {
377 octant |= TOPBIT;
378 }
379
380 if (pt.z() == min().z())
381 {
382 octant |= BACKBIT;
383 }
384 else if (pt.z() == max().z())
385 {
386 octant |= FRONTBIT;
387 }
388
389 return octant;
390}
391
392
394{
395 direction octant = 0;
396
397 if (pt.x() < min().x())
398 {
399 octant |= LEFTBIT;
400 }
401 else if (pt.x() > max().x())
402 {
403 octant |= RIGHTBIT;
404 }
405
406 if (pt.y() < min().y())
407 {
408 octant |= BOTTOMBIT;
409 }
410 else if (pt.y() > max().y())
411 {
412 octant |= TOPBIT;
413 }
414
415 if (pt.z() < min().z())
416 {
417 octant |= BACKBIT;
418 }
419 else if (pt.z() > max().z())
420 {
421 octant |= FRONTBIT;
422 }
423
424 return octant;
425}
426
427
429(
430 const point& pt,
431 point& nearest,
432 point& furthest
433) const
434{
435 for (direction cmpt=0; cmpt < point::nComponents; ++cmpt)
436 {
437 if
438 (
439 Foam::mag(min()[cmpt] - pt[cmpt])
440 < Foam::mag(max()[cmpt] - pt[cmpt])
441 )
442 {
443 nearest[cmpt] = min()[cmpt];
444 furthest[cmpt] = max()[cmpt];
445 }
446 else
447 {
448 nearest[cmpt] = max()[cmpt];
449 furthest[cmpt] = min()[cmpt];
450 }
451 }
452}
453
454
455Foam::scalar Foam::treeBoundBox::maxDist(const point& pt) const
456{
457 point near, far;
458 calcExtremities(pt, near, far);
459
460 return Foam::mag(far - pt);
461}
462
463
465(
466 const point& pt,
467 const treeBoundBox& other
468) const
469{
470 //
471 // Distance point <-> nearest and furthest away vertex of this
472 //
473
474 point nearThis, farThis;
475
476 // get nearest and furthest away vertex
477 calcExtremities(pt, nearThis, farThis);
478
479 const scalar minDistThis =
480 sqr(nearThis.x() - pt.x())
481 + sqr(nearThis.y() - pt.y())
482 + sqr(nearThis.z() - pt.z());
483 const scalar maxDistThis =
484 sqr(farThis.x() - pt.x())
485 + sqr(farThis.y() - pt.y())
486 + sqr(farThis.z() - pt.z());
487
488 //
489 // Distance point <-> other
490 //
491
492 point nearOther, farOther;
493
494 // get nearest and furthest away vertex
495 other.calcExtremities(pt, nearOther, farOther);
496
497 const scalar minDistOther =
498 sqr(nearOther.x() - pt.x())
499 + sqr(nearOther.y() - pt.y())
500 + sqr(nearOther.z() - pt.z());
501 const scalar maxDistOther =
502 sqr(farOther.x() - pt.x())
503 + sqr(farOther.y() - pt.y())
504 + sqr(farOther.z() - pt.z());
505
506 //
507 // Categorize
508 //
509 if (maxDistThis < minDistOther)
510 {
511 // All vertices of this are nearer to point than any vertex of other
512 return -1;
513 }
514 else if (minDistThis > maxDistOther)
515 {
516 // All vertices of this are further from point than any vertex of other
517 return 1;
518 }
519 else
520 {
521 // Mixed bag
522 return 0;
523 }
524}
525
526
527// * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
528
530{
531 return os << static_cast<const boundBox&>(bb);
532}
533
534
536{
537 return is >> static_cast<boundBox&>(bb);
538}
539
540
541// ************************************************************************* //
scalar y
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:94
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
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
static constexpr direction nComponents
Number of components in bool is 1.
Definition: bool.H:98
A class for managing temporary objects.
Definition: tmp.H:65
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:89
void calcExtremities(const point &pt, point &nearest, point &furthest) const
Calculate nearest and furthest (to point) vertex coords of.
Definition: treeBoundBox.C:429
label distanceCmp(const point &pt, const treeBoundBox &other) const
Compare distance to point with other bounding box.
Definition: treeBoundBox.C:465
static const edgeList edges
Edge to point addressing.
Definition: treeBoundBox.H:154
treeBoundBox subBbox(const direction octant) const
Sub-box of given octant. Midpoint calculated.
Definition: treeBoundBox.C:106
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:162
direction faceBits(const point &pt) const
Code position of point on bounding box faces.
Definition: treeBoundBox.C:358
bool contains(const vector &dir, const point &) const
Contains point (inside or on edge) and moving in direction.
Definition: treeBoundBox.C:320
direction posBits(const point &pt) const
Position of point relative to bounding box.
Definition: treeBoundBox.C:393
@ 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
static const faceList faces
Face to point addressing.
Definition: treeBoundBox.H:151
treeBoundBox()
Construct without any points - an inverted bounding box.
Definition: treeBoundBoxI.H:34
scalar maxDist(const point &pt) const
Returns distance point to furthest away corner.
Definition: treeBoundBox.C:455
tmp< pointField > points() const
Vertex coordinates. In octant coding.
Definition: treeBoundBox.C:92
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
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))
#define WarningInFunction
Report a warning using Foam::Warning.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Istream & operator>>(Istream &, directionInfo &)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
errorManip< error > abort(error &err)
Definition: errorManip.H:144
uint8_t direction
Definition: direction.H:56
error FatalError
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333