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