boundBox.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) 2016-2019 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 #include "PstreamReduceOps.H"
31 #include "tmp.H"
32 #include "plane.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
37 (
38  point::uniform(-ROOTVGREAT),
39  point::uniform(ROOTVGREAT)
40 );
41 
43 (
44  point::uniform(ROOTVGREAT),
45  point::uniform(-ROOTVGREAT)
46 );
47 
49 ({
50  // Point and face order as per hex cellmodel
51  face({0, 4, 7, 3}), // 0: x-min, left
52  face({1, 2, 6, 5}), // 1: x-max, right
53  face({0, 1, 5, 4}), // 2: y-min, bottom
54  face({3, 7, 6, 2}), // 3: y-max, top
55  face({0, 3, 2, 1}), // 4: z-min, back
56  face({4, 5, 6, 7}) // 5: z-max, front
57 });
58 
60 ({
61  vector(-1, 0, 0), // 0: x-min, left
62  vector( 1, 0, 0), // 1: x-max, right
63  vector( 0, -1, 0), // 2: y-min, bottom
64  vector( 0, 1, 0), // 3: y-max, top
65  vector( 0, 0, -1), // 4: z-min, back
66  vector( 0, 0, 1) // 5: z-max, front
67 });
68 
69 
70 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
71 
73 :
74  boundBox()
75 {
76  add(points);
77 
78  if (doReduce)
79  {
80  reduce();
81  }
82 }
83 
84 
85 Foam::boundBox::boundBox(const tmp<pointField>& tpoints, bool doReduce)
86 :
87  boundBox()
88 {
89  add(tpoints);
90 
91  if (doReduce)
92  {
93  reduce();
94  }
95 }
96 
97 
99 (
100  const UList<point>& points,
101  const labelUList& indices,
102  bool doReduce
103 )
104 :
105  boundBox()
106 {
107  add(points, indices);
108 
109  if (doReduce)
110  {
111  reduce();
112  }
113 }
114 
115 
116 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
117 
119 {
120  auto tpt = tmp<pointField>::New(8);
121  auto& pt = tpt.ref();
122 
123  pt[0] = min_; // min-x, min-y, min-z
124  pt[1] = point(max_.x(), min_.y(), min_.z()); // max-x, min-y, min-z
125  pt[2] = point(max_.x(), max_.y(), min_.z()); // max-x, max-y, min-z
126  pt[3] = point(min_.x(), max_.y(), min_.z()); // min-x, max-y, min-z
127  pt[4] = point(min_.x(), min_.y(), max_.z()); // min-x, min-y, max-z
128  pt[5] = point(max_.x(), min_.y(), max_.z()); // max-x, min-y, max-z
129  pt[6] = max_; // max-x, max-y, max-z
130  pt[7] = point(min_.x(), max_.y(), max_.z()); // min-x, max-y, max-z
131 
132  return tpt;
133 }
134 
135 
137 {
138  auto tpts = tmp<pointField>::New(6);
139  auto& pts = tpts.ref();
140 
141  forAll(pts, facei)
142  {
143  pts[facei] = faceCentre(facei);
144  }
145 
146  return tpts;
147 }
148 
149 
151 {
152  point pt = boundBox::centre();
153 
154  if (facei > 5)
155  {
157  << "face should be [0..5]"
158  << abort(FatalError);
159  }
160 
161  switch (facei)
162  {
163  case 0: pt.x() = min().x(); break; // 0: x-min, left
164  case 1: pt.x() = max().x(); break; // 1: x-max, right
165  case 2: pt.y() = min().y(); break; // 2: y-min, bottom
166  case 3: pt.y() = max().y(); break; // 3: y-max, top
167  case 4: pt.z() = min().z(); break; // 4: z-min, back
168  case 5: pt.z() = max().z(); break; // 5: z-max, front
169  }
170 
171  return pt;
172 }
173 
174 
175 void Foam::boundBox::inflate(const scalar s)
176 {
177  const vector ext = vector::one*s*mag();
178 
179  min_ -= ext;
180  max_ += ext;
181 }
182 
183 
185 {
186  Foam::reduce(min_, minOp<point>());
187  Foam::reduce(max_, maxOp<point>());
188 }
189 
190 
192 {
193  min_ = ::Foam::max(min_, bb.min_);
194  max_ = ::Foam::min(max_, bb.max_);
195 
196  return valid();
197 }
198 
199 
200 bool Foam::boundBox::intersects(const plane& pln) const
201 {
202  // Require a full 3D box
203  if (nDim() != 3)
204  {
205  return false;
206  }
207 
208  bool above = false;
209  bool below = false;
210 
211  tmp<pointField> tpts(points());
212  const auto& pts = tpts();
213 
214  for (const point& p : pts)
215  {
216  if (pln.sideOfPlane(p) == plane::FRONT)
217  {
218  above = true;
219  }
220  else
221  {
222  below = true;
223  }
224  }
225 
226  return (above && below);
227 }
228 
229 
231 {
232  if (points.empty())
233  {
234  return true;
235  }
236 
237  for (const point& p : points)
238  {
239  if (!contains(p))
240  {
241  return false;
242  }
243  }
244 
245  return true;
246 }
247 
248 
250 {
251  if (points.empty())
252  {
253  return true;
254  }
255 
256  for (const point& p : points)
257  {
258  if (contains(p))
259  {
260  return true;
261  }
262  }
263 
264  return false;
265 }
266 
267 
269 {
270  // Clip the point to the range of the bounding box
271  const scalar surfPtx = Foam::max(Foam::min(pt.x(), max_.x()), min_.x());
272  const scalar surfPty = Foam::max(Foam::min(pt.y(), max_.y()), min_.y());
273  const scalar surfPtz = Foam::max(Foam::min(pt.z(), max_.z()), min_.z());
274 
275  return point(surfPtx, surfPty, surfPtz);
276 }
277 
278 
279 // * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
280 
282 {
283  if (os.format() == IOstream::ASCII)
284  {
285  os << bb.min_ << token::SPACE << bb.max_;
286  }
287  else
288  {
289  os.write
290  (
291  reinterpret_cast<const char*>(&bb.min_),
292  sizeof(boundBox)
293  );
294  }
295 
297  return os;
298 }
299 
300 
302 {
303  if (is.format() == IOstream::ASCII)
304  {
305  is >> bb.min_ >> bb.max_;
306  }
307  else
308  {
309  Detail::readContiguous<boundBox>
310  (
311  is,
312  reinterpret_cast<char*>(&bb.min_),
313  sizeof(boundBox)
314  );
315  }
316 
317  is.check(FUNCTION_NAME);
318  return is;
319 }
320 
321 
322 // ************************************************************************* //
Foam::maxOp
Definition: ops.H:223
Foam::Vector::x
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
Foam::boundBox::reduce
void reduce()
Parallel reduction of min/max values.
Definition: boundBox.C:184
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::one
static const Vector< scalar > one
Definition: VectorSpace.H:116
s
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))
Definition: gmvOutputSpray.H:25
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::boundBox::faces
static const faceList faces
Faces to point addressing, as per a 'hex' cell.
Definition: boundBox.H:89
Foam::minOp
Definition: ops.H:224
Foam::boundBox::faceCentre
point faceCentre(const direction facei) const
Face centre of given face index.
Definition: boundBox.C:150
Foam::boundBox::containsAny
bool containsAny(const UList< point > &points) const
Contains any of the points? (inside or on edge)
Definition: boundBox.C:249
Foam::boundBox::inflate
void inflate(const scalar s)
Inflate box by factor*mag(span) in all dimensions.
Definition: boundBox.C:175
Foam::boundBox::invertedBox
static const boundBox invertedBox
A large inverted boundBox: min/max == +/- ROOTVGREAT.
Definition: boundBox.H:86
Foam::IOstreamOption::format
streamFormat format() const noexcept
Get the current stream format.
Definition: IOstreamOption.H:286
Foam::operator>>
Istream & operator>>(Istream &, directionInfo &)
Definition: directionInfo.C:230
Foam::OBJstream::write
virtual Ostream & write(const char c)
Write character.
Definition: OBJstream.C:78
Foam::Vector::z
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::boundBox::intersect
bool intersect(const boundBox &bb)
Intersection (union) with the second box.
Definition: boundBox.C:191
Foam::plane
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition: plane.H:89
Foam::operator<<
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
plane.H
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::boundBox::centre
point centre() const
The centre (midpoint) of the bounding box.
Definition: boundBoxI.H:115
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::IOstream::check
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:58
Foam::FatalError
error FatalError
Foam::boundBox::boundBox
boundBox()
Construct without any points - an inverted bounding box.
Definition: boundBoxI.H:33
Foam::add
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
Definition: FieldFieldFunctions.C:939
os
OBJstream os(runTime.globalPath()/outputName)
Foam::boundBox::nearest
point nearest(const point &pt) const
Return the nearest point on the boundBox to the supplied point.
Definition: boundBox.C:268
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
PstreamReduceOps.H
Inter-processor communication reduction functions.
Foam::boundBox::contains
bool contains(const point &pt) const
Contains point? (inside or on edge)
Definition: boundBoxI.H:271
boundBox.H
Foam::Vector::y
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
tmp.H
Foam::IOstreamOption::ASCII
"ascii" (normal default)
Definition: IOstreamOption.H:72
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::boundBox::greatBox
static const boundBox greatBox
A large boundBox: min/max == -/+ ROOTVGREAT.
Definition: boundBox.H:83
Foam::Vector< scalar >
Foam::boundBox::intersects
bool intersects(const plane &pln) const
Does plane intersect this bounding box.
Definition: boundBox.C:200
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::FixedList
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:104
Foam::token::SPACE
Space [isspace].
Definition: token.H:125
Foam::UList
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:103
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::direction
uint8_t direction
Definition: direction.H:52
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:295
Foam::boundBox
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
Foam::boundBox::faceNormals
static const FixedList< vector, 6 > faceNormals
The unit normal per face.
Definition: boundBox.H:92
Foam::plane::FRONT
The front (positive normal) side of the plane.
Definition: plane.H:96
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::plane::sideOfPlane
side sideOfPlane(const point &p) const
Return the side of the plane that the point is on.
Definition: planeI.H:87
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::boundBox::points
tmp< pointField > points() const
Corner points in an order corresponding to a 'hex' cell.
Definition: boundBox.C:118
Foam::boundBox::faceCentres
tmp< pointField > faceCentres() const
Face midpoints.
Definition: boundBox.C:136
Foam::boundBox::add
void add(const boundBox &bb)
Extend to include the second box.
Definition: boundBoxI.H:191