searchableBox.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) 2018-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 "searchableBox.H"
31 #include "SortableList.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(searchableBox, 0);
39  (
40  searchableSurface,
41  searchableBox,
42  dict
43  );
45  (
46  searchableSurface,
47  searchableBox,
48  dict,
49  box
50  );
51 }
52 
53 
54 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55 
56 void Foam::searchableBox::projectOntoCoordPlane
57 (
58  const direction dir,
59  const point& planePt,
60  pointIndexHit& info
61 ) const
62 {
63  // Set point
64  info.rawPoint()[dir] = planePt[dir];
65 
66  // Set face
67  if (planePt[dir] == min()[dir])
68  {
69  info.setIndex(dir*2);
70  }
71  else if (planePt[dir] == max()[dir])
72  {
73  info.setIndex(dir*2+1);
74  }
75  else
76  {
78  << "Point on plane " << planePt
79  << " is not on coordinate " << min()[dir]
80  << " nor " << max()[dir] << abort(FatalError);
81  }
82 }
83 
84 
85 // Returns miss or hit with face (0..5) and region(always 0)
86 Foam::pointIndexHit Foam::searchableBox::findNearest
87 (
88  const point& bbMid,
89  const point& sample,
90  const scalar nearestDistSqr
91 ) const
92 {
93  // Point can be inside or outside. For every component direction can be
94  // left of min, right of max or inbetween.
95  // - outside points: project first one x plane (either min().x()
96  // or max().x()), then onto y plane and finally z. You should be left
97  // with intersection point
98  // - inside point: find nearest side (compare to mid point). Project onto
99  // that.
100 
101  // The face is set to the last projected face.
102 
103 
104  // Outside point projected onto cube. Assume faces 0..5.
105  pointIndexHit info(true, sample, -1);
106  bool outside = false;
107 
108  // (for internal points) per direction what nearest cube side is
109  point near;
110 
111  for (direction dir = 0; dir < vector::nComponents; ++dir)
112  {
113  if (info.rawPoint()[dir] < min()[dir])
114  {
115  projectOntoCoordPlane(dir, min(), info);
116  outside = true;
117  }
118  else if (info.rawPoint()[dir] > max()[dir])
119  {
120  projectOntoCoordPlane(dir, max(), info);
121  outside = true;
122  }
123  else if (info.rawPoint()[dir] > bbMid[dir])
124  {
125  near[dir] = max()[dir];
126  }
127  else
128  {
129  near[dir] = min()[dir];
130  }
131  }
132 
133 
134  // For outside points the info will be correct now. Handle inside points
135  // using the three near distances. Project onto the nearest plane.
136  if (!outside)
137  {
138  vector dist(cmptMag(info.rawPoint() - near));
139 
140  if (dist.x() < dist.y())
141  {
142  if (dist.x() < dist.z())
143  {
144  // Project onto x plane
145  projectOntoCoordPlane(vector::X, near, info);
146  }
147  else
148  {
149  projectOntoCoordPlane(vector::Z, near, info);
150  }
151  }
152  else
153  {
154  if (dist.y() < dist.z())
155  {
156  projectOntoCoordPlane(vector::Y, near, info);
157  }
158  else
159  {
160  projectOntoCoordPlane(vector::Z, near, info);
161  }
162  }
163  }
164 
165 
166  // Check if outside. Optimisation: could do some checks on distance already
167  // on components above
168  if (magSqr(info.rawPoint() - sample) > nearestDistSqr)
169  {
170  info.setMiss();
171  info.setIndex(-1);
172  }
173 
174  return info;
175 }
176 
177 
178 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
179 
180 Foam::searchableBox::searchableBox
181 (
182  const IOobject& io,
183  const treeBoundBox& bb
184 )
185 :
186  searchableSurface(io),
187  treeBoundBox(bb)
188 {
189  if (!treeBoundBox::valid())
190  {
192  << "Illegal bounding box specification : "
193  << static_cast<const treeBoundBox>(*this) << nl
194  << exit(FatalError);
195  }
196 
197  bounds() = static_cast<boundBox>(*this);
198 }
199 
200 
201 Foam::searchableBox::searchableBox
202 (
203  const IOobject& io,
204  const dictionary& dict
205 )
206 :
207  searchableSurface(io),
208  treeBoundBox(dict.get<point>("min"), dict.get<point>("max"))
209 {
210  if (!treeBoundBox::valid())
211  {
213  << "Illegal bounding box specification : "
214  << static_cast<const treeBoundBox>(*this) << nl
215  << exit(FatalError);
216  }
217 
218  bounds() = static_cast<boundBox>(*this);
219 }
220 
221 
222 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
223 
225 {
226  if (regions_.empty())
227  {
228  regions_.setSize(1);
229  regions_[0] = "region0";
230  }
231  return regions_;
232 }
233 
234 
236 {
237  auto tctrs = tmp<pointField>::New(6);
238  auto& ctrs = tctrs.ref();
239 
240  const pointField pts(treeBoundBox::points());
241  const faceList& fcs = treeBoundBox::faces;
242 
243  forAll(fcs, i)
244  {
245  ctrs[i] = fcs[i].centre(pts);
246  }
247 
248  return tctrs;
249 }
250 
251 
253 (
254  pointField& centres,
255  scalarField& radiusSqr
256 ) const
257 {
258  centres.setSize(size());
259  radiusSqr.setSize(size());
260  radiusSqr = 0.0;
261 
262  const pointField pts(treeBoundBox::points());
263  const faceList& fcs = treeBoundBox::faces;
264 
265  forAll(fcs, i)
266  {
267  const face& f = fcs[i];
268 
269  centres[i] = f.centre(pts);
270  for (const label pointi : f)
271  {
272  const point& pt = pts[pointi];
273 
274  radiusSqr[i] = Foam::max
275  (
276  radiusSqr[i],
277  Foam::magSqr(pt-centres[i])
278  );
279  }
280  }
281 
282  // Add a bit to make sure all points are tested inside
283  radiusSqr += Foam::sqr(SMALL);
284 }
285 
286 
288 {
289  return treeBoundBox::points();
290 }
291 
292 
293 Foam::pointIndexHit Foam::searchableBox::findNearest
294 (
295  const point& sample,
296  const scalar nearestDistSqr
297 ) const
298 {
299  return findNearest(centre(), sample, nearestDistSqr);
300 }
301 
302 
304 (
305  const point& sample,
306  const scalar nearestDistSqr
307 ) const
308 {
309  const point bbMid(centre());
310 
311  // Outside point projected onto cube. Assume faces 0..5.
312  pointIndexHit info(true, sample, -1);
313  bool outside = false;
314 
315  // (for internal points) per direction what nearest cube side is
316  point near;
317 
318  for (direction dir = 0; dir < vector::nComponents; ++dir)
319  {
320  if (info.rawPoint()[dir] < min()[dir])
321  {
322  projectOntoCoordPlane(dir, min(), info);
323  outside = true;
324  }
325  else if (info.rawPoint()[dir] > max()[dir])
326  {
327  projectOntoCoordPlane(dir, max(), info);
328  outside = true;
329  }
330  else if (info.rawPoint()[dir] > bbMid[dir])
331  {
332  near[dir] = max()[dir];
333  }
334  else
335  {
336  near[dir] = min()[dir];
337  }
338  }
339 
340 
341  // For outside points the info will be correct now. Handle inside points
342  // using the three near distances. Project onto the nearest two planes.
343  if (!outside)
344  {
345  // Get the per-component distance to nearest wall
346  vector dist(cmptMag(info.rawPoint() - near));
347 
348  SortableList<scalar> sortedDist(3);
349  sortedDist[0] = dist[0];
350  sortedDist[1] = dist[1];
351  sortedDist[2] = dist[2];
352  sortedDist.sort();
353 
354  // Project onto nearest
355  projectOntoCoordPlane(sortedDist.indices()[0], near, info);
356  // Project onto second nearest
357  projectOntoCoordPlane(sortedDist.indices()[1], near, info);
358  }
359 
360 
361  // Check if outside. Optimisation: could do some checks on distance already
362  // on components above
363  if (magSqr(info.rawPoint() - sample) > nearestDistSqr)
364  {
365  info.setMiss();
366  info.setIndex(-1);
367  }
368 
369  return info;
370 }
371 
372 
373 Foam::pointIndexHit Foam::searchableBox::findNearest
374 (
375  const linePointRef& ln,
376  treeBoundBox& tightest,
377  point& linePoint
378 ) const
379 {
381  return pointIndexHit();
382 }
383 
384 
386 (
387  const point& start,
388  const point& end
389 ) const
390 {
391  pointIndexHit info(false, start, -1);
392 
393  bool foundInter;
394 
395  if (posBits(start) == 0)
396  {
397  if (posBits(end) == 0)
398  {
399  // Both start and end inside.
400  foundInter = false;
401  }
402  else
403  {
404  // end is outside. Clip to bounding box.
405  foundInter = intersects(end, start, info.rawPoint());
406  }
407  }
408  else
409  {
410  // start is outside. Clip to bounding box.
411  foundInter = intersects(start, end, info.rawPoint());
412  }
413 
414 
415  // Classify point
416  if (foundInter)
417  {
418  info.setHit();
419 
420  for (direction dir = 0; dir < vector::nComponents; ++dir)
421  {
422  if (info.rawPoint()[dir] == min()[dir])
423  {
424  info.setIndex(2*dir);
425  break;
426  }
427  else if (info.rawPoint()[dir] == max()[dir])
428  {
429  info.setIndex(2*dir+1);
430  break;
431  }
432  }
433 
434  if (info.index() == -1)
435  {
437  << "point " << info.rawPoint()
438  << " on segment " << start << end
439  << " should be on face of " << *this
440  << " but it isn't." << abort(FatalError);
441  }
442  }
443 
444  return info;
445 }
446 
447 
449 (
450  const point& start,
451  const point& end
452 ) const
453 {
454  return findLine(start, end);
455 }
456 
457 
458 void Foam::searchableBox::findNearest
459 (
460  const pointField& samples,
461  const scalarField& nearestDistSqr,
462  List<pointIndexHit>& info
463 ) const
464 {
465  info.setSize(samples.size());
466 
467  const point bbMid(centre());
468 
469  forAll(samples, i)
470  {
471  info[i] = findNearest(bbMid, samples[i], nearestDistSqr[i]);
472  }
473 }
474 
475 
477 (
478  const pointField& start,
479  const pointField& end,
480  List<pointIndexHit>& info
481 ) const
482 {
483  info.setSize(start.size());
484 
485  forAll(start, i)
486  {
487  info[i] = findLine(start[i], end[i]);
488  }
489 }
490 
491 
493 (
494  const pointField& start,
495  const pointField& end,
496  List<pointIndexHit>& info
497 ) const
498 {
499  info.setSize(start.size());
500 
501  forAll(start, i)
502  {
503  info[i] = findLineAny(start[i], end[i]);
504  }
505 }
506 
507 
509 (
510  const pointField& start,
511  const pointField& end,
513 ) const
514 {
515  info.setSize(start.size());
516 
517  // Work array
519 
520  // Tolerances:
521  // To find all intersections we add a small vector to the last intersection
522  // This is chosen such that
523  // - it is significant (SMALL is smallest representative relative tolerance;
524  // we need something bigger since we're doing calculations)
525  // - if the start-end vector is zero we still progress
526  const vectorField dirVec(end-start);
527  const scalarField magSqrDirVec(magSqr(dirVec));
528  const vectorField smallVec
529  (
530  ROOTSMALL*dirVec
531  + vector(ROOTVSMALL,ROOTVSMALL,ROOTVSMALL)
532  );
533 
534  forAll(start, pointi)
535  {
536  // See if any intersection between pt and end
537  pointIndexHit inter = findLine(start[pointi], end[pointi]);
538 
539  if (inter.hit())
540  {
541  hits.clear();
542  hits.append(inter);
543 
544  point pt = inter.hitPoint() + smallVec[pointi];
545 
546  while (((pt-start[pointi])&dirVec[pointi]) <= magSqrDirVec[pointi])
547  {
548  // See if any intersection between pt and end
549  pointIndexHit inter = findLine(pt, end[pointi]);
550 
551  // Check for not hit or hit same face as before (can happen
552  // if vector along surface of face)
553  if
554  (
555  !inter.hit()
556  || (inter.index() == hits.last().index())
557  )
558  {
559  break;
560  }
561  hits.append(inter);
562 
563  pt = inter.hitPoint() + smallVec[pointi];
564  }
565 
566  info[pointi].transfer(hits);
567  }
568  else
569  {
570  info[pointi].clear();
571  }
572  }
573 }
574 
575 
577 (
578  const List<pointIndexHit>& info,
579  labelList& region
580 ) const
581 {
582  region.setSize(info.size());
583  region = 0;
584 }
585 
586 
588 (
589  const List<pointIndexHit>& info,
590  vectorField& normal
591 ) const
592 {
593  normal.setSize(info.size());
594  normal = Zero;
595 
596  forAll(info, i)
597  {
598  if (info[i].hit())
599  {
600  normal[i] = treeBoundBox::faceNormals[info[i].index()];
601  }
602  else
603  {
604  // Set to what?
605  }
606  }
607 }
608 
609 
611 (
612  const pointField& points,
613  List<volumeType>& volType
614 ) const
615 {
616  volType.setSize(points.size());
617 
618  forAll(points, pointi)
619  {
620  const point& pt = points[pointi];
621 
623 
624  for (direction dir=0; dir < vector::nComponents; ++dir)
625  {
626  if (pt[dir] < min()[dir] || pt[dir] > max()[dir])
627  {
628  vt = volumeType::OUTSIDE;
629  break;
630  }
631  }
632 
633  volType[pointi] = vt;
634  }
635 }
636 
637 
638 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::PointIndexHit::index
label index() const
Return index.
Definition: PointIndexHit.H:113
Foam::Vector< scalar >::Z
Definition: Vector.H:80
Foam::searchableBox::findLine
pointIndexHit findLine(const point &start, const point &end) const
Find nearest intersection of line between start and end.
Definition: searchableBox.C:386
Foam::Vector< scalar >::Y
Definition: Vector.H:80
Foam::PointIndexHit::setIndex
void setIndex(const label index)
Definition: PointIndexHit.H:170
Foam::SortableList::sort
void sort()
Forward (stable) sort the list (if changed after construction).
Definition: SortableList.C:138
Foam::treeBoundBox::points
tmp< pointField > points() const
Vertex coordinates. In octant coding.
Definition: treeBoundBox.C:92
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:57
searchableBox.H
Foam::searchableBox::findLineAny
pointIndexHit findLineAny(const point &start, const point &end) const
Find any intersection of line between start and end.
Definition: searchableBox.C:449
Foam::treeBoundBox
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:87
Foam::PointIndexHit::hit
bool hit() const
Is there a hit.
Definition: PointIndexHit.H:107
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
Foam::PointIndexHit::setMiss
void setMiss()
Definition: PointIndexHit.H:160
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::IOobject::info
InfoProxy< IOobject > info() const
Return info proxy.
Definition: IOobject.H:519
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
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::cmptMag
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:400
SortableList.H
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:419
Foam::searchableBox::getRegion
virtual void getRegion(const List< pointIndexHit > &, labelList &region) const
From a set of points and indices get the region.
Definition: searchableBox.C:577
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:55
Foam::PointIndexHit::rawPoint
const Point & rawPoint() const
Return point with no checking.
Definition: PointIndexHit.H:145
Foam::searchableBox::findNearestOnEdge
pointIndexHit findNearestOnEdge(const point &sample, const scalar nearestDistSqr) const
Calculate nearest point on edge.
Definition: searchableBox.C:304
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::volumeType
An enumeration wrapper for classification of a location as being inside/outside of a volume.
Definition: volumeType.H:60
Foam::Field< vector >
Foam::PointIndexHit::hitPoint
const Point & hitPoint() const
Return hit point.
Definition: PointIndexHit.H:119
Foam::treeBoundBox::faces
static const faceList faces
Face to point addressing.
Definition: treeBoundBox.H:150
Foam::Vector< scalar >::X
Definition: Vector.H:80
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:472
Foam::searchableSurface
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
Definition: searchableSurface.H:69
Foam::searchableBox::coordinates
virtual tmp< pointField > coordinates() const
Get representative set of element coordinates.
Definition: searchableBox.C:235
Foam::DynamicList::clear
void clear()
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:348
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:436
Foam::addNamedToRunTimeSelectionTable
addNamedToRunTimeSelectionTable(topoSetCellSource, badQualityToCell, word, badQuality)
samples
scalarField samples(nIntervals, Zero)
Foam::PointIndexHit::setHit
void setHit()
Definition: PointIndexHit.H:155
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::searchableBox::points
virtual tmp< pointField > points() const
Get the points that define the surface.
Definition: searchableBox.C:287
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::SortableList
A list that is sorted upon construction or when explicitly requested with the sort() method.
Definition: List.H:66
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:115
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::searchableBox::getVolumeType
virtual void getVolumeType(const pointField &points, List< volumeType > &volType) const
Determine type (inside/outside) for points.
Definition: searchableBox.C:611
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::pointIndexHit
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:45
Foam::nl
constexpr char nl
Definition: Ostream.H:372
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List< word >
Foam::start
label ListType::const_reference const label start
Definition: ListOps.H:408
Foam::searchableBox::regions
virtual const wordList & regions() const
Names of regions.
Definition: searchableBox.C:224
Foam::searchableBox::findLineAll
virtual void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit >> &) const
Get all intersections in order from start to end.
Definition: searchableBox.C:509
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:115
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::direction
uint8_t direction
Definition: direction.H:47
Foam::volumeType::INSIDE
A location inside the volume.
Definition: volumeType.H:68
Foam::line
A line primitive.
Definition: line.H:59
Foam::SortableList::indices
const labelList & indices() const
Return the list of sorted indices. Updated every sort.
Definition: SortableList.H:113
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::ln
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: MSwindows.C:915
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:74
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::boundBox::valid
bool valid() const
Bounding box is non-inverted.
Definition: boundBoxI.H:76
Foam::searchableBox::getNormal
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
Definition: searchableBox.C:588
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::searchableBox::boundingSpheres
virtual void boundingSpheres(pointField &centres, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared), one per element.
Definition: searchableBox.C:253
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::nComponents
static constexpr direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:101
Foam::volumeType::OUTSIDE
A location outside the volume.
Definition: volumeType.H:69