treeBoundBox.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-2024 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 "treeBoundBox.H"
27 #include "ListOps.H"
28 #include "OFstream.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 const Foam::scalar Foam::treeBoundBox::great(great);
33 
35 (
36  vector(-great, -great, -great),
37  vector(great, great, great)
38 );
39 
40 
42 (
43  vector(great, great, great),
44  vector(-great, -great, -great)
45 );
46 
47 
48 const Foam::label facesArray[6][4] =
49 {
50  {0, 4, 6, 2}, // left
51  {1, 3, 7, 5}, // right
52  {0, 1, 5, 4}, // bottom
53  {2, 6, 7, 3}, // top
54  {0, 2, 3, 1}, // back
55  {4, 5, 7, 6} // front
56 };
57 
58 
60 (
61  initListList<face, label, 6, 4>(facesArray)
62 );
63 
64 
65 const Foam::label edgesArray[12][2] =
66 {
67  {0, 1}, // 0
68  {1, 3},
69  {2, 3}, // 2
70  {0, 2},
71  {4, 5}, // 4
72  {5, 7},
73  {6, 7}, // 6
74  {4, 6},
75  {0, 4}, // 8
76  {1, 5},
77  {3, 7}, // 10
78  {2, 6}
79 };
80 
81 
83 (
84  calcEdges(edgesArray)
85 );
86 
87 
89 (
90  calcFaceNormals()
91 );
92 
93 
94 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
95 
96 Foam::edgeList Foam::treeBoundBox::calcEdges(const label edgesArray[12][2])
97 {
98  edgeList edges(12);
99  forAll(edges, edgeI)
100  {
101  edges[edgeI][0] = edgesArray[edgeI][0];
102  edges[edgeI][1] = edgesArray[edgeI][1];
103  }
104  return edges;
105 }
106 
107 
108 Foam::FixedList<Foam::vector, 6> Foam::treeBoundBox::calcFaceNormals()
109 {
110  FixedList<vector, 6> normals;
111  normals[faceId::left] = vector(-1, 0, 0);
112  normals[faceId::right] = vector( 1, 0, 0);
113  normals[faceId::bottom] = vector( 0, -1, 0);
114  normals[faceId::top] = vector( 0, 1, 0);
115  normals[faceId::back] = vector( 0, 0, -1);
116  normals[faceId::front] = vector( 0, 0, 1);
117  return normals;
118 }
119 
120 
121 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
122 
124 :
125  boundBox(points, false)
126 {
127  if (points.empty())
128  {
130  << "cannot find bounding box for zero-sized pointField, "
131  << "returning zero" << endl;
132 
133  return;
134  }
135 }
136 
137 
139 (
140  const UList<point>& points,
141  const labelUList& indices
142 )
143 :
144  boundBox(points, indices, false)
145 {
146  if (points.empty() || indices.empty())
147  {
149  << "cannot find bounding box for zero-sized pointField, "
150  << "returning zero" << endl;
151 
152  return;
153  }
154 }
155 
156 
157 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
158 
160 {
162  pointField& points = tPts.ref();
163 
164  forAll(points, octant)
165  {
166  points[octant] = corner(octant);
167  }
168 
169  return tPts;
170 }
171 
172 
174 {
175  return subBbox(midpoint(), octant);
176 }
177 
178 
180 (
181  const point& mid,
182  const direction octant
183 ) const
184 {
185  if (octant > 7)
186  {
188  << "octant should be [0..7]"
189  << abort(FatalError);
190  }
191 
192  // start with a copy of this bounding box and adjust limits accordingly
193  treeBoundBox subBb(*this);
194  point& bbMin = subBb.min();
195  point& bbMax = subBb.max();
196 
197  if (octant & octantBit::rightHalf)
198  {
199  bbMin.x() = mid.x(); // mid -> max
200  }
201  else
202  {
203  bbMax.x() = mid.x(); // min -> mid
204  }
205 
206  if (octant & octantBit::topHalf)
207  {
208  bbMin.y() = mid.y(); // mid -> max
209  }
210  else
211  {
212  bbMax.y() = mid.y(); // min -> mid
213  }
214 
215  if (octant & octantBit::frontHalf)
216  {
217  bbMin.z() = mid.z(); // mid -> max
218  }
219  else
220  {
221  bbMax.z() = mid.z(); // min -> mid
222  }
223 
224  return subBb;
225 }
226 
227 
229 (
230  const point& overallStart,
231  const vector& overallVec,
232  const point& start,
233  const point& end,
234  point& pt,
235  direction& ptOnFaces
236 ) const
237 {
238  // Sutherlands algorithm:
239  // loop
240  // - start = intersection of line with one of the planes bounding
241  // the bounding box
242  // - stop if start inside bb (return true)
243  // - stop if start and end in same 'half' (e.g. both above bb)
244  // (return false)
245  //
246  // Uses posBits to efficiently determine 'half' in which start and end
247  // point are.
248  //
249  // Note:
250  // - sets coordinate to exact position: e.g. pt.x() = min().x()
251  // since plane intersect routine might have truncation error.
252  // This makes sure that posBits tests 'inside'
253 
254  const direction endBits = posBits(end);
255  pt = start;
256 
257  // Allow maximum of 3 clips.
258  for (label i = 0; i < 4; ++i)
259  {
260  direction ptBits = posBits(pt);
261 
262  if (ptBits == 0)
263  {
264  // pt inside bb
265  ptOnFaces = faceBits(pt);
266  return true;
267  }
268 
269  if ((ptBits & endBits) != 0)
270  {
271  // pt and end in same block outside of bb
272  ptOnFaces = faceBits(pt);
273  return false;
274  }
275 
276  if (ptBits & faceBit::left)
277  {
278  // Intersect with plane V=min, n=-1,0,0
279  if (Foam::mag(overallVec.x()) > vSmall)
280  {
281  scalar s = (min().x() - overallStart.x())/overallVec.x();
282  pt.x() = min().x();
283  pt.y() = overallStart.y() + overallVec.y()*s;
284  pt.z() = overallStart.z() + overallVec.z()*s;
285  }
286  else
287  {
288  // Vector not in x-direction. But still intersecting bb planes.
289  // So must be close - just snap to bb.
290  pt.x() = min().x();
291  }
292  }
293  else if (ptBits & faceBit::right)
294  {
295  // Intersect with plane V=max, n=1,0,0
296  if (Foam::mag(overallVec.x()) > vSmall)
297  {
298  scalar s = (max().x() - overallStart.x())/overallVec.x();
299  pt.x() = max().x();
300  pt.y() = overallStart.y() + overallVec.y()*s;
301  pt.z() = overallStart.z() + overallVec.z()*s;
302  }
303  else
304  {
305  pt.x() = max().x();
306  }
307  }
308  else if (ptBits & faceBit::bottom)
309  {
310  // Intersect with plane V=min, n=0,-1,0
311  if (Foam::mag(overallVec.y()) > vSmall)
312  {
313  scalar s = (min().y() - overallStart.y())/overallVec.y();
314  pt.x() = overallStart.x() + overallVec.x()*s;
315  pt.y() = min().y();
316  pt.z() = overallStart.z() + overallVec.z()*s;
317  }
318  else
319  {
320  pt.x() = min().y();
321  }
322  }
323  else if (ptBits & faceBit::top)
324  {
325  // Intersect with plane V=max, n=0,1,0
326  if (Foam::mag(overallVec.y()) > vSmall)
327  {
328  scalar s = (max().y() - overallStart.y())/overallVec.y();
329  pt.x() = overallStart.x() + overallVec.x()*s;
330  pt.y() = max().y();
331  pt.z() = overallStart.z() + overallVec.z()*s;
332  }
333  else
334  {
335  pt.y() = max().y();
336  }
337  }
338  else if (ptBits & faceBit::back)
339  {
340  // Intersect with plane V=min, n=0,0,-1
341  if (Foam::mag(overallVec.z()) > vSmall)
342  {
343  scalar s = (min().z() - overallStart.z())/overallVec.z();
344  pt.x() = overallStart.x() + overallVec.x()*s;
345  pt.y() = overallStart.y() + overallVec.y()*s;
346  pt.z() = min().z();
347  }
348  else
349  {
350  pt.z() = min().z();
351  }
352  }
353  else if (ptBits & faceBit::front)
354  {
355  // Intersect with plane V=max, n=0,0,1
356  if (Foam::mag(overallVec.z()) > vSmall)
357  {
358  scalar s = (max().z() - overallStart.z())/overallVec.z();
359  pt.x() = overallStart.x() + overallVec.x()*s;
360  pt.y() = overallStart.y() + overallVec.y()*s;
361  pt.z() = max().z();
362  }
363  else
364  {
365  pt.z() = max().z();
366  }
367  }
368  }
369 
370  // Can end up here if the end point is on the edge of the boundBox
371  return true;
372 }
373 
374 
376 (
377  const point& start,
378  const point& end,
379  point& pt
380 ) const
381 {
382  direction ptBits;
383  return intersects(start, end-start, start, end, pt, ptBits);
384 }
385 
386 
387 bool Foam::treeBoundBox::contains(const vector& dir, const point& pt) const
388 {
389  // Compare all components against min and max of bb
390 
391  for (direction cmpt=0; cmpt<3; cmpt++)
392  {
393  if (pt[cmpt] < min()[cmpt])
394  {
395  return false;
396  }
397  else if (pt[cmpt] == min()[cmpt])
398  {
399  // On edge. Outside if direction points outwards.
400  if (dir[cmpt] < 0)
401  {
402  return false;
403  }
404  }
405 
406  if (pt[cmpt] > max()[cmpt])
407  {
408  return false;
409  }
410  else if (pt[cmpt] == max()[cmpt])
411  {
412  // On edge. Outside if direction points outwards.
413  if (dir[cmpt] > 0)
414  {
415  return false;
416  }
417  }
418  }
419 
420  // All components inside bb
421  return true;
422 }
423 
424 
426 {
427  direction faceBits = 0;
428  if (pt.x() == min().x())
429  {
430  faceBits |= faceBit::left;
431  }
432  else if (pt.x() == max().x())
433  {
434  faceBits |= faceBit::right;
435  }
436 
437  if (pt.y() == min().y())
438  {
439  faceBits |= faceBit::bottom;
440  }
441  else if (pt.y() == max().y())
442  {
443  faceBits |= faceBit::top;
444  }
445 
446  if (pt.z() == min().z())
447  {
448  faceBits |= faceBit::back;
449  }
450  else if (pt.z() == max().z())
451  {
452  faceBits |= faceBit::front;
453  }
454  return faceBits;
455 }
456 
457 
459 {
460  direction posBits = 0;
461 
462  if (pt.x() < min().x())
463  {
464  posBits |= faceBit::left;
465  }
466  else if (pt.x() > max().x())
467  {
468  posBits |= faceBit::right;
469  }
470 
471  if (pt.y() < min().y())
472  {
473  posBits |= faceBit::bottom;
474  }
475  else if (pt.y() > max().y())
476  {
477  posBits |= faceBit::top;
478  }
479 
480  if (pt.z() < min().z())
481  {
482  posBits |= faceBit::back;
483  }
484  else if (pt.z() > max().z())
485  {
486  posBits |= faceBit::front;
487  }
488  return posBits;
489 }
490 
491 
493 (
494  const point& pt,
495  point& nearest,
496  point& furthest
497 ) const
498 {
499  scalar nearX, nearY, nearZ;
500  scalar farX, farY, farZ;
501 
502  if (Foam::mag(min().x() - pt.x()) < Foam::mag(max().x() - pt.x()))
503  {
504  nearX = min().x();
505  farX = max().x();
506  }
507  else
508  {
509  nearX = max().x();
510  farX = min().x();
511  }
512 
513  if (Foam::mag(min().y() - pt.y()) < Foam::mag(max().y() - pt.y()))
514  {
515  nearY = min().y();
516  farY = max().y();
517  }
518  else
519  {
520  nearY = max().y();
521  farY = min().y();
522  }
523 
524  if (Foam::mag(min().z() - pt.z()) < Foam::mag(max().z() - pt.z()))
525  {
526  nearZ = min().z();
527  farZ = max().z();
528  }
529  else
530  {
531  nearZ = max().z();
532  farZ = min().z();
533  }
534 
535  nearest = point(nearX, nearY, nearZ);
536  furthest = point(farX, farY, farZ);
537 }
538 
539 
540 Foam::scalar Foam::treeBoundBox::maxDist(const point& pt) const
541 {
542  point near, far;
543  calcExtremities(pt, near, far);
544 
545  return Foam::mag(far - pt);
546 }
547 
548 
550 (
551  const point& pt,
552  const treeBoundBox& other
553 ) const
554 {
555  //
556  // Distance point <-> nearest and furthest away vertex of this
557  //
558 
559  point nearThis, farThis;
560 
561  // get nearest and furthest away vertex
562  calcExtremities(pt, nearThis, farThis);
563 
564  const scalar minDistThis =
565  sqr(nearThis.x() - pt.x())
566  + sqr(nearThis.y() - pt.y())
567  + sqr(nearThis.z() - pt.z());
568  const scalar maxDistThis =
569  sqr(farThis.x() - pt.x())
570  + sqr(farThis.y() - pt.y())
571  + sqr(farThis.z() - pt.z());
572 
573  //
574  // Distance point <-> other
575  //
576 
577  point nearOther, farOther;
578 
579  // get nearest and furthest away vertex
580  other.calcExtremities(pt, nearOther, farOther);
581 
582  const scalar minDistOther =
583  sqr(nearOther.x() - pt.x())
584  + sqr(nearOther.y() - pt.y())
585  + sqr(nearOther.z() - pt.z());
586  const scalar maxDistOther =
587  sqr(farOther.x() - pt.x())
588  + sqr(farOther.y() - pt.y())
589  + sqr(farOther.z() - pt.z());
590 
591  //
592  // Categorise
593  //
594  if (maxDistThis < minDistOther)
595  {
596  // All vertices of this are nearer to point than any vertex of other
597  return -1;
598  }
599  else if (minDistThis > maxDistOther)
600  {
601  // All vertices of this are further from point than any vertex of other
602  return 1;
603  }
604  else
605  {
606  // Mixed bag
607  return 0;
608  }
609 }
610 
611 
612 void Foam::treeBoundBox::writeOBJ(const fileName& fName) const
613 {
614  OFstream str(fName);
615 
616  Info<< "Dumping bounding box " << *this << " as lines to obj file "
617  << str.name() << endl;
618 
619  pointField bbPoints(points());
620 
621  forAll(bbPoints, i)
622  {
623  const point& pt = bbPoints[i];
624 
625  str<< "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << endl;
626  }
627 
629  {
630  const edge& e = treeBoundBox::edges[i];
631 
632  str<< "l " << e[0]+1 << ' ' << e[1]+1 << nl;
633  }
634 }
635 
636 
637 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
638 
640 {
641  return operator==
642  (
643  static_cast<const boundBox&>(a),
644  static_cast<const boundBox&>(b)
645  );
646 }
647 
648 
650 {
651  return !(a == b);
652 }
653 
654 
655 // * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
656 
658 {
659  return os << static_cast<const boundBox&>(bb);
660 }
661 
662 
664 {
665  return is >> static_cast<boundBox&>(bb);
666 }
667 
668 
669 // ************************************************************************* //
scalar y
Various functions to operate on Lists.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:78
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
Output to file stream.
Definition: OFstream.H:86
const fileName & name() const
Return the name of the stream.
Definition: OFstream.H:120
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:74
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:325
const Cmpt & z() const
Definition: VectorI.H:87
const Cmpt & y() const
Definition: VectorI.H:81
const Cmpt & x() const
Definition: VectorI.H:75
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
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
A class for handling file names.
Definition: fileName.H:82
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:181
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:90
void writeOBJ(const fileName &fName) const
Definition: treeBoundBox.C:612
void calcExtremities(const point &pt, point &nearest, point &furthest) const
Calculate nearest and furthest (to point) vertex coords of.
Definition: treeBoundBox.C:493
scalar maxDist(const point &) const
Returns distance point to furthest away corner.
Definition: treeBoundBox.C:540
label distanceCmp(const point &, const treeBoundBox &other) const
Compare distance to point with other bounding box.
Definition: treeBoundBox.C:550
direction faceBits(const point &) const
Code position of point on bounding box faces.
Definition: treeBoundBox.C:425
treeBoundBox subBbox(const direction) const
Sub box given by octant number. Midpoint calculated.
Definition: treeBoundBox.C:173
static const scalar great
The great value used for greatBox and invertedBox.
Definition: treeBoundBox.H:102
static const edgeList edges
Edge to point addressing.
Definition: treeBoundBox.H:178
bool contains(const point &) const
Contains point or other bounding box?
Definition: boundBoxI.H:176
bool intersects(const point &overallStart, const vector &overallVec, const point &start, const point &end, point &pt, direction &ptBits) const
Intersects segment; set point to intersection position and face,.
Definition: treeBoundBox.C:229
static const treeBoundBox invertedBox
As per boundBox::invertedBox, but with great instead of vGreat.
Definition: treeBoundBox.H:108
static const treeBoundBox greatBox
As per boundBox::greatBox, but with great instead of vGreat.
Definition: treeBoundBox.H:105
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
direction posBits(const point &) const
Position of point relative to bounding box.
Definition: treeBoundBox.C:458
treeBoundBox()
Construct null setting points to zero.
Definition: treeBoundBoxI.H:31
tmp< pointField > points() const
Vertex coordinates. In octant coding.
Definition: treeBoundBox.C:159
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const pointField & points
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
volScalarField & b
Definition: createFields.H:25
#define WarningInFunction
Report a warning using Foam::Warning.
bool operator!=(const particle &, const particle &)
Definition: particle.C:1210
const doubleScalar e
Definition: doubleScalar.H:106
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:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
errorManip< error > abort(error &err)
Definition: errorManip.H:131
messageStream Info
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
Istream & operator>>(Istream &, pistonPointEdgeData &)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
error FatalError
List< edge > edgeList
Definition: edgeList.H:38
static const char nl
Definition: Ostream.H:266
uint8_t direction
Definition: direction.H:45
const Foam::label edgesArray[12][2]
Definition: treeBoundBox.C:65
const Foam::label facesArray[6][4]
Definition: treeBoundBox.C:48