boundBox.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-2019 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 Class
25  Foam::boundBox
26 
27 Description
28  A bounding box defined in terms of the points at its extremities.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #ifndef boundBox_H
33 #define boundBox_H
34 
35 #include "pointField.H"
36 #include "faceList.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 
43 // Forward declaration of friend functions and operators
44 
45 class boundBox;
46 template<class T> class tmp;
47 
48 bool operator==(const boundBox&, const boundBox&);
49 bool operator!=(const boundBox&, const boundBox&);
50 
51 Istream& operator>>(Istream&, boundBox&);
52 Ostream& operator<<(Ostream&, const boundBox&);
53 
54 
55 /*---------------------------------------------------------------------------*\
56  Class boundBox Declaration
57 \*---------------------------------------------------------------------------*/
58 
59 class boundBox
60 {
61  // Private Data
62 
63  //- Minimum and maximum describing the bounding box
64  point min_, max_;
65 
66  // Private Member Functions
67 
68  //- Calculate the bounding box from the given points.
69  // Does parallel communication (doReduce = true)
70  void calculate(const UList<point>&, const bool doReduce = true);
71 
72 public:
73 
74  // Static Data Members
75 
76  //- The great value used for greatBox and invertedBox
77  static const scalar great;
78 
79  //- A very large boundBox: min/max == -/+ vGreat
80  static const boundBox greatBox;
81 
82  //- A very large inverted boundBox: min/max == +/- vGreat
83  static const boundBox invertedBox;
84 
85 
86  // Constructors
87 
88  //- Construct null, setting points to zero
89  inline boundBox();
90 
91  //- Construct from components
92  inline boundBox(const point& min, const point& max);
93 
94  //- Construct as the bounding box of the given points
95  // Does parallel communication (doReduce = true)
96  boundBox(const UList<point>&, const bool doReduce = true);
97 
98  //- Construct as the bounding box of the given temporary pointField.
99  // Does parallel communication (doReduce = true)
100  boundBox(const tmp<pointField>&, const bool doReduce = true);
101 
102  //- Construct bounding box as subset of the pointField.
103  // The indices could be from cell/face etc.
104  // Does parallel communication (doReduce = true)
105  boundBox
106  (
107  const UList<point>&,
108  const labelUList& indices,
109  const bool doReduce = true
110  );
111 
112  //- Construct bounding box as subset of the pointField.
113  // The indices could be from edge/triFace etc.
114  // Does parallel communication (doReduce = true)
115  template<unsigned Size>
116  boundBox
117  (
118  const UList<point>&,
119  const FixedList<label, Size>& indices,
120  const bool doReduce = true
121  );
122 
123  //- Construct from Istream
124  inline boundBox(Istream&);
125 
126 
127  // Member Functions
128 
129  // Access
130 
131  //- Minimum describing the bounding box
132  inline const point& min() const;
133 
134  //- Maximum describing the bounding box
135  inline const point& max() const;
136 
137  //- Minimum describing the bounding box, non-const access
138  inline point& min();
139 
140  //- Maximum describing the bounding box, non-const access
141  inline point& max();
142 
143  //- The midpoint of the bounding box
144  inline point midpoint() const;
145 
146  //- The bounding box span (from minimum to maximum)
147  inline vector span() const;
148 
149  //- The magnitude of the bounding box span
150  inline scalar mag() const;
151 
152  //- The volume of the bound box
153  inline scalar volume() const;
154 
155  //- Smallest length/height/width dimension
156  inline scalar minDim() const;
157 
158  //- Largest length/height/width dimension
159  inline scalar maxDim() const;
160 
161  //- Average length/height/width dimension
162  inline scalar avgDim() const;
163 
164  //- Return corner points in an order corresponding to a 'hex' cell
165  tmp<pointField> points() const;
166 
167  //- Return faces with correct point order
168  static faceList faces();
169 
170 
171  // Manipulate
172 
173  //- Inflate box by factor*mag(span) in all dimensions
174  void inflate(const scalar s);
175 
176 
177  // Query
178 
179  //- Overlaps/touches boundingBox?
180  inline bool overlaps(const boundBox&) const;
181 
182  //- Overlaps boundingSphere (centre + sqr(radius))?
183  inline bool overlaps(const point&, const scalar radiusSqr) const;
184 
185  //- Contains point? (inside or on edge)
186  inline bool contains(const point&) const;
187 
188  //- Fully contains other boundingBox?
189  inline bool contains(const boundBox&) const;
190 
191  //- Contains point? (inside only)
192  inline bool containsInside(const point&) const;
193 
194  //- Contains all of the points? (inside or on edge)
195  bool contains(const UList<point>&) const;
196 
197  //- Contains all of the points? (inside or on edge)
198  bool contains
199  (
200  const UList<point>&,
201  const labelUList& indices
202  ) const;
203 
204  //- Contains all of the points? (inside or on edge)
205  template<unsigned Size>
206  bool contains
207  (
208  const UList<point>&,
209  const FixedList<label, Size>& indices
210  ) const;
211 
212 
213  //- Contains any of the points? (inside or on edge)
214  bool containsAny(const UList<point>&) const;
215 
216  //- Contains any of the points? (inside or on edge)
217  bool containsAny
218  (
219  const UList<point>&,
220  const labelUList& indices
221  ) const;
222 
223  //- Contains any of the points? (inside or on edge)
224  template<unsigned Size>
225  bool containsAny
226  (
227  const UList<point>&,
228  const FixedList<label, Size>& indices
229  ) const;
230 
231  //- Return the nearest point on the boundBox to the supplied point.
232  // If point is inside the boundBox then the point is returned
233  // unchanged.
234  point nearest(const point&) const;
235 
236 
237  // Friend Operators
238 
239  inline friend bool operator==(const boundBox&, const boundBox&);
240  inline friend bool operator!=(const boundBox&, const boundBox&);
241 
242 
243  // IOstream operator
244 
245  friend Istream& operator>>(Istream&, boundBox&);
246  friend Ostream& operator<<(Ostream&, const boundBox&);
247 };
248 
249 
250 //- Data associated with boundBox type are contiguous
251 template<>
252 inline bool contiguous<boundBox>() {return contiguous<point>();}
253 
254 
255 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
256 
257 } // End namespace Foam
258 
259 #include "boundBoxI.H"
260 
261 #ifdef NoRepository
262  #include "boundBoxTemplates.C"
263 #endif
264 
265 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
266 
267 #endif
268 
269 // ************************************************************************* //
boundBox()
Construct null, setting points to zero.
Definition: boundBoxI.H:32
bool overlaps(const boundBox &) const
Overlaps/touches boundingBox?
Definition: boundBoxI.H:120
point nearest(const point &) const
Return the nearest point on the boundBox to the supplied point.
Definition: boundBox.C:303
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:54
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
static const scalar great
The great value used for greatBox and invertedBox.
Definition: boundBox.H:76
static faceList faces()
Return faces with correct point order.
Definition: boundBox.C:167
bool contains(const point &) const
Contains point? (inside or on edge)
Definition: boundBoxI.H:170
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
static const boundBox greatBox
A very large boundBox: min/max == -/+ vGreat.
Definition: boundBox.H:79
friend bool operator!=(const boundBox &, const boundBox &)
scalar avgDim() const
Average length/height/width dimension.
Definition: boundBoxI.H:114
tmp< pointField > points() const
Return corner points in an order corresponding to a &#39;hex&#39; cell.
Definition: boundBox.C:149
friend Ostream & operator<<(Ostream &, const boundBox &)
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().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))
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Istream & operator>>(Istream &, directionInfo &)
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:90
point midpoint() const
The midpoint of the bounding box.
Definition: boundBoxI.H:78
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:60
void inflate(const scalar s)
Inflate box by factor*mag(span) in all dimensions.
Definition: boundBox.C:210
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
static const boundBox invertedBox
A very large inverted boundBox: min/max == +/- vGreat.
Definition: boundBox.H:82
Ostream & operator<<(Ostream &, const ensightPart &)
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:84
scalar volume() const
The volume of the bound box.
Definition: boundBoxI.H:96
bool contiguous< boundBox >()
Data associated with boundBox type are contiguous.
Definition: boundBox.H:251
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
scalar maxDim() const
Largest length/height/width dimension.
Definition: boundBoxI.H:108
A class for managing temporary objects.
Definition: PtrList.H:53
bool operator!=(const particle &, const particle &)
Definition: particle.C:1204
friend Istream & operator>>(Istream &, boundBox &)
scalar minDim() const
Smallest length/height/width dimension.
Definition: boundBoxI.H:102
bool containsInside(const point &) const
Contains point? (inside only)
Definition: boundBoxI.H:188
friend bool operator==(const boundBox &, const boundBox &)
Namespace for OpenFOAM.
bool containsAny(const UList< point > &) const
Contains any of the points? (inside or on edge)
Definition: boundBox.C:261