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