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-2020 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] << nl
81  << abort(FatalError);
82  }
83 }
84 
85 
86 // Returns miss or hit with face (0..5) and region(always 0)
87 Foam::pointIndexHit Foam::searchableBox::findNearest
88 (
89  const point& bbMid,
90  const point& sample,
91  const scalar nearestDistSqr
92 ) const
93 {
94  // Point can be inside or outside. For every component direction can be
95  // left of min, right of max or inbetween.
96  // - outside points: project first one x plane (either min().x()
97  // or max().x()), then onto y plane and finally z. You should be left
98  // with intersection point
99  // - inside point: find nearest side (compare to mid point). Project onto
100  // that.
101 
102  // The face is set to the last projected face.
103 
104 
105  // Outside point projected onto cube. Assume faces 0..5.
106  pointIndexHit info(true, sample, -1);
107  bool outside = false;
108 
109  // (for internal points) per direction what nearest cube side is
110  point near;
111 
112  for (direction dir = 0; dir < vector::nComponents; ++dir)
113  {
114  if (info.rawPoint()[dir] < min()[dir])
115  {
116  projectOntoCoordPlane(dir, min(), info);
117  outside = true;
118  }
119  else if (info.rawPoint()[dir] > max()[dir])
120  {
121  projectOntoCoordPlane(dir, max(), info);
122  outside = true;
123  }
124  else if (info.rawPoint()[dir] > bbMid[dir])
125  {
126  near[dir] = max()[dir];
127  }
128  else
129  {
130  near[dir] = min()[dir];
131  }
132  }
133 
134 
135  // For outside points the info will be correct now. Handle inside points
136  // using the three near distances. Project onto the nearest plane.
137  if (!outside)
138  {
139  const vector dist(cmptMag(info.point() - near));
140 
141  direction projNorm(vector::Z);
142 
143  if (dist.x() < dist.y())
144  {
145  if (dist.x() < dist.z())
146  {
147  projNorm = vector::X;
148  }
149  }
150  else
151  {
152  if (dist.y() < dist.z())
153  {
154  projNorm = vector::Y;
155  }
156  }
157 
158  projectOntoCoordPlane(projNorm, near, info);
159  }
160 
161 
162  // Check if outside. Optimisation: could do some checks on distance already
163  // on components above
164  if (magSqr(info.rawPoint() - sample) > nearestDistSqr)
165  {
166  info.setMiss();
167  info.setIndex(-1);
168  }
169 
170  return info;
171 }
172 
173 
174 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
175 
176 Foam::searchableBox::searchableBox
177 (
178  const IOobject& io,
179  const treeBoundBox& bb
180 )
181 :
182  searchableSurface(io),
183  treeBoundBox(bb)
184 {
185  if (!treeBoundBox::valid())
186  {
188  << "Illegal bounding box specification : "
189  << static_cast<const treeBoundBox>(*this) << nl
190  << exit(FatalError);
191  }
192 
193  bounds() = static_cast<treeBoundBox>(*this);
194 }
195 
196 
197 Foam::searchableBox::searchableBox
198 (
199  const IOobject& io,
200  const dictionary& dict
201 )
202 :
204  (
205  io,
206  treeBoundBox(dict.get<point>("min"), dict.get<point>("max"))
207  )
208 {}
209 
210 
211 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
212 
214 {
215  if (regions_.empty())
216  {
217  regions_.setSize(1);
218  regions_[0] = "region0";
219  }
220  return regions_;
221 }
222 
223 
225 {
226  auto tctrs = tmp<pointField>::New(6);
227  auto& ctrs = tctrs.ref();
228 
229  const pointField pts(treeBoundBox::points());
230  const faceList& fcs = treeBoundBox::faces;
231 
232  forAll(fcs, i)
233  {
234  ctrs[i] = fcs[i].centre(pts);
235  }
236 
237  return tctrs;
238 }
239 
240 
242 (
243  pointField& centres,
244  scalarField& radiusSqr
245 ) const
246 {
247  centres.setSize(size());
248  radiusSqr.setSize(size());
249  radiusSqr = Zero;
250 
251  const pointField pts(treeBoundBox::points());
252  const faceList& fcs = treeBoundBox::faces;
253 
254  forAll(fcs, i)
255  {
256  const face& f = fcs[i];
257 
258  centres[i] = f.centre(pts);
259  for (const label pointi : f)
260  {
261  const point& pt = pts[pointi];
262 
263  radiusSqr[i] = Foam::max
264  (
265  radiusSqr[i],
266  Foam::magSqr(pt-centres[i])
267  );
268  }
269  }
270 
271  // Add a bit to make sure all points are tested inside
272  radiusSqr += Foam::sqr(SMALL);
273 }
274 
275 
277 {
278  return treeBoundBox::points();
279 }
280 
281 
282 Foam::pointIndexHit Foam::searchableBox::findNearest
283 (
284  const point& sample,
285  const scalar nearestDistSqr
286 ) const
287 {
288  return findNearest(centre(), sample, nearestDistSqr);
289 }
290 
291 
293 (
294  const point& sample,
295  const scalar nearestDistSqr
296 ) const
297 {
298  const point bbMid(centre());
299 
300  // Outside point projected onto cube. Assume faces 0..5.
301  pointIndexHit info(true, sample, -1);
302  bool outside = false;
303 
304  // (for internal points) per direction what nearest cube side is
305  point near;
306 
307  for (direction dir = 0; dir < vector::nComponents; ++dir)
308  {
309  if (info.rawPoint()[dir] < min()[dir])
310  {
311  projectOntoCoordPlane(dir, min(), info);
312  outside = true;
313  }
314  else if (info.rawPoint()[dir] > max()[dir])
315  {
316  projectOntoCoordPlane(dir, max(), info);
317  outside = true;
318  }
319  else if (info.rawPoint()[dir] > bbMid[dir])
320  {
321  near[dir] = max()[dir];
322  }
323  else
324  {
325  near[dir] = min()[dir];
326  }
327  }
328 
329 
330  // For outside points the info will be correct now. Handle inside points
331  // using the three near distances. Project onto the nearest two planes.
332  if (!outside)
333  {
334  // Get the per-component distance to nearest wall
335  vector dist(cmptMag(info.rawPoint() - near));
336 
337  SortableList<scalar> sortedDist(3);
338  sortedDist[0] = dist[0];
339  sortedDist[1] = dist[1];
340  sortedDist[2] = dist[2];
341  sortedDist.sort();
342 
343  // Project onto nearest
344  projectOntoCoordPlane(sortedDist.indices()[0], near, info);
345  // Project onto second nearest
346  projectOntoCoordPlane(sortedDist.indices()[1], near, info);
347  }
348 
349 
350  // Check if outside. Optimisation: could do some checks on distance already
351  // on components above
352  if (magSqr(info.rawPoint() - sample) > nearestDistSqr)
353  {
354  info.setMiss();
355  info.setIndex(-1);
356  }
357 
358  return info;
359 }
360 
361 
362 Foam::pointIndexHit Foam::searchableBox::findNearest
363 (
364  const linePointRef& ln,
365  treeBoundBox& tightest,
366  point& linePoint
367 ) const
368 {
370  return pointIndexHit();
371 }
372 
373 
375 (
376  const point& start,
377  const point& end
378 ) const
379 {
380  pointIndexHit info(false, start, -1);
381 
382  bool foundInter;
383 
384  if (posBits(start) == 0)
385  {
386  if (posBits(end) == 0)
387  {
388  // Both start and end inside.
389  foundInter = false;
390  }
391  else
392  {
393  // end is outside. Clip to bounding box.
394  foundInter = intersects(end, start, info.rawPoint());
395  }
396  }
397  else
398  {
399  // start is outside. Clip to bounding box.
400  foundInter = intersects(start, end, info.rawPoint());
401  }
402 
403 
404  // Classify point
405  if (foundInter)
406  {
407  info.setHit();
408 
409  for (direction dir = 0; dir < vector::nComponents; ++dir)
410  {
411  if (info.rawPoint()[dir] == min()[dir])
412  {
413  info.setIndex(2*dir);
414  break;
415  }
416  else if (info.rawPoint()[dir] == max()[dir])
417  {
418  info.setIndex(2*dir+1);
419  break;
420  }
421  }
422 
423  if (info.index() == -1)
424  {
426  << "point " << info.rawPoint()
427  << " on segment " << start << end
428  << " should be on face of " << *this
429  << " but it isn't." << abort(FatalError);
430  }
431  }
432 
433  return info;
434 }
435 
436 
438 (
439  const point& start,
440  const point& end
441 ) const
442 {
443  return findLine(start, end);
444 }
445 
446 
447 void Foam::searchableBox::findNearest
448 (
449  const pointField& samples,
450  const scalarField& nearestDistSqr,
451  List<pointIndexHit>& info
452 ) const
453 {
454  info.setSize(samples.size());
455 
456  const point bbMid(centre());
457 
458  forAll(samples, i)
459  {
460  info[i] = findNearest(bbMid, samples[i], nearestDistSqr[i]);
461  }
462 }
463 
464 
466 (
467  const pointField& start,
468  const pointField& end,
469  List<pointIndexHit>& info
470 ) const
471 {
472  info.setSize(start.size());
473 
474  forAll(start, i)
475  {
476  info[i] = findLine(start[i], end[i]);
477  }
478 }
479 
480 
482 (
483  const pointField& start,
484  const pointField& end,
485  List<pointIndexHit>& info
486 ) const
487 {
488  info.setSize(start.size());
489 
490  forAll(start, i)
491  {
492  info[i] = findLineAny(start[i], end[i]);
493  }
494 }
495 
496 
498 (
499  const pointField& start,
500  const pointField& end,
502 ) const
503 {
504  info.setSize(start.size());
505 
506  // Work array
508 
509  // Tolerances:
510  // To find all intersections we add a small vector to the last intersection
511  // This is chosen such that
512  // - it is significant (SMALL is smallest representative relative tolerance;
513  // we need something bigger since we're doing calculations)
514  // - if the start-end vector is zero we still progress
515  const vectorField dirVec(end-start);
516  const scalarField magSqrDirVec(magSqr(dirVec));
517  const vectorField smallVec
518  (
519  ROOTSMALL*dirVec + vector::uniform(ROOTVSMALL)
520  );
521 
522  forAll(start, pointi)
523  {
524  // See if any intersection between pt and end
525  pointIndexHit inter = findLine(start[pointi], end[pointi]);
526 
527  if (inter.hit())
528  {
529  hits.clear();
530  hits.append(inter);
531 
532  point pt = inter.hitPoint() + smallVec[pointi];
533 
534  while (((pt-start[pointi])&dirVec[pointi]) <= magSqrDirVec[pointi])
535  {
536  // See if any intersection between pt and end
537  pointIndexHit inter = findLine(pt, end[pointi]);
538 
539  // Check for not hit or hit same face as before (can happen
540  // if vector along surface of face)
541  if
542  (
543  !inter.hit()
544  || (inter.index() == hits.last().index())
545  )
546  {
547  break;
548  }
549  hits.append(inter);
550 
551  pt = inter.hitPoint() + smallVec[pointi];
552  }
553 
554  info[pointi].transfer(hits);
555  }
556  else
557  {
558  info[pointi].clear();
559  }
560  }
561 }
562 
563 
565 (
566  const List<pointIndexHit>& info,
567  labelList& region
568 ) const
569 {
570  region.setSize(info.size());
571  region = 0;
572 }
573 
574 
576 (
577  const List<pointIndexHit>& info,
578  vectorField& normal
579 ) const
580 {
581  normal.setSize(info.size());
582  normal = Zero;
583 
584  forAll(info, i)
585  {
586  if (info[i].hit())
587  {
588  normal[i] = treeBoundBox::faceNormals[info[i].index()];
589  }
590  else
591  {
592  // Set to what?
593  }
594  }
595 }
596 
597 
599 (
600  const pointField& points,
601  List<volumeType>& volType
602 ) const
603 {
604  volType.setSize(points.size());
605 
606  forAll(points, pointi)
607  {
608  const point& pt = points[pointi];
609 
611 
612  for (direction dir=0; dir < vector::nComponents; ++dir)
613  {
614  if (pt[dir] < min()[dir] || pt[dir] > max()[dir])
615  {
616  vt = volumeType::OUTSIDE;
617  break;
618  }
619  }
620 
621  volType[pointi] = vt;
622  }
623 }
624 
625 
626 // ************************************************************************* //
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::Vector::Z
Definition: Vector.H:81
Foam::searchableBox::findLine
pointIndexHit findLine(const point &start, const point &end) const
Find nearest intersection of line between start and end.
Definition: searchableBox.C:375
Foam::Vector::Y
Definition: Vector.H:81
Foam::VectorSpace< Vector< Cmpt >, Cmpt, 3 >::uniform
static Vector< Cmpt > uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:164
Foam::PointIndexHit::setIndex
void setIndex(const label index) noexcept
Set the index.
Definition: PointIndexHit.H:211
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:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:55
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:438
Foam::treeBoundBox
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:86
Foam::PointIndexHit::setMiss
void setMiss() noexcept
Set the hit status off.
Definition: PointIndexHit.H:199
Foam::PointIndexHit::hitPoint
const point_type & hitPoint() const
Return hit point. Fatal if not hit.
Definition: PointIndexHit.H:154
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
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:540
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:296
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:445
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:565
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:52
Foam::PointIndexHit::hit
bool hit() const noexcept
Is there a hit?
Definition: PointIndexHit.H:130
Foam::searchableBox::findNearestOnEdge
pointIndexHit findNearestOnEdge(const point &sample, const scalar nearestDistSqr) const
Calculate nearest point on edge.
Definition: searchableBox.C:293
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::treeBoundBox::faces
static const faceList faces
Face to point addressing.
Definition: treeBoundBox.H:151
Foam::Vector::X
Definition: Vector.H:81
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:474
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:224
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:459
Foam::addNamedToRunTimeSelectionTable
addNamedToRunTimeSelectionTable(topoSetCellSource, badQualityToCell, word, badQuality)
samples
scalarField samples(nIntervals, Zero)
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:276
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::PointIndexHit::rawPoint
const point_type & rawPoint() const noexcept
The point, no checks. Same as point()
Definition: PointIndexHit.H:178
Foam::SortableList
A list that is sorted upon construction or when explicitly requested with the sort() method.
Definition: List.H:63
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:121
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
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
Foam::searchableBox::getVolumeType
virtual void getVolumeType(const pointField &points, List< volumeType > &volType) const
Determine type (inside/outside) for points.
Definition: searchableBox.C:599
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:381
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::pointIndexHit
PointIndexHit< point > pointIndexHit
A PointIndexHit for 3D points.
Definition: pointIndexHit.H:46
Foam::nl
constexpr char nl
Definition: Ostream.H:385
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List< word >
Foam::PointIndexHit::index
label index() const noexcept
Return the hit index.
Definition: PointIndexHit.H:136
Foam::PointIndexHit::setHit
void setHit() noexcept
Set the hit status on.
Definition: PointIndexHit.H:193
Foam::searchableBox::regions
virtual const wordList & regions() const
Names of regions.
Definition: searchableBox.C:213
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:498
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:52
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::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:72
Foam::searchableBox
Searching on bounding box.
Definition: searchableBox.H:83
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:576
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:242
Foam::VectorSpace< Vector< Cmpt >, Cmpt, 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