Loading...
Searching...
No Matches
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-2025 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
27Class
28 Foam::boundBox
29
30Description
31 A bounding box defined in terms of min/max extrema points.
32
33Note
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
38SeeAlso
39 Foam::treeBoundBox
40
41SourceFiles
42 boundBoxI.H
43 boundBox.cxx
44 boundBox.txx
45
46\*---------------------------------------------------------------------------*/
47
48#ifndef Foam_boundBox_H
49#define Foam_boundBox_H
50
51#include "pointField.H"
52#include "faceList.H"
53#include "Pair.H"
54#include "triangleFwd.H"
55
56// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57
58namespace Foam
59{
60
61// Forward Declarations
63class plane;
64class Random;
65template<class T> class MinMax;
66
69
71/*---------------------------------------------------------------------------*\
72 Class boundBox Declaration
73\*---------------------------------------------------------------------------*/
74
75class boundBox
76{
77 // Private Data
78
79 //- Minimum and maximum points describing the bounding box
80 point min_, max_;
81
82
83protected:
84
85 // Protected Member Functions
86
87 //- Test for overlap of box and box (inclusive check)
88 inline static bool box_box_overlaps
89 (
90 const point& minA, // boxA (min)
91 const point& maxA, // boxA (max)
92 const point& minB, // boxB (min)
93 const point& maxB // boxB (max)
94 );
95
96 //- Test for overlap of box and boundingSphere (centre + sqr(radius))
97 // Note: ordering of corners is irrelevant
98 inline static bool box_sphere_overlaps
99 (
100 const point& corner0, // box corner
101 const point& corner1, // box corner
102 const point& centre, // sphere centre
103 const scalar radiusSqr // sqr(radius)
104 );
105
106
107public:
108
109 // Static Data Members
110
111 //- Bits used for (x/y/z) directional encoding.
113 {
114 XDIR = 1,
115 YDIR = 2,
116 ZDIR = 4
117 };
119 //- A large boundBox: min/max == -/+ ROOTVGREAT
120 static const boundBox greatBox;
121
122 //- A large inverted boundBox: min/max == +/- ROOTVGREAT
123 static const boundBox invertedBox;
124
125 //- The unit normal per face
127
128
129 // Static Methods
130
131 //- The null boundBox is the same as an inverted box
132 static const boundBox& null() noexcept
133 {
134 return invertedBox;
135 }
137 //- The boundBox faces as a hexCell, using hexCorner points.
138 //- Same as hexCell::modelFaces()
139 static const Foam::faceList& hexFaces();
140
141 //- Number of points for boundBox and HEX
142 static constexpr label nPoints() noexcept { return 8; }
143
144 //- Number of edges for boundBox and HEX
145 static constexpr label nEdges() noexcept { return 12; }
146
147 //- Number of faces for boundBox and HEX
148 static constexpr label nFaces() noexcept { return 6; }
149
150
151 // Standard (Generated) Methods
152
153 //- Default construct: an inverted bounding box
154 inline boundBox();
155
156 //- Copy construct
157 boundBox(const boundBox&) = default;
159 //- Copy assignment
160 boundBox& operator=(const boundBox&) = default;
161
162
163 // Constructors
164
165 //- Copy construct with specified global reduction
166 boundBox(const boundBox& bb, bool doReduce);
167
168 //- Construct a bounding box containing a single initial point
169 inline explicit boundBox(const point& p);
170
171 //- Construct a 0/1 unit bounding box
172 inline explicit boundBox(Foam::zero_one);
173
174 //- Construct from bound box min/max points
175 inline boundBox(const point& min, const point& max);
176
179
180 //- Construct from bound box min/max points
181 inline explicit boundBox(const Pair<point>& bb);
182
183 //- Construct as the bounding box of the given points
184 // Does parallel communication (doReduce = true)
185 explicit boundBox(const UList<point>& points, bool doReduce = true);
187 //- Construct as the bounding box of the given temporary pointField.
188 // Does parallel communication (doReduce = true)
189 explicit boundBox(const tmp<pointField>& tpoints, bool doReduce = true);
190
191 //- Construct bounding box as an indirect subset of the points.
192 // The indices could be from cell/face etc.
193 // Does parallel communication (doReduce = true)
195 (
196 const UList<point>& points,
197 const labelUList& indices,
198 bool doReduce = true
199 );
200
201 //- Construct bounding box as an indirect subset of the points.
202 // The indices could be from edge/triFace etc.
203 // Does parallel communication (doReduce = true)
204 template<unsigned N>
206 (
207 const UList<point>& points,
208 const FixedList<label, N>& indices,
209 bool doReduce = true
210 );
211
212 //- Construct from Istream
213 inline explicit boundBox(Istream& is);
214
215
216 // Member Functions
217
218 // Access
219
220 //- Bounding box is inverted, contains no points.
221 inline bool empty() const;
222
223 //- Bounding box is non-inverted.
224 inline bool good() const;
225
226 //- Bounding box is non-inverted - same as good().
227 bool valid() const { return good(); }
228
229 //- Minimum describing the bounding box
230 inline const point& min() const noexcept;
232 //- Maximum describing the bounding box
233 inline const point& max() const noexcept;
234
235 //- Minimum describing the bounding box, non-const access
236 inline point& min() noexcept;
237
238 //- Maximum describing the bounding box, non-const access
239 inline point& max() noexcept;
240
241 //- The centre (midpoint) of the bounding box
242 inline point centre() const;
243
244 //- The bounding box span (from minimum to maximum)
245 inline vector span() const;
246
247 //- The magnitude/length of the bounding box diagonal
248 inline scalar mag() const;
249
250 //- The magnitude/length squared of bounding box diagonal
251 inline scalar magSqr() const;
252
253 //- The volume of the bound box
254 inline scalar volume() const;
255
256 //- Smallest length/height/width dimension
257 inline scalar minDim() const;
258
259 //- Largest length/height/width dimension
260 inline scalar maxDim() const;
261
262 //- Average length/height/width dimension
263 inline scalar avgDim() const;
264
265 //- Direction (X/Y/Z) of the smallest span (for empty box: 0)
266 inline direction minDir() const;
267
268 //- Direction (X/Y/Z) of the largest span (for empty box: 0)
269 inline direction maxDir() const;
270
271 //- Count the number of positive, non-zero dimensions.
272 // \return
273 // - -1 : if any dimensions are negative
274 // - 0 : 0D (point)
275 // - 1 : 1D (line aligned with an axis)
276 // - 2 : 2D (plane aligned with an axis)
277 // - 3 : 3D (box)
278 inline int nDim() const;
279
280 //- Return corner point [0..7] corresponding to a 'hex' cell
281 template<direction CornerNumber>
282 inline point hexCorner() const;
284 //- Corner points in an order corresponding to a 'hex' cell
285 tmp<pointField> hexCorners() const;
286
287 //- Corner points in an order corresponding to a 'hex' cell
288 tmp<pointField> points() const { return hexCorners(); }
289
290 //- Face midpoints
292
293 //- Face centre of given face index
294 point faceCentre(const direction facei) const;
295
296
297 // Manipulate
298
299 //- Reset to an inverted box
300 inline void reset();
301
302 //- Reset to a 0/1 unit bounding box
303 inline void reset(Foam::zero_one);
304
305 //- Reset min/max to be identical to the specified point
306 inline void reset(const point& pt);
307
308 //- Reset min/max to specified values
309 inline void reset(const point& min, const point& max);
310
311 //- Same as reset() - reset to an inverted box
312 void clear() { reset(); }
313
314 //- Extend to include the second box.
315 inline void add(const boundBox& bb);
316
317 //- Extend to include the point.
318 inline void add(const point& pt);
319
320 //- Extend to include two additional points.
321 inline void add(const point& pt0, const point& pt1);
322
323 //- Extend to include two additional points.
324 inline void add(const Pair<point>& points);
325
326 //- Extend to include the points.
327 inline void add(const UList<point>& points);
328
329 //- Extend to include the points from the temporary point field.
330 inline void add(const tmp<pointField>& tpoints);
331
332 //- Extend to include the points.
333 template<unsigned N>
334 void add(const FixedList<point, N>& points);
335
336 //- Extend to include a (subsetted) point field.
337 // The indices could be from edge/triFace etc.
338 template<unsigned N>
339 void add
340 (
341 const UList<point>& points,
342 const FixedList<label, N>& indices
343 );
344
345 //- Extend to include a (subsetted) point field.
346 //
347 // \tparam IntContainer A container with an iterator that
348 // dereferences to an label
349 template<class IntContainer>
350 void add
351 (
352 const UList<point>& points,
353 const IntContainer& indices
354 );
355
356 //- Expand box by adjusting min/max by specified amount
357 //- in each dimension
358 inline void grow(const scalar delta);
359
360 //- Expand box by adjusting min/max by specified amounts
361 inline void grow(const vector& delta);
362
363 //- Expand box by factor*mag(span) in all dimensions
364 inline void inflate(const scalar factor);
365
366 //- Expand box slightly by expanding all dimensions with
367 //- factor*span*(random 0-1) and guarantees
368 //- factor*mag(span) minimum width in any direction
369 void inflate(Random& rndGen, const scalar factor);
370
371 //- As per two parameter version but with additional
372 //- grow() by given amount in each dimension
373 void inflate(Random& r, const scalar factor, const scalar delta);
374
375
376 // Parallel
377
378 //- Inplace parallel reduction of min/max values,
379 //- using UPstream::worldComm
380 void reduce();
382 //- Inplace parallel reduction of min/max values,
383 //- using the specified communicator
384 void reduce(int communicator);
385
386 //- Perform a reduction on a copy and return the result,
387 //- using UPstream::worldComm
388 static boundBox returnReduce(const boundBox& bb);
389
390 //- Perform a reduction on a copy and return the result,
391 //- using the specified communicator
392 static boundBox returnReduce(const boundBox& bb, int communicator);
393
394
395 // Query
396
397 //- Does plane intersect this bounding box.
398 // There is an intersection if the plane segments the corner points
399 // \note the results are unreliable when plane coincides almost
400 // exactly with a box face
401 bool intersects(const plane& pln) const;
402
403 //- Does triangle intersect this bounding box
404 //- or is contained within this bounding box.
405 // \note results may be unreliable when it coincides almost
406 // exactly with a box face
407 bool intersects(const triPointRef& tri) const;
408
409 //- Overlaps/touches boundingBox?
410 inline bool overlaps(const boundBox& bb) const;
411
412 //- Overlaps boundingSphere (centre + sqr(radius))?
413 inline bool overlaps
414 (
415 const point& centre,
416 const scalar radiusSqr
417 ) const;
418
419 //- Contains point? (inside or on edge)
420 inline bool contains(const point& pt) const;
421
422 //- Fully contains other boundingBox?
423 inline bool contains(const boundBox& bb) const;
424
425 //- Contains point? (inside only)
426 inline bool containsInside(const point& pt) const;
427
428 //- Contains all points? (inside or on edge)
429 bool contains(const UList<point>& points) const;
430
431 //- Contains all of the (subsetted) points? (inside or on edge)
432 template<unsigned N>
433 bool contains
434 (
435 const UList<point>& points,
436 const FixedList<label, N>& indices
437 ) const;
438
439
440 //- Contains all of the (subsetted) points? (inside or on edge)
441 //
442 // \tparam IntContainer A container with an iterator that
443 // dereferences to an label
444 template<class IntContainer>
445 bool contains
446 (
447 const UList<point>& points,
448 const IntContainer& indices
449 ) const;
450
451
452 //- Contains any of the points? (inside or on edge)
453 bool containsAny(const UList<point>& points) const;
454
455 //- Contains any of the (subsetted) points? (inside or on edge)
456 template<unsigned N>
457 bool containsAny
458 (
459 const UList<point>& points,
460 const FixedList<label, N>& indices
461 ) const;
462
463 //- Contains any of the (subsetted) points? (inside or on edge)
464 //
465 // \tparam IntContainer A container with an iterator that
466 // dereferences to an label
467 template<class IntContainer>
468 bool containsAny
469 (
470 const UList<point>& points,
471 const IntContainer& indices
472 ) const;
473
474
475 //- Return the nearest point on the boundBox to the supplied point.
476 // If point is inside the boundBox then the point is returned
477 // unchanged.
478 point nearest(const point& p) const;
479
480
481 // Member Operators
482
483 //- Restrict min/max to union with other box.
484 void operator&=(const boundBox& bb);
485
486 //- Extend box to include the second box, as per the add() method.
487 // Can be used in a reduction operation.
488 void operator+=(const boundBox& bb) { add(bb); }
489
490
491 // IOstream Operators
492
493 friend Istream& operator>>(Istream& is, boundBox& bb);
494 friend Ostream& operator<<(Ostream& os, const boundBox& bb);
495
496
497 // Housekeeping
498
499 //- Deprecated(2022-10) - use 'operator&=' to avoid confusion with
500 //- other intersects() methods.
501 // \deprecated(2022-10) - use 'operator&='
502 FOAM_DEPRECATED_FOR(2022-10, "boundBox::operator&=(const boundBox&)")
503 bool intersect(const boundBox& bb) { *this &= bb; return good(); }
505 //- Identical to centre()
506 point midpoint() const { return centre(); }
507};
508
509
510// * * * * * * * * * * * * * * * * * Traits * * * * * * * * * * * * * * * * //
511
512//- Contiguous data for boundBox
513template<> struct is_contiguous<boundBox> : is_contiguous<point> {};
514
515//- Contiguous scalar data for boundBox
516template<> struct is_contiguous_scalar<boundBox>
517:
519{};
520
521
522// * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
523
524inline bool operator==(const boundBox& a, const boundBox& b)
526 return (a.min() == b.min() && a.max() == b.max());
527}
528
529
530inline bool operator!=(const boundBox& a, const boundBox& b)
532 return !(a == b);
533}
534
535
536inline bool operator<(const boundBox& a, const boundBox& b)
538 return
539 (
540 a.min() < b.min()
541 || (!(b.min() < a.min()) && a.max() < b.max())
542 );
543}
544
545
546inline bool operator<=(const boundBox& a, const boundBox& b)
547{
548 return !(b < a);
550
551
552inline bool operator>(const boundBox& a, const boundBox& b)
553{
554 return (b < a);
555}
556
557
558inline bool operator>=(const boundBox& a, const boundBox& b)
559{
560 return !(a < b);
561}
562
563
564// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
565
566} // End namespace Foam
567
568#include "boundBoxI.H"
569
570#ifdef NoRepository
571 #include "boundBox.txx"
572#endif
573
574// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
575
576#endif
577
578// ************************************************************************* //
scalar delta
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
A min/max value pair with additional methods. In addition to conveniently storing values,...
Definition MinMax.H:106
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
An ordered pair of two objects of type <T> with first() and second() elements.
Definition Pair.H:66
Random number generator.
Definition Random.H:56
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:89
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
boundBox(const boundBox &bb, bool doReduce)
Copy construct with specified global reduction.
bool containsInside(const point &pt) const
Contains point? (inside only).
Definition boundBoxI.H:472
int nDim() const
Count the number of positive, non-zero dimensions.
Definition boundBoxI.H:274
static boundBox returnReduce(const boundBox &bb)
Perform a reduction on a copy and return the result, using UPstream::worldComm.
point faceCentre(const direction facei) const
Face centre of given face index.
static const Foam::faceList & hexFaces()
The boundBox faces as a hexCell, using hexCorner points. Same as hexCell::modelFaces().
const point & max() const noexcept
Maximum describing the bounding box.
Definition boundBoxI.H:168
friend Ostream & operator<<(Ostream &os, const boundBox &bb)
void grow(const scalar delta)
Expand box by adjusting min/max by specified amount in each dimension.
Definition boundBoxI.H:367
bool contains(const UList< point > &points, const FixedList< label, N > &indices) const
Contains all of the (subsetted) points? (inside or on edge).
bool valid() const
Bounding box is non-inverted - same as good().
Definition boundBox.H:283
void operator+=(const boundBox &bb)
Extend box to include the second box, as per the add() method.
Definition boundBox.H:671
bool overlaps(const boundBox &bb) const
Overlaps/touches boundingBox?
Definition boundBoxI.H:439
void reduce()
Inplace parallel reduction of min/max values, using UPstream::worldComm.
directionBit
Bits used for (x/y/z) directional encoding.
Definition boundBox.H:117
@ XDIR
1: x-direction. Same as (1 << vector::X)
Definition boundBox.H:118
@ ZDIR
4: z-direction. Same as (1 << vector::Z)
Definition boundBox.H:120
@ YDIR
2: y-direction. Same as (1 << vector::Y)
Definition boundBox.H:119
static const boundBox invertedBox
A large inverted boundBox: min/max == +/- ROOTVGREAT.
Definition boundBox.H:131
tmp< pointField > points() const
Corner points in an order corresponding to a 'hex' cell.
Definition boundBox.H:381
static bool box_sphere_overlaps(const point &corner0, const point &corner1, const point &centre, const scalar radiusSqr)
Test for overlap of box and boundingSphere (centre + sqr(radius)).
Definition boundBoxI.H:403
static boundBox returnReduce(const boundBox &bb, int communicator)
Perform a reduction on a copy and return the result, using the specified communicator.
direction minDir() const
Direction (X/Y/Z) of the smallest span (for empty box: 0).
Definition boundBoxI.H:234
bool containsAny(const UList< point > &points, const FixedList< label, N > &indices) const
Contains any of the (subsetted) points? (inside or on edge).
const point & min() const noexcept
Minimum describing the bounding box.
Definition boundBoxI.H:162
boundBox & operator=(const boundBox &)=default
Copy assignment.
bool contains(const point &pt) const
Contains point? (inside or on edge).
Definition boundBoxI.H:455
scalar volume() const
The volume of the bound box.
Definition boundBoxI.H:210
static const boundBox & null() noexcept
The null boundBox is the same as an inverted box.
Definition boundBox.H:144
void add(const UList< point > &points, const FixedList< label, N > &indices)
Extend to include a (subsetted) point field.
bool empty() const
Bounding box is inverted, contains no points.
Definition boundBoxI.H:144
void add(const UList< point > &points, const IntContainer &indices)
Extend to include a (subsetted) point field.
boundBox(const UList< point > &points, const labelUList &indices, bool doReduce=true)
Construct bounding box as an indirect subset of the points.
void add(const boundBox &bb)
Extend to include the second box.
Definition boundBoxI.H:323
tmp< pointField > hexCorners() const
Corner points in an order corresponding to a 'hex' cell.
tmp< pointField > faceCentres() const
Face midpoints.
direction maxDir() const
Direction (X/Y/Z) of the largest span (for empty box: 0).
Definition boundBoxI.H:254
void add(const FixedList< point, N > &points)
Extend to include the points.
void reduce(int communicator)
Inplace parallel reduction of min/max values, using the specified communicator.
point midpoint() const
Identical to centre().
Definition boundBox.H:694
friend Istream & operator>>(Istream &is, boundBox &bb)
bool containsAny(const UList< point > &points, const IntContainer &indices) const
Contains any of the (subsetted) points? (inside or on edge).
bool contains(const UList< point > &points) const
Contains all points? (inside or on edge).
void inflate(Random &rndGen, const scalar factor)
Expand box slightly by expanding all dimensions with factor*span*(random 0-1) and guarantees factor*m...
void inflate(Random &r, const scalar factor, const scalar delta)
As per two parameter version but with additional grow() by given amount in each dimension.
scalar minDim() const
Smallest length/height/width dimension.
Definition boundBoxI.H:216
void operator&=(const boundBox &bb)
Restrict min/max to union with other box.
scalar mag() const
The magnitude/length of the bounding box diagonal.
Definition boundBoxI.H:198
boundBox(const tmp< pointField > &tpoints, bool doReduce=true)
Construct as the bounding box of the given temporary pointField.
static const boundBox greatBox
A large boundBox: min/max == -/+ ROOTVGREAT.
Definition boundBox.H:126
point nearest(const point &p) const
Return the nearest point on the boundBox to the supplied point.
scalar avgDim() const
Average length/height/width dimension.
Definition boundBoxI.H:228
boundBox()
Default construct: an inverted bounding box.
Definition boundBoxI.H:101
boundBox(const UList< point > &points, bool doReduce=true)
Construct as the bounding box of the given points.
bool good() const
Bounding box is non-inverted.
Definition boundBoxI.H:156
void clear()
Same as reset() - reset to an inverted box.
Definition boundBox.H:419
static constexpr label nFaces() noexcept
Number of faces for boundBox and HEX.
Definition boundBox.H:168
boundBox(const boundBox &)=default
Copy construct.
bool contains(const UList< point > &points, const IntContainer &indices) const
Contains all of the (subsetted) points? (inside or on edge).
bool intersects(const plane &pln) const
Does plane intersect this bounding box.
void reset()
Reset to an inverted box.
Definition boundBoxI.H:295
boundBox(const UList< point > &points, const FixedList< label, N > &indices, bool doReduce=true)
Construct bounding box as an indirect subset of the points.
scalar magSqr() const
The magnitude/length squared of bounding box diagonal.
Definition boundBoxI.H:204
static constexpr label nEdges() noexcept
Number of edges for boundBox and HEX.
Definition boundBox.H:163
void inflate(const scalar factor)
Expand box by factor*mag(span) in all dimensions.
Definition boundBoxI.H:381
bool intersects(const triPointRef &tri) const
Does triangle intersect this bounding box or is contained within this bounding box.
static constexpr label nPoints() noexcept
Number of points for boundBox and HEX.
Definition boundBox.H:158
bool intersect(const boundBox &bb)
Deprecated(2022-10) - use 'operator&=' to avoid confusion with other intersects() methods.
Definition boundBox.H:689
point hexCorner() const
Return corner point [0..7] corresponding to a 'hex' cell.
vector span() const
The bounding box span (from minimum to maximum).
Definition boundBoxI.H:192
static bool box_box_overlaps(const point &minA, const point &maxA, const point &minB, const point &maxB)
Test for overlap of box and box (inclusive check).
Definition boundBoxI.H:388
point centre() const
The centre (midpoint) of the bounding box.
Definition boundBoxI.H:186
bool containsAny(const UList< point > &points) const
Contains any of the points? (inside or on edge).
scalar maxDim() const
Largest length/height/width dimension.
Definition boundBoxI.H:222
static const FixedList< vector, 6 > faceNormals
The unit normal per face.
Definition boundBox.H:136
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition plane.H:91
A class for managing temporary objects.
Definition tmp.H:75
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
Definition pTraits.H:51
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
Namespace for OpenFOAM.
bool operator!=(const eddy &a, const eddy &b)
Definition eddy.H:297
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
bool operator<=(const IOstreamOption::versionNumber &a, const IOstreamOption::versionNumber &b) noexcept
Version A same or older than B.
bool operator>=(const IOstreamOption::versionNumber &a, const IOstreamOption::versionNumber &b) noexcept
Version A same or newer than B.
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces).
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Istream & operator>>(Istream &, directionInfo &)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
uint8_t direction
Definition direction.H:49
vector point
Point is a vector.
Definition point.H:37
const direction noexcept
Definition scalarImpl.H:265
bool operator<(const IOstreamOption::versionNumber &a, const IOstreamOption::versionNumber &b) noexcept
Version A older than B.
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition triangleFwd.H:46
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
UList< label > labelUList
A UList of labels.
Definition UList.H:75
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
bool operator>(const IOstreamOption::versionNumber &a, const IOstreamOption::versionNumber &b) noexcept
Version A newer than B.
volScalarField & b
#define FOAM_DEPRECATED_FOR(since, replacement)
Definition stdFoam.H:43
A template class to specify if a data type is composed solely of Foam::scalar elements.
Definition contiguous.H:87
A template class to specify that a data type can be considered as being contiguous in memory.
Definition contiguous.H:70
Random rndGen