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