boundBox.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 -------------------------------------------------------------------------------
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 Class
28  Foam::boundBox
29 
30 Description
31  A bounding box defined in terms of min/max extrema points.
32 
33 Note
34  When a bounding box is created without any points, it creates an inverted
35  bounding box. Points can be added later and the bounding box will grow to
36  include them.
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef boundBox_H
41 #define boundBox_H
42 
43 #include "pointField.H"
44 #include "faceList.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 // Forward Declarations
52 class boundBox;
53 class plane;
54 template<class T> class tmp;
55 
56 Istream& operator>>(Istream& is, boundBox& bb);
57 Ostream& operator<<(Ostream& os, const boundBox& bb);
58 
59 
60 /*---------------------------------------------------------------------------*\
61  Class boundBox Declaration
62 \*---------------------------------------------------------------------------*/
63 
64 class boundBox
65 {
66  // Private data
67 
68  //- Minimum and maximum points describing the bounding box
69  point min_, max_;
70 
71 public:
72 
73  // Static Data Members
74 
75  //- Bits used for (x/y/z) direction encoding.
76  enum directionBit : direction
77  {
78  XDIR = 0x1,
79  YDIR = 0x2,
80  ZDIR = 0x4
81  };
82 
83  //- A large boundBox: min/max == -/+ ROOTVGREAT
84  static const boundBox greatBox;
85 
86  //- A large inverted boundBox: min/max == +/- ROOTVGREAT
87  static const boundBox invertedBox;
88 
89  //- Faces to point addressing, as per a 'hex' cell
90  static const faceList faces;
91 
92  //- The unit normal per face
93  static const FixedList<vector, 6> faceNormals;
94 
95 
96  // Constructors
97 
98  //- Construct without any points - an inverted bounding box
99  inline boundBox();
100 
101  //- Construct a bounding box containing a single initial point
102  inline explicit boundBox(const point& pt);
103 
104  //- Construct from components
105  inline boundBox(const point& min, const point& max);
106 
107  //- Construct as the bounding box of the given points
108  // Does parallel communication (doReduce = true)
109  explicit boundBox(const UList<point>& points, bool doReduce = true);
110 
111  //- Construct as the bounding box of the given temporary pointField.
112  // Does parallel communication (doReduce = true)
113  explicit boundBox(const tmp<pointField>& tpoints, bool doReduce = true);
114 
115  //- Construct bounding box as an indirect subset of the points.
116  // The indices could be from cell/face etc.
117  // Does parallel communication (doReduce = true)
118  boundBox
119  (
120  const UList<point>& points,
121  const labelUList& indices,
122  bool doReduce = true
123  );
124 
125  //- Construct bounding box as an indirect subset of the points.
126  // The indices could be from edge/triFace etc.
127  // Does parallel communication (doReduce = true)
128  template<unsigned N>
129  boundBox
130  (
131  const UList<point>& points,
132  const FixedList<label, N>& indices,
133  bool doReduce = true
134  );
135 
136  //- Construct from Istream
137  inline boundBox(Istream& is);
138 
139 
140  // Member Functions
141 
142  // Access
143 
144  //- Bounding box is inverted, contains no points.
145  inline bool empty() const;
146 
147  //- Bounding box is non-inverted.
148  inline bool valid() const;
149 
150  //- Minimum describing the bounding box
151  inline const point& min() const;
152 
153  //- Maximum describing the bounding box
154  inline const point& max() const;
155 
156  //- Minimum describing the bounding box, non-const access
157  inline point& min();
158 
159  //- Maximum describing the bounding box, non-const access
160  inline point& max();
161 
162  //- The centre (midpoint) of the bounding box
163  inline point centre() const;
164 
165  //- The midpoint (centre) of the bounding box. Identical to centre()
166  inline point midpoint() const;
167 
168  //- The bounding box span (from minimum to maximum)
169  inline vector span() const;
170 
171  //- The magnitude of the bounding box span
172  inline scalar mag() const;
173 
174  //- The volume of the bound box
175  inline scalar volume() const;
176 
177  //- Smallest length/height/width dimension
178  inline scalar minDim() const;
179 
180  //- Largest length/height/width dimension
181  inline scalar maxDim() const;
182 
183  //- Average length/height/width dimension
184  inline scalar avgDim() const;
185 
186  //- Count the number of positive, non-zero dimensions.
187  // \return -1 if any dimensions are negative,
188  // 0 = 0D (point),
189  // 1 = 1D (line aligned with an axis),
190  // 2 = 2D (plane aligned with an axis),
191  // 3 = 3D (box)
192  inline label nDim() const;
193 
194  //- Corner points in an order corresponding to a 'hex' cell
195  tmp<pointField> points() const;
196 
197  //- Face midpoints
199 
200  //- Face centre of given face index
201  point faceCentre(const direction facei) const;
202 
203 
204  // Manipulate
205 
206  //- Clear bounding box and make it an inverted box
207  inline void clear();
208 
209  //- Extend to include the second box.
210  inline void add(const boundBox& bb);
211 
212  //- Extend to include the point.
213  inline void add(const point& pt);
214 
215  //- Extend to include the points.
216  inline void add(const UList<point>& points);
217 
218  //- Extend to include the points from the temporary point field.
219  inline void add(const tmp<pointField>& tpoints);
220 
221  //- Extend to include the points.
222  template<unsigned N>
223  void add(const FixedList<point, N>& points);
224 
225  //- Extend to include a (subsetted) point field.
226  // The indices could be from edge/triFace etc.
227  template<unsigned N>
228  void add
229  (
230  const UList<point>& points,
231  const FixedList<label, N>& indices
232  );
233 
234  //- Extend to include a (subsetted) point field.
235  //
236  // \tparam IntContainer A container with an iterator that
237  // dereferences to an label
238  template<class IntContainer>
239  void add
240  (
241  const UList<point>& points,
242  const IntContainer& indices
243  );
244 
245  //- Inflate box by factor*mag(span) in all dimensions
246  void inflate(const scalar s);
247 
248  //- Parallel reduction of min/max values
249  void reduce();
250 
251 
252  // Query
253 
254  //- Intersection (union) with the second box.
255  // The return value is true if the intersection is non-empty.
256  bool intersect(const boundBox& bb);
257 
258  //- Does plane intersect this bounding box.
259  // There is an intersection if the plane segments the corner points
260  // \note the results are unreliable when plane coincides almost
261  // exactly with a box face
262  bool intersects(const plane& pln) const;
263 
264  //- Overlaps/touches boundingBox?
265  inline bool overlaps(const boundBox& bb) const;
266 
267  //- Overlaps boundingSphere (centre + sqr(radius))?
268  inline bool overlaps
269  (
270  const point& centre,
271  const scalar radiusSqr
272  ) const;
273 
274  //- Contains point? (inside or on edge)
275  inline bool contains(const point& pt) const;
276 
277  //- Fully contains other boundingBox?
278  inline bool contains(const boundBox& bb) const;
279 
280  //- Contains point? (inside only)
281  inline bool containsInside(const point& pt) const;
282 
283  //- Contains all points? (inside or on edge)
284  bool contains(const UList<point>& points) const;
285 
286  //- Contains all of the (subsetted) points? (inside or on edge)
287  template<unsigned N>
288  bool contains
289  (
290  const UList<point>& points,
291  const FixedList<label, N>& indices
292  ) const;
293 
294 
295  //- Contains all of the (subsetted) points? (inside or on edge)
296  //
297  // \tparam IntContainer A container with an iterator that
298  // dereferences to an label
299  template<class IntContainer>
300  bool contains
301  (
302  const UList<point>& points,
303  const IntContainer& indices
304  ) const;
305 
306 
307  //- Contains any of the points? (inside or on edge)
308  bool containsAny(const UList<point>& points) const;
309 
310  //- Contains any of the (subsetted) points? (inside or on edge)
311  template<unsigned N>
312  bool containsAny
313  (
314  const UList<point>& points,
315  const FixedList<label, N>& indices
316  ) const;
317 
318  //- Contains any of the (subsetted) points? (inside or on edge)
319  //
320  // \tparam IntContainer A container with an iterator that
321  // dereferences to an label
322  template<class IntContainer>
323  bool containsAny
324  (
325  const UList<point>& points,
326  const IntContainer& indices
327  ) const;
328 
329 
330  //- Return the nearest point on the boundBox to the supplied point.
331  // If point is inside the boundBox then the point is returned
332  // unchanged.
333  point nearest(const point& pt) const;
334 
335 
336  // Member Operators
337 
338  //- Extend box to include the second box, as per the add() method.
339  inline void operator+=(const boundBox& bb);
340 
341 
342  // IOstream operator
343 
344  friend Istream& operator>>(Istream& is, boundBox& bb);
345  friend Ostream& operator<<(Ostream& os, const boundBox& bb);
346 };
347 
348 
349 // * * * * * * * * * * * * * * * * * Traits * * * * * * * * * * * * * * * * //
350 
351 //- Contiguous data for boundBox
352 template<> struct is_contiguous<boundBox> : is_contiguous<point> {};
353 
354 //- Contiguous scalar data for boundBox
355 template<> struct is_contiguous_scalar<boundBox>
356 :
357  is_contiguous_scalar<point>
358 {};
359 
360 
361 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
362 
363 inline bool operator==(const boundBox& a, const boundBox& b);
364 inline bool operator!=(const boundBox& a, const boundBox& b);
365 
366 
367 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
368 
369 } // End namespace Foam
370 
371 #include "boundBoxI.H"
372 
373 #ifdef NoRepository
374  #include "boundBoxTemplates.C"
375 #endif
376 
377 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
378 
379 #endif
380 
381 // ************************************************************************* //
Foam::boundBox::reduce
void reduce()
Parallel reduction of min/max values.
Definition: boundBox.C:184
Foam::boundBox::midpoint
point midpoint() const
The midpoint (centre) of the bounding box. Identical to centre()
Definition: boundBoxI.H:121
Foam::boundBox::mag
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:133
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::boundBox::XDIR
1: x-direction (vector component 0)
Definition: boundBox.H:77
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::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::boundBox::YDIR
2: y-direction (vector component 1)
Definition: boundBox.H:78
Foam::operator>>
Istream & operator>>(Istream &, directionInfo &)
Definition: directionInfo.C:230
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
Foam::boundBox::volume
scalar volume() const
The volume of the bound box.
Definition: boundBoxI.H:139
faceList.H
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
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::boundBox::overlaps
bool overlaps(const boundBox &bb) const
Overlaps/touches boundingBox?
Definition: boundBoxI.H:221
Foam::boundBox::span
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:127
Foam::operator!=
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:239
Foam::operator<<
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
Foam::boundBox::nDim
label nDim() const
Count the number of positive, non-zero dimensions.
Definition: boundBoxI.H:163
Foam::boundBox::maxDim
scalar maxDim() const
Largest length/height/width dimension.
Definition: boundBoxI.H:151
boundBoxI.H
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::boundBox::operator>>
friend Istream & operator>>(Istream &is, boundBox &bb)
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::boundBox::containsInside
bool containsInside(const point &pt) const
Contains point? (inside only)
Definition: boundBoxI.H:288
Foam::boundBox::minDim
scalar minDim() const
Smallest length/height/width dimension.
Definition: boundBoxI.H:145
Foam::boundBox::boundBox
boundBox()
Construct without any points - an inverted bounding box.
Definition: boundBoxI.H:33
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
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::boundBox::empty
bool empty() const
Bounding box is inverted, contains no points.
Definition: boundBoxI.H:62
Foam::boundBox::directionBit
directionBit
Bits used for (x/y/z) direction encoding.
Definition: boundBox.H:75
Foam::is_contiguous_scalar
A template class to specify if a data type is composed solely of Foam::scalar elements.
Definition: contiguous.H:91
Foam::boundBox::clear
void clear()
Clear bounding box and make it an inverted box.
Definition: boundBoxI.H:184
Foam::boundBox::contains
bool contains(const point &pt) const
Contains point? (inside or on edge)
Definition: boundBoxI.H:271
pointField.H
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< face >
Foam::FixedList
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:104
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
Foam::boundBox::avgDim
scalar avgDim() const
Average length/height/width dimension.
Definition: boundBoxI.H:157
Foam::boundBox::ZDIR
4: z-direction (vector component 2)
Definition: boundBox.H:79
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::boundBox::operator+=
void operator+=(const boundBox &bb)
Extend box to include the second box, as per the add() method.
Definition: boundBoxI.H:301
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::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::boundBox::points
tmp< pointField > points() const
Corner points in an order corresponding to a 'hex' cell.
Definition: boundBox.C:118
Foam::boundBox::operator<<
friend Ostream & operator<<(Ostream &os, const boundBox &bb)
Foam::boundBox::valid
bool valid() const
Bounding box is non-inverted.
Definition: boundBoxI.H:76
boundBoxTemplates.C
Foam::is_contiguous
A template class to specify that a data type can be considered as being contiguous in memory.
Definition: contiguous.H:75
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