All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
boundBoxI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "boundBox.H"
27 #include "pointField.H"
28 
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
33 :
34  min_(Zero),
35  max_(Zero)
36 {}
37 
38 
39 inline Foam::boundBox::boundBox(const point& min, const point& max)
40 :
41  min_(min),
42  max_(max)
43 {}
44 
45 
47 {
48  operator>>(is, *this);
49 }
50 
51 
52 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
53 
54 inline const Foam::point& Foam::boundBox::min() const
55 {
56  return min_;
57 }
58 
59 
60 inline const Foam::point& Foam::boundBox::max() const
61 {
62  return max_;
63 }
64 
65 
67 {
68  return min_;
69 }
70 
71 
73 {
74  return max_;
75 }
76 
77 
79 {
80  return 0.5 * (max_ + min_);
81 }
82 
83 
85 {
86  return (max_ - min_);
87 }
88 
89 
90 inline Foam::scalar Foam::boundBox::mag() const
91 {
92  return ::Foam::mag(max_ - min_);
93 }
94 
95 
96 inline Foam::scalar Foam::boundBox::volume() const
97 {
98  return cmptProduct(span());
99 }
100 
101 
102 inline Foam::scalar Foam::boundBox::minDim() const
103 {
104  return cmptMin(span());
105 }
106 
107 
108 inline Foam::scalar Foam::boundBox::maxDim() const
109 {
110  return cmptMax(span());
111 }
112 
113 
114 inline Foam::scalar Foam::boundBox::avgDim() const
115 {
116  return cmptAv(span());
117 }
118 
119 
120 inline bool Foam::boundBox::overlaps(const boundBox& bb) const
121 {
122  return
123  (
124  bb.max_.x() >= min_.x() && bb.min_.x() <= max_.x()
125  && bb.max_.y() >= min_.y() && bb.min_.y() <= max_.y()
126  && bb.max_.z() >= min_.z() && bb.min_.z() <= max_.z()
127  );
128 }
129 
130 
132 (
133  const point& centre,
134  const scalar radiusSqr
135 ) const
136 {
137  // Find out where centre is in relation to bb.
138  // Find nearest point on bb.
139  scalar distSqr = 0;
140 
141  for (direction dir = 0; dir < vector::nComponents; dir++)
142  {
143  scalar d0 = min_[dir] - centre[dir];
144  scalar d1 = max_[dir] - centre[dir];
145 
146  if ((d0 > 0) != (d1 > 0))
147  {
148  // centre inside both extrema. This component does not add any
149  // distance.
150  }
151  else if (Foam::mag(d0) < Foam::mag(d1))
152  {
153  distSqr += d0*d0;
154  }
155  else
156  {
157  distSqr += d1*d1;
158  }
159 
160  if (distSqr > radiusSqr)
161  {
162  return false;
163  }
164  }
165 
166  return true;
167 }
168 
169 
170 inline bool Foam::boundBox::contains(const point& pt) const
171 {
172  return
173  (
174  pt.x() >= min_.x() && pt.x() <= max_.x()
175  && pt.y() >= min_.y() && pt.y() <= max_.y()
176  && pt.z() >= min_.z() && pt.z() <= max_.z()
177  );
178 }
179 
180 
181 // this.bb fully contains bb
182 inline bool Foam::boundBox::contains(const boundBox& bb) const
183 {
184  return contains(bb.min()) && contains(bb.max());
185 }
186 
187 
188 inline bool Foam::boundBox::containsInside(const point& pt) const
189 {
190  return
191  (
192  pt.x() > min_.x() && pt.x() < max_.x()
193  && pt.y() > min_.y() && pt.y() < max_.y()
194  && pt.z() > min_.z() && pt.z() < max_.z()
195  );
196 }
197 
198 
199 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
200 
201 inline bool Foam::operator==(const boundBox& a, const boundBox& b)
202 {
203  return (a.min_ == b.min_) && (a.max_ == b.max_);
204 }
205 
206 
207 inline bool Foam::operator!=(const boundBox& a, const boundBox& b)
208 {
209  return !(a == b);
210 }
211 
212 
213 // ************************************************************************* //
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
static const direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:99
const Cmpt & z() const
Definition: VectorI.H:87
const Cmpt & y() const
Definition: VectorI.H:81
const Cmpt & x() const
Definition: VectorI.H:75
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:59
const point & min() const
Minimum point defining the bounding box.
Definition: boundBoxI.H:54
const point & max() const
Maximum point defining the bounding box.
Definition: boundBoxI.H:60
scalar volume() const
The volume of the bound box.
Definition: boundBoxI.H:96
bool overlaps(const boundBox &) const
Overlaps/touches boundingBox?
Definition: boundBoxI.H:120
bool contains(const point &) const
Contains point? (inside or on edge)
Definition: boundBoxI.H:170
scalar minDim() const
Smallest length/height/width dimension.
Definition: boundBoxI.H:102
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:90
scalar avgDim() const
Average length/height/width dimension.
Definition: boundBoxI.H:114
boundBox()
Construct null, setting points to zero.
Definition: boundBoxI.H:32
point midpoint() const
The midpoint of the bounding box.
Definition: boundBoxI.H:78
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:84
scalar maxDim() const
Largest length/height/width dimension.
Definition: boundBoxI.H:108
bool containsInside(const point &) const
Contains point? (inside only)
Definition: boundBoxI.H:188
volScalarField & b
Definition: createFields.H:27
static const zero Zero
Definition: zero.H:97
bool operator!=(const particle &, const particle &)
Definition: particle.C:1257
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Cmpt cmptProduct(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:521
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh >> cmptAv(const DimensionedField< Type, GeoMesh > &df)
Istream & operator>>(Istream &, directionInfo &)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
uint8_t direction
Definition: direction.H:45