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 | Copyright (C) 2011-2015 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 "searchableBox.H"
28 #include "SortableList.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34  defineTypeNameAndDebug(searchableBox, 0);
35  addToRunTimeSelectionTable(searchableSurface, searchableBox, dict);
36 }
37 
38 
39 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
40 
41 void Foam::searchableBox::projectOntoCoordPlane
42 (
43  const direction dir,
44  const point& planePt,
45  pointIndexHit& info
46 ) const
47 {
48  // Set point
49  info.rawPoint()[dir] = planePt[dir];
50  // Set face
51  if (planePt[dir] == min()[dir])
52  {
53  info.setIndex(dir*2);
54  }
55  else if (planePt[dir] == max()[dir])
56  {
57  info.setIndex(dir*2+1);
58  }
59  else
60  {
61  FatalErrorIn("searchableBox::projectOntoCoordPlane(..)")
62  << "Point on plane " << planePt
63  << " is not on coordinate " << min()[dir]
64  << " nor " << max()[dir] << abort(FatalError);
65  }
66 }
67 
68 
69 // Returns miss or hit with face (0..5) and region(always 0)
70 Foam::pointIndexHit Foam::searchableBox::findNearest
71 (
72  const point& bbMid,
73  const point& sample,
74  const scalar nearestDistSqr
75 ) const
76 {
77  // Point can be inside or outside. For every component direction can be
78  // left of min, right of max or inbetween.
79  // - outside points: project first one x plane (either min().x()
80  // or max().x()), then onto y plane and finally z. You should be left
81  // with intersection point
82  // - inside point: find nearest side (compare to mid point). Project onto
83  // that.
84 
85  // The face is set to the last projected face.
86 
87 
88  // Outside point projected onto cube. Assume faces 0..5.
89  pointIndexHit info(true, sample, -1);
90  bool outside = false;
91 
92  // (for internal points) per direction what nearest cube side is
93  point near;
94 
95  for (direction dir = 0; dir < vector::nComponents; dir++)
96  {
97  if (info.rawPoint()[dir] < min()[dir])
98  {
99  projectOntoCoordPlane(dir, min(), info);
100  outside = true;
101  }
102  else if (info.rawPoint()[dir] > max()[dir])
103  {
104  projectOntoCoordPlane(dir, max(), info);
105  outside = true;
106  }
107  else if (info.rawPoint()[dir] > bbMid[dir])
108  {
109  near[dir] = max()[dir];
110  }
111  else
112  {
113  near[dir] = min()[dir];
114  }
115  }
116 
117 
118  // For outside points the info will be correct now. Handle inside points
119  // using the three near distances. Project onto the nearest plane.
120  if (!outside)
121  {
122  vector dist(cmptMag(info.rawPoint() - near));
123 
124  if (dist.x() < dist.y())
125  {
126  if (dist.x() < dist.z())
127  {
128  // Project onto x plane
129  projectOntoCoordPlane(vector::X, near, info);
130  }
131  else
132  {
133  projectOntoCoordPlane(vector::Z, near, info);
134  }
135  }
136  else
137  {
138  if (dist.y() < dist.z())
139  {
140  projectOntoCoordPlane(vector::Y, near, info);
141  }
142  else
143  {
144  projectOntoCoordPlane(vector::Z, near, info);
145  }
146  }
147  }
148 
149 
150  // Check if outside. Optimisation: could do some checks on distance already
151  // on components above
152  if (magSqr(info.rawPoint() - sample) > nearestDistSqr)
153  {
154  info.setMiss();
155  info.setIndex(-1);
156  }
157 
158  return info;
159 }
160 
161 
162 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
163 
164 Foam::searchableBox::searchableBox
165 (
166  const IOobject& io,
167  const treeBoundBox& bb
168 )
169 :
170  searchableSurface(io),
171  treeBoundBox(bb)
172 {
173  if (!contains(midpoint()))
174  {
176  (
177  "Foam::searchableBox::searchableBox\n"
178  "(\n"
179  " const IOobject& io,\n"
180  " const treeBoundBox& bb\n"
181  ")\n"
182  ) << "Illegal bounding box specification : "
183  << static_cast<const treeBoundBox>(*this) << exit(FatalError);
184  }
185 
186  bounds() = static_cast<boundBox>(*this);
187 }
188 
189 
190 Foam::searchableBox::searchableBox
191 (
192  const IOobject& io,
193  const dictionary& dict
194 )
195 :
196  searchableSurface(io),
197  treeBoundBox(dict.lookup("min"), dict.lookup("max"))
198 {
199  if (!contains(midpoint()))
200  {
202  (
203  "Foam::searchableBox::searchableBox\n"
204  "(\n"
205  " const IOobject& io,\n"
206  " const treeBoundBox& bb\n"
207  ")\n"
208  ) << "Illegal bounding box specification : "
209  << static_cast<const treeBoundBox>(*this) << exit(FatalError);
210  }
211 
212  bounds() = static_cast<boundBox>(*this);
213 }
214 
215 
216 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
217 
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 {
238  pointField& ctrs = tCtrs();
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  forAll(f, fp)
271  {
272  const point& pt = pts[f[fp]];
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(midpoint(), sample, nearestDistSqr);
300 }
301 
302 
304 (
305  const point& sample,
306  const scalar nearestDistSqr
307 ) const
308 {
309  const point bbMid(midpoint());
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  (
382  "searchableBox::findNearest"
383  "(const linePointRef&, treeBoundBox&, point&)"
384  );
385  return pointIndexHit();
386 }
387 
388 
390 (
391  const point& start,
392  const point& end
393 ) const
394 {
395  pointIndexHit info(false, start, -1);
396 
397  bool foundInter;
398 
399  if (posBits(start) == 0)
400  {
401  if (posBits(end) == 0)
402  {
403  // Both start and end inside.
404  foundInter = false;
405  }
406  else
407  {
408  // end is outside. Clip to bounding box.
409  foundInter = intersects(end, start, info.rawPoint());
410  }
411  }
412  else
413  {
414  // start is outside. Clip to bounding box.
415  foundInter = intersects(start, end, info.rawPoint());
416  }
417 
418 
419  // Classify point
420  if (foundInter)
421  {
422  info.setHit();
423 
424  for (direction dir = 0; dir < vector::nComponents; dir++)
425  {
426  if (info.rawPoint()[dir] == min()[dir])
427  {
428  info.setIndex(2*dir);
429  break;
430  }
431  else if (info.rawPoint()[dir] == max()[dir])
432  {
433  info.setIndex(2*dir+1);
434  break;
435  }
436  }
437 
438  if (info.index() == -1)
439  {
440  FatalErrorIn("searchableBox::findLine(const point&, const point&)")
441  << "point " << info.rawPoint()
442  << " on segment " << start << end
443  << " should be on face of " << *this
444  << " but it isn't." << abort(FatalError);
445  }
446  }
447 
448  return info;
449 }
450 
451 
453 (
454  const point& start,
455  const point& end
456 ) const
457 {
458  return findLine(start, end);
459 }
460 
461 
462 void Foam::searchableBox::findNearest
463 (
464  const pointField& samples,
465  const scalarField& nearestDistSqr,
466  List<pointIndexHit>& info
467 ) const
468 {
469  info.setSize(samples.size());
470 
471  const point bbMid(midpoint());
472 
473  forAll(samples, i)
474  {
475  info[i] = findNearest(bbMid, samples[i], nearestDistSqr[i]);
476  }
477 }
478 
479 
481 (
482  const pointField& start,
483  const pointField& end,
484  List<pointIndexHit>& info
485 ) const
486 {
487  info.setSize(start.size());
488 
489  forAll(start, i)
490  {
491  info[i] = findLine(start[i], end[i]);
492  }
493 }
494 
495 
497 (
498  const pointField& start,
499  const pointField& end,
500  List<pointIndexHit>& info
501 ) const
502 {
503  info.setSize(start.size());
504 
505  forAll(start, i)
506  {
507  info[i] = findLineAny(start[i], end[i]);
508  }
509 }
510 
511 
513 (
514  const pointField& start,
515  const pointField& end,
516  List<List<pointIndexHit> >& info
517 ) const
518 {
519  info.setSize(start.size());
520 
521  // Work array
523 
524  // Tolerances:
525  // To find all intersections we add a small vector to the last intersection
526  // This is chosen such that
527  // - it is significant (SMALL is smallest representative relative tolerance;
528  // we need something bigger since we're doing calculations)
529  // - if the start-end vector is zero we still progress
530  const vectorField dirVec(end-start);
531  const scalarField magSqrDirVec(magSqr(dirVec));
532  const vectorField smallVec
533  (
534  ROOTSMALL*dirVec
535  + vector(ROOTVSMALL,ROOTVSMALL,ROOTVSMALL)
536  );
537 
538  forAll(start, pointI)
539  {
540  // See if any intersection between pt and end
541  pointIndexHit inter = findLine(start[pointI], end[pointI]);
542 
543  if (inter.hit())
544  {
545  hits.clear();
546  hits.append(inter);
547 
548  point pt = inter.hitPoint() + smallVec[pointI];
549 
550  while (((pt-start[pointI])&dirVec[pointI]) <= magSqrDirVec[pointI])
551  {
552  // See if any intersection between pt and end
553  pointIndexHit inter = findLine(pt, end[pointI]);
554 
555  // Check for not hit or hit same face as before (can happen
556  // if vector along surface of face)
557  if
558  (
559  !inter.hit()
560  || (inter.index() == hits.last().index())
561  )
562  {
563  break;
564  }
565  hits.append(inter);
566 
567  pt = inter.hitPoint() + smallVec[pointI];
568  }
569 
570  info[pointI].transfer(hits);
571  }
572  else
573  {
574  info[pointI].clear();
575  }
576  }
577 }
578 
579 
581 (
582  const List<pointIndexHit>& info,
583  labelList& region
584 ) const
585 {
586  region.setSize(info.size());
587  region = 0;
588 }
589 
590 
592 (
593  const List<pointIndexHit>& info,
594  vectorField& normal
595 ) const
596 {
597  normal.setSize(info.size());
598  normal = vector::zero;
599 
600  forAll(info, i)
601  {
602  if (info[i].hit())
603  {
604  normal[i] = treeBoundBox::faceNormals[info[i].index()];
605  }
606  else
607  {
608  // Set to what?
609  }
610  }
611 }
612 
613 
615 (
616  const pointField& points,
617  List<volumeType>& volType
618 ) const
619 {
620  volType.setSize(points.size());
621  volType = volumeType::INSIDE;
622 
623  forAll(points, pointI)
624  {
625  const point& pt = points[pointI];
626 
627  for (direction dir = 0; dir < vector::nComponents; dir++)
628  {
629  if (pt[dir] < min()[dir] || pt[dir] > max()[dir])
630  {
631  volType[pointI] = volumeType::OUTSIDE;
632  break;
633  }
634  }
635  }
636 }
637 
638 
639 // ************************************************************************* //
unsigned char direction
Definition: direction.H:43
const labelList & indices() const
Return the list of sorted indices. Updated every sort.
Definition: SortableList.H:90
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
vector point
Point is a vector.
Definition: point.H:41
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:390
Base class of (analytical or triangulated) surface. Encapsulates all the search routines. WIP.
labelList f(nPoints)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
virtual ~searchableBox()
Destructor.
const Point & rawPoint() const
Return point with no checking.
T & last()
Return the last element of the list.
Definition: UListI.H:131
pointIndexHit findNearestOnEdge(const point &sample, const scalar nearestDistSqr) const
Calculate nearest point on edge.
void clear()
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:242
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:75
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:310
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
virtual void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit > > &) const
Get all intersections in order from start to end.
bool hit() const
Is there a hit.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
tmp< pointField > points() const
Vertex coordinates. In octant coding.
Definition: treeBoundBox.C:168
void sort()
(stable) sort the list (if changed after construction time)
Definition: SortableList.C:95
label index() const
Return index.
point centre(const pointField &) const
Centre point of face.
Definition: face.C:493
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Namespace for OpenFOAM.
A line primitive.
Definition: line.H:56
const Point & hitPoint() const
Return hit point.
static faceList faces()
Return faces with correct point order.
Definition: boundBox.C:167
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Number of components in this vector space.
Definition: VectorSpace.H:88
void clear()
Clear the list, i.e. set size to zero.
Definition: List.C:379
void setSize(const label)
Reset size of List.
Definition: List.C:318
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
virtual tmp< pointField > coordinates() const
Get representative set of element coordinates.
A list that is sorted upon construction or when explicitly requested with the sort() method...
Definition: List.H:65
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
pointIndexHit findLine(const point &start, const point &end) const
Find nearest intersection of line between start and end.
virtual const wordList & regions() const
Names of regions.
#define forAll(list, i)
Definition: UList.H:421
static const FixedList< vector, 6 > faceNormals
Per face the unit normal.
Definition: treeBoundBox.H:160
Macros for easy insertion into run-time selection tables.
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void setIndex(const label index)
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:55
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
error FatalError
static const Vector zero
Definition: Vector.H:80
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
virtual tmp< pointField > points() const
Get the points that define the surface.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
virtual void boundingSpheres(pointField &centres, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared), one per element.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
#define notImplemented(functionName)
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:356
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
virtual void getRegion(const List< pointIndexHit > &, labelList &region) const
From a set of points and indices get the region.
A class for managing temporary objects.
Definition: PtrList.H:118
virtual void getVolumeType(const pointField &, List< volumeType > &) const
Determine type (inside/outside/mixed) for point. unknown if.
defineTypeNameAndDebug(combustionModel, 0)
pointIndexHit findLineAny(const point &start, const point &end) const
Find any intersection of line between start and end.