Loading...
Searching...
No Matches
boundBoxI.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-2016 OpenFOAM Foundation
9 Copyright (C) 2016-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 "boundBox.H"
30
31// * * * * * * * * * * * * * Geometrical Information * * * * * * * * * * * * //
33namespace Foam
34{
35
36// Box corners as per hex cellmodel
37
38template<>
40{
41 // == octCorner<0>()
42 return min_;
43}
44
45template<>
47{
48 // == octCorner<1>()
49 return point(max_.x(), min_.y(), min_.z());
50}
51
52template<>
54{
55 // == octCorner<3>()
56 return point(max_.x(), max_.y(), min_.z());
57}
58
59template<>
61{
62 // == octCorner<2>()
63 return point(min_.x(), max_.y(), min_.z());
64}
65
66template<>
68{
69 // == octCorner<4>()
70 return point(min_.x(), min_.y(), max_.z());
71}
72
73template<>
75{
76 // == octCorner<5>()
77 return point(max_.x(), min_.y(), max_.z());
78}
79
80template<>
82{
83 // == octCorner<7>()
84 return max_;
85}
86
87template<>
89{
90 // == octCorner<6>()
91 return point(min_.x(), max_.y(), max_.z());
93
94} // End namespace Foam
95
96
97// Non-specialized version is compile-time disabled
98template<Foam::direction CornerNumber>
100{
101 static_assert(CornerNumber < 8, "Corner index [0..7]");
102 return point();
103}
104
105
106// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
107
109:
110 min_(invertedBox.min()),
111 max_(invertedBox.max())
112{}
113
114
116:
117 min_(point::zero),
118 max_(point::one)
119{}
120
121
123:
124 min_(p),
125 max_(p)
126{}
127
128
130:
131 min_(min),
132 max_(max)
133{}
134
135
137:
138 min_(bb.first()),
139 max_(bb.second())
140{}
141
142
145 operator>>(is, *this);
146}
147
148
149// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150
151inline bool Foam::boundBox::empty() const
152{
153 // Is empty/invalid if any component has (max < min), ie, !(min <= max)
154 return
155 (
156 (max_.x() < min_.x())
157 || (max_.y() < min_.y())
158 || (max_.z() < min_.z())
159 );
160}
161
163inline bool Foam::boundBox::good() const
164{
165 return !empty();
166}
167
169inline const Foam::point& Foam::boundBox::min() const noexcept
170{
171 return min_;
172}
173
175inline const Foam::point& Foam::boundBox::max() const noexcept
176{
177 return max_;
178}
179
182{
183 return min_;
184}
185
188{
189 return max_;
190}
191
194{
195 return 0.5 * (min_ + max_);
196}
197
200{
201 return (max_ - min_);
202}
203
205inline Foam::scalar Foam::boundBox::mag() const
206{
207 return min_.dist(max_);
208}
209
211inline Foam::scalar Foam::boundBox::magSqr() const
212{
213 return min_.distSqr(max_);
214}
215
217inline Foam::scalar Foam::boundBox::volume() const
218{
219 return cmptProduct(span());
220}
221
223inline Foam::scalar Foam::boundBox::minDim() const
224{
225 return cmptMin(span());
226}
227
229inline Foam::scalar Foam::boundBox::maxDim() const
230{
231 return cmptMax(span());
232}
233
235inline Foam::scalar Foam::boundBox::avgDim() const
236{
237 return cmptAv(span());
238}
239
240
242{
243 direction cmpt = 0;
244
245 scalar best = ROOTVGREAT;
246
247 for (direction dir = 0; dir < vector::nComponents; ++dir)
248 {
249 const scalar dist = (max_[dir] - min_[dir]);
250 if (dist < best && dist > 0)
251 {
252 best = dist;
253 cmpt = dir;
255 }
256
257 return cmpt;
258}
259
260
262{
263 direction cmpt = 0;
264
265 scalar best = 0;
266
267 for (direction dir = 0; dir < vector::nComponents; ++dir)
268 {
269 const scalar dist = (max_[dir] - min_[dir]);
270 if (dist > best)
271 {
272 best = dist;
273 cmpt = dir;
275 }
276
277 return cmpt;
278}
279
280
281inline int Foam::boundBox::nDim() const
282{
283 int ngood = 0;
284
285 for (direction dir = 0; dir < vector::nComponents; ++dir)
286 {
287 const scalar dist = (max_[dir] - min_[dir]);
288 if (dist < 0)
289 {
290 return -1;
291 }
292 else if (dist > 0)
293 {
294 ++ngood;
296 }
297
298 return ngood;
299}
300
301
303{
304 min_ = invertedBox.min();
305 max_ = invertedBox.max();
306}
307
308
310{
311 min_ = point::zero;
312 max_ = point::one;
313}
314
315
316inline void Foam::boundBox::reset(const point& pt)
317{
318 min_ = pt;
319 max_ = pt;
320}
321
322
323inline void Foam::boundBox::reset(const point& min, const point& max)
324{
325 min_ = min;
326 max_ = max;
327}
328
329
330inline void Foam::boundBox::add(const boundBox& bb)
331{
332 min_ = ::Foam::min(min_, bb.min_);
333 max_ = ::Foam::max(max_, bb.max_);
334}
335
336
337inline void Foam::boundBox::add(const point& pt)
338{
339 min_ = ::Foam::min(min_, pt);
340 max_ = ::Foam::max(max_, pt);
341}
342
343
344inline void Foam::boundBox::add(const point& pt0, const point& pt1)
345{
346 add(pt0);
347 add(pt1);
348}
349
350
352{
353 add(points.first());
354 add(points.second());
355}
356
357
358inline void Foam::boundBox::add(const UList<point>& points)
359{
360 for (const point& p : points)
361 {
362 add(p);
363 }
364}
365
366
367inline void Foam::boundBox::add(const tmp<pointField>& tpoints)
368{
369 add(tpoints());
370 tpoints.clear();
371}
372
373
374inline void Foam::boundBox::grow(const scalar delta)
375{
376 min_.x() -= delta; min_.y() -= delta; min_.z() -= delta;
377 max_.x() += delta; max_.y() += delta; max_.z() += delta;
378}
379
380
382{
383 min_ -= delta;
384 max_ += delta;
385}
386
388inline void Foam::boundBox::inflate(const scalar factor)
389{
390 grow(factor*mag());
391}
392
393
395(
396 const point& minA, const point& maxA, // boxA
397 const point& minB, const point& maxB // boxB
398)
399{
400 return
401 (
402 minA.x() <= maxB.x() && minB.x() <= maxA.x()
403 && minA.y() <= maxB.y() && minB.y() <= maxA.y()
404 && minA.z() <= maxB.z() && minB.z() <= maxA.z()
405 );
406}
407
408
410(
411 const point& corner0,
412 const point& corner1,
413 const point& centre,
414 const scalar radiusSqr
415)
416{
417 // Find out where centre is in relation to bb.
418 // Find nearest point on bb.
419 scalar distSqr = 0;
420
421 for (direction dir = 0; dir < vector::nComponents; ++dir)
422 {
423 const scalar d0 = corner0[dir] - centre[dir];
424 const scalar d1 = corner1[dir] - centre[dir];
425
426 if ((d0 > 0) != (d1 > 0))
427 {
428 // centre inside both extrema. This component does not add any
429 // distance.
430 }
431 else
432 {
433 distSqr += ::Foam::min(Foam::magSqr(d0), Foam::magSqr(d1));
434
435 if (distSqr > radiusSqr)
436 {
437 return false;
438 }
440 }
441
442 return true;
443}
444
446inline bool Foam::boundBox::overlaps(const boundBox& bb) const
447{
448 return box_box_overlaps(min_, max_, bb.min(), bb.max());
449}
450
451
452inline bool Foam::boundBox::overlaps
453(
454 const point& centre,
455 const scalar radiusSqr
456) const
457{
458 return box_sphere_overlaps(min_, max_, centre, radiusSqr);
459}
460
461
462inline bool Foam::boundBox::contains(const point& pt) const
463{
464 return
465 (
466 min_.x() <= pt.x() && pt.x() <= max_.x()
467 && min_.y() <= pt.y() && pt.y() <= max_.y()
468 && min_.z() <= pt.z() && pt.z() <= max_.z()
469 );
470}
471
473inline bool Foam::boundBox::contains(const boundBox& bb) const
474{
475 return contains(bb.min()) && contains(bb.max());
476}
477
478
479inline bool Foam::boundBox::containsInside(const point& pt) const
480{
481 return
482 (
483 min_.x() < pt.x() && pt.x() < max_.x()
484 && min_.y() < pt.y() && pt.y() < max_.y()
485 && min_.z() < pt.z() && pt.z() < max_.z()
486 );
487}
488
489
490// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
scalar delta
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
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:89
static const Form zero
static const Form one
static constexpr direction nComponents
Number of components in this vector space.
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector.H:135
scalar dist(const Vector< Cmpt > &v2) const
The L2-norm distance from another vector. The mag() of the difference.
Definition VectorI.H:107
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
bool containsInside(const point &pt) const
Contains point? (inside only).
Definition boundBoxI.H:472
int nDim() const
Count the number of positive, non-zero dimensions.
Definition boundBoxI.H:274
const point & max() const noexcept
Maximum describing the bounding box.
Definition boundBoxI.H:168
void grow(const scalar delta)
Expand box by adjusting min/max by specified amount in each dimension.
Definition boundBoxI.H:367
bool overlaps(const boundBox &bb) const
Overlaps/touches boundingBox?
Definition boundBoxI.H:439
static const boundBox invertedBox
A large inverted boundBox: min/max == +/- ROOTVGREAT.
Definition boundBox.H:131
tmp< pointField > points() const
Corner points in an order corresponding to a 'hex' cell.
Definition boundBox.H:381
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
direction minDir() const
Direction (X/Y/Z) of the smallest span (for empty box: 0).
Definition boundBoxI.H:234
const point & min() const noexcept
Minimum describing the bounding box.
Definition boundBoxI.H:162
bool contains(const point &pt) const
Contains point? (inside or on edge).
Definition boundBoxI.H:455
scalar volume() const
The volume of the bound box.
Definition boundBoxI.H:210
bool empty() const
Bounding box is inverted, contains no points.
Definition boundBoxI.H:144
void add(const boundBox &bb)
Extend to include the second box.
Definition boundBoxI.H:323
direction maxDir() const
Direction (X/Y/Z) of the largest span (for empty box: 0).
Definition boundBoxI.H:254
scalar minDim() const
Smallest length/height/width dimension.
Definition boundBoxI.H:216
scalar mag() const
The magnitude/length of the bounding box diagonal.
Definition boundBoxI.H:198
scalar avgDim() const
Average length/height/width dimension.
Definition boundBoxI.H:228
boundBox()
Default construct: an inverted bounding box.
Definition boundBoxI.H:101
bool good() const
Bounding box is non-inverted.
Definition boundBoxI.H:156
void reset()
Reset to an inverted box.
Definition boundBoxI.H:295
scalar magSqr() const
The magnitude/length squared of bounding box diagonal.
Definition boundBoxI.H:204
void inflate(const scalar factor)
Expand box by factor*mag(span) in all dimensions.
Definition boundBoxI.H:381
point hexCorner() const
Return corner point [0..7] corresponding to a 'hex' cell.
vector span() const
The bounding box span (from minimum to maximum).
Definition boundBoxI.H:192
static bool box_box_overlaps(const point &minA, const point &maxA, const point &minB, const point &maxB)
Test for overlap of box and box (inclusive check).
Definition boundBoxI.H:388
point centre() const
The centre (midpoint) of the bounding box.
Definition boundBoxI.H:186
scalar maxDim() const
Largest length/height/width dimension.
Definition boundBoxI.H:222
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition one.H:57
A class for managing temporary objects.
Definition tmp.H:75
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition tmpI.H:289
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
Definition pTraits.H:51
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
volScalarField & p
const pointField & points
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
Cmpt cmptProduct(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
void add(DimensionedField< scalar, GeoMesh > &result, const dimensioned< scalar > &dt1, const DimensionedField< scalar, GeoMesh > &f2)
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh > > cmptAv(const DimensionedField< Type, GeoMesh > &f1)
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:26
uint8_t direction
Definition direction.H:49
vector point
Point is a vector.
Definition point.H:37
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
const direction noexcept
Definition scalarImpl.H:265
Vector< scalar > vector
Definition vector.H:57
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)