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-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 "boundBox.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
34:
35 min_(invertedBox.min()),
36 max_(invertedBox.max())
37{}
38
39
41:
42 min_(pt),
43 max_(pt)
44{}
45
46
48:
49 min_(min),
50 max_(max)
51{}
52
53
55{
56 operator>>(is, *this);
57}
58
59
60// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
61
62inline bool Foam::boundBox::empty() const
63{
64 // Check each component for max < min
65 for (direction dir = 0; dir < vector::nComponents; ++dir)
66 {
67 if (max_[dir] < min_[dir])
68 {
69 return true;
70 }
71 }
72 return false;
73}
74
75
76inline bool Foam::boundBox::valid() const
77{
78 // Check each component for max < min
79 for (direction dir = 0; dir < vector::nComponents; ++dir)
80 {
81 if (max_[dir] < min_[dir])
82 {
83 return false;
84 }
85 }
86
87 return true;
88}
89
90
91inline const Foam::point& Foam::boundBox::min() const
92{
93 return min_;
94}
95
96
97inline const Foam::point& Foam::boundBox::max() const
98{
99 return max_;
100}
101
102
104{
105 return min_;
106}
107
108
110{
111 return max_;
112}
113
114
116{
117 return 0.5 * (min_ + max_);
118}
119
120
122{
123 return this->centre();
124}
125
126
128{
129 return (max_ - min_);
130}
131
132
133inline Foam::scalar Foam::boundBox::mag() const
134{
135 return ::Foam::mag(max_ - min_);
136}
137
138
139inline Foam::scalar Foam::boundBox::volume() const
140{
141 return cmptProduct(span());
142}
143
144
145inline Foam::scalar Foam::boundBox::minDim() const
146{
147 return cmptMin(span());
148}
149
150
151inline Foam::scalar Foam::boundBox::maxDim() const
152{
153 return cmptMax(span());
154}
155
156
157inline Foam::scalar Foam::boundBox::avgDim() const
158{
159 return cmptAv(span());
160}
161
162
163inline Foam::label Foam::boundBox::nDim() const
164{
165 label ngood = 0;
166
167 for (direction dir = 0; dir < vector::nComponents; ++dir)
168 {
169 const scalar diff = (max_[dir] - min_[dir]);
170 if (diff < 0)
171 {
172 return -1;
173 }
174 else if (diff > 0)
175 {
176 ++ngood;
177 }
178 }
179
180 return ngood;
181}
182
183
185{
186 min_ = invertedBox.min();
187 max_ = invertedBox.max();
188}
189
190
191inline void Foam::boundBox::add(const boundBox& bb)
192{
193 min_ = ::Foam::min(min_, bb.min_);
194 max_ = ::Foam::max(max_, bb.max_);
195}
196
197
198inline void Foam::boundBox::add(const point& pt)
199{
200 min_ = ::Foam::min(min_, pt);
201 max_ = ::Foam::max(max_, pt);
202}
203
204
206{
207 for (const point& p : points)
208 {
209 add(p);
210 }
211}
212
213
214inline void Foam::boundBox::add(const tmp<pointField>& tpoints)
215{
216 add(tpoints());
217 tpoints.clear();
218}
219
220
221inline bool Foam::boundBox::overlaps(const boundBox& bb) const
222{
223 return
224 (
225 bb.max_.x() >= min_.x() && bb.min_.x() <= max_.x()
226 && bb.max_.y() >= min_.y() && bb.min_.y() <= max_.y()
227 && bb.max_.z() >= min_.z() && bb.min_.z() <= max_.z()
228 );
229}
230
231
233(
234 const point& centre,
235 const scalar radiusSqr
236) const
237{
238 // Find out where centre is in relation to bb.
239 // Find nearest point on bb.
240 scalar distSqr = 0;
241
242 for (direction dir = 0; dir < vector::nComponents; ++dir)
243 {
244 const scalar d0 = min_[dir] - centre[dir];
245 const scalar d1 = max_[dir] - centre[dir];
246
247 if ((d0 > 0) != (d1 > 0))
248 {
249 // centre inside both extrema. This component does not add any
250 // distance.
251 }
252 else if (Foam::mag(d0) < Foam::mag(d1))
253 {
254 distSqr += d0*d0;
255 }
256 else
257 {
258 distSqr += d1*d1;
259 }
260
261 if (distSqr > radiusSqr)
262 {
263 return false;
265 }
266
267 return true;
268}
269
270
271inline bool Foam::boundBox::contains(const point& pt) const
272{
273 return
275 min_.x() <= pt.x() && pt.x() <= max_.x()
276 && min_.y() <= pt.y() && pt.y() <= max_.y()
277 && min_.z() <= pt.z() && pt.z() <= max_.z()
278 );
279}
280
281
282inline bool Foam::boundBox::contains(const boundBox& bb) const
283{
284 return contains(bb.min()) && contains(bb.max());
285}
286
287
288inline bool Foam::boundBox::containsInside(const point& pt) const
289{
290 return
291 (
292 min_.x() < pt.x() && pt.x() < max_.x()
293 && min_.y() < pt.y() && pt.y() < max_.y()
294 && min_.z() < pt.z() && pt.z() < max_.z()
295 );
296}
297
298
299// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
300
302{
303 add(bb);
304}
305
306
307// * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * * //
308
309inline bool Foam::operator==(const boundBox& a, const boundBox& b)
310{
311 return (a.min() == b.min()) && (a.max() == b.max());
312}
313
314
315inline bool Foam::operator!=(const boundBox& a, const boundBox& b)
316{
317 return !(a == b);
318}
319
320
321// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Y[inertIndex] max(0.0)
void max(const dimensioned< Type > &dt)
Use the maximum of the field and specified value.
void min(const dimensioned< Type > &dt)
Use the minimum of the field and specified value.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
int overlaps
Flag to control which overlap calculations are performed.
Definition: PDRparams.H:97
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
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
bool containsInside(const point &pt) const
Contains point? (inside only)
Definition: boundBoxI.H:288
bool valid() const
Bounding box is non-inverted.
Definition: boundBoxI.H:76
void operator+=(const boundBox &bb)
Extend box to include the second box, as per the add() method.
Definition: boundBoxI.H:301
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
bool contains(const point &pt) const
Contains point? (inside or on edge)
Definition: boundBoxI.H:271
scalar volume() const
The volume of the bound box.
Definition: boundBoxI.H:139
bool empty() const
Bounding box is inverted, contains no points.
Definition: boundBoxI.H:62
scalar minDim() const
Smallest length/height/width dimension.
Definition: boundBoxI.H:145
label nDim() const
Count the number of positive, non-zero dimensions.
Definition: boundBoxI.H:163
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:133
scalar avgDim() const
Average length/height/width dimension.
Definition: boundBoxI.H:157
boundBox()
Construct without any points - an inverted bounding box.
Definition: boundBoxI.H:33
void clear()
Clear bounding box and make it an inverted box.
Definition: boundBoxI.H:184
point midpoint() const
The midpoint (centre) of the bounding box. Identical to centre()
Definition: boundBoxI.H:121
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:127
point centre() const
The centre (midpoint) of the bounding box.
Definition: boundBoxI.H:115
scalar maxDim() const
Largest length/height/width dimension.
Definition: boundBoxI.H:151
Sums a given list of (at least two or more) fields and outputs the result into a new field,...
Definition: add.H:161
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
void clear() const noexcept
Definition: tmpI.H:287
volScalarField & p
const pointField & points
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
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
Cmpt cmptProduct(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:598
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
Istream & operator>>(Istream &, directionInfo &)
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:378
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
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh > > cmptAv(const DimensionedField< Type, GeoMesh > &df)
volScalarField & b
Definition: createFields.H:27