indexedOctree.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) 2016-2022 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 "indexedOctree.H"
30 #include "line.H"
31 #include "OFstream.H"
32 #include "ListOps.H"
33 #include "memInfo.H"
34 
35 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 
37 template<class Type>
39 (
40  const labelUList& indices,
41  const treeBoundBox& bb,
42  FixedList<labelList, 8>& dividedIndices
43 ) const
44 {
45  // Scratch array
46  DynamicList<label> contains(indices.size());
47 
48  for (direction octant = 0; octant < 8; ++octant)
49  {
50  const treeBoundBox subBbs(bb.subBbox(octant));
51 
52  contains.clear();
53 
54  for (const label index : indices)
55  {
56  if (shapes_.overlaps(index, subBbs))
57  {
58  contains.push_back(index);
59  }
60  }
61 
62  // The sub-divided indices:
63  dividedIndices[octant] = contains;
64  }
65 }
66 
67 
68 template<class Type>
71 (
72  const treeBoundBox& bb,
73  DynamicList<labelList>& contents,
74  label contentIndex
75 ) const
76 {
77  if
78  (
79  bb.min().x() >= bb.max().x()
80  || bb.min().y() >= bb.max().y()
81  || bb.min().z() >= bb.max().z()
82  )
83  {
85  << "Badly formed bounding box:" << bb
86  << abort(FatalError);
87  }
88 
89  // Divide the indices into 8 (possibly empty) subsets.
90  // Replace current contentIndex with the first (non-empty) subset.
91  // Append the rest.
92 
93  const labelList& indices = contents[contentIndex];
94 
95  FixedList<labelList, 8> dividedIndices;
96  divide(indices, bb, dividedIndices);
97 
98  node nod;
99  nod.bb_ = bb;
100  nod.parent_ = -1;
101 
102  bool replaceNode = true;
103 
104  for (direction octant = 0; octant < 8; ++octant)
105  {
106  auto& subIndices = dividedIndices[octant];
107 
108  if (subIndices.size())
109  {
110  if (replaceNode)
111  {
112  // Replace existing
113  contents[contentIndex] = std::move(subIndices);
114  replaceNode = false;
115  }
116  else
117  {
118  // Append to contents
119  contentIndex = contents.size();
120  contents.push_back(std::move(subIndices));
121  }
122 
123  nod.subNodes_[octant] = contentPlusOctant(contentIndex, octant);
124  }
125  else
126  {
127  // Mark node as empty
128  nod.subNodes_[octant] = emptyPlusOctant(octant);
129  }
130  }
131 
132  return nod;
133 }
134 
135 
136 template<class Type>
138 (
139  const label minSize,
140  DynamicList<indexedOctreeBase::node>& nodes,
141  DynamicList<labelList>& contents
142 ) const
143 {
144  label currentSize = nodes.size();
145 
146  // Take care to loop only over old nodes.
147  // Also we loop over the same DynamicList which gets modified and
148  // moved so make sure not to keep any references!
149  for (label nodeI = 0; nodeI < currentSize; nodeI++)
150  {
151  for (direction octant = 0; octant < node::nChildren; ++octant)
152  {
153  labelBits index = nodes[nodeI].subNodes_[octant];
154 
155  if (isNode(index))
156  {
157  // tree node. Leave intact.
158  }
159  else if (isContent(index))
160  {
161  label contentI = getContent(index);
162 
163  if (contents[contentI].size() > minSize)
164  {
165  // Create node for content.
166 
167  // Find the bounding box for the subnode
168  const node& nod = nodes[nodeI];
169  const treeBoundBox bb(nod.bb_.subBbox(octant));
170 
171  node subNode(divide(bb, contents, contentI));
172  subNode.parent_ = nodeI;
173  label sz = nodes.size();
174  nodes.push_back(subNode);
175  nodes[nodeI].subNodes_[octant] = nodePlusOctant(sz, octant);
176  }
177  }
178  }
179  }
180 }
181 
182 
183 template<class Type>
185 (
186  DynamicList<node>& nodes,
187  DynamicList<labelList>& contents,
188  const label compactLevel,
189  const label nodeI,
190  const label level,
191 
192  List<labelList>& compactedContents,
193  label& compactI
194 )
195 {
196  const node& nod = nodes[nodeI];
197 
198  label nNodes = 0;
199 
200  if (level < compactLevel)
201  {
202  for (direction octant = 0; octant < node::nChildren; ++octant)
203  {
204  labelBits index = nod.subNodes_[octant];
205 
206  if (isNode(index))
207  {
208  nNodes += compactContents
209  (
210  nodes,
211  contents,
212  compactLevel,
213  getNode(index),
214  level+1,
215  compactedContents,
216  compactI
217  );
218  }
219  }
220  }
221  else if (level == compactLevel)
222  {
223  // Compact all content on this level
224  for (direction octant = 0; octant < node::nChildren; ++octant)
225  {
226  labelBits index = nod.subNodes_[octant];
227 
228  if (isContent(index))
229  {
230  label contentI = getContent(index);
231 
232  compactedContents[compactI].transfer(contents[contentI]);
233 
234  // Subnode is at compactI. Adapt nodeI to point to it
235  nodes[nodeI].subNodes_[octant] =
236  contentPlusOctant(compactI, octant);
237 
238  compactI++;
239  }
240  else if (isNode(index))
241  {
242  nNodes++;
243  }
244  }
245  }
246  return nNodes;
247 }
248 
249 
250 template<class Type>
252 (
253  const label nodeI
254 ) const
255 {
256  // Pre-calculates wherever possible the volume status per node/subnode.
257  // Recurses to determine status of lowest level boxes. Level above is
258  // combination of octants below.
259 
260  const node& nod = nodes_[nodeI];
261 
262  volumeType myType = volumeType::UNKNOWN;
263 
264  for (direction octant = 0; octant < node::nChildren; ++octant)
265  {
266  volumeType subType;
267 
268  labelBits index = nod.subNodes_[octant];
269 
270  if (isNode(index))
271  {
272  // tree node. Recurse.
273  subType = calcVolumeType(getNode(index));
274  }
275  else if (isContent(index))
276  {
277  // Contents. Depending on position in box might be on either
278  // side.
279  subType = volumeType::MIXED;
280  }
281  else
282  {
283  // No data in this octant. Set type for octant acc. to the mid
284  // of its bounding box.
285  const treeBoundBox subBb = nod.bb_.subBbox(octant);
286 
287  subType = shapes_.getVolumeType(*this, subBb.centre());
288  }
289 
290  // Store octant type
291  nodeTypes_.set(labelBits::pack(nodeI, octant), subType);
292 
293  // Combine sub node types into type for treeNode. Result is 'mixed' if
294  // types differ among subnodes.
295  if (myType == volumeType::UNKNOWN)
296  {
297  myType = subType;
298  }
299  else if (subType != myType)
300  {
301  myType = volumeType::MIXED;
302  }
303  }
304  return myType;
305 }
306 
307 
308 template<class Type>
310 (
311  const label nodeI,
312  const point& sample
313 ) const
314 {
315  const node& nod = nodes_[nodeI];
316 
317  direction octant = nod.bb_.subOctant(sample);
318 
319  volumeType octantType = volumeType::type
320  (
321  nodeTypes_.get(labelBits::pack(nodeI, octant))
322  );
323 
324  if (octantType == volumeType::INSIDE)
325  {
326  return octantType;
327  }
328  else if (octantType == volumeType::OUTSIDE)
329  {
330  return octantType;
331  }
332  else if (octantType == volumeType::UNKNOWN)
333  {
334  // Can happen for e.g. non-manifold surfaces.
335  return octantType;
336  }
337  else if (octantType == volumeType::MIXED)
338  {
339  labelBits index = nod.subNodes_[octant];
340 
341  if (isNode(index))
342  {
343  // Recurse
344  volumeType subType = getVolumeType(getNode(index), sample);
345 
346  return subType;
347  }
348  else if (isContent(index))
349  {
350  // Content. Defer to shapes.
351  return volumeType(shapes_.getVolumeType(*this, sample));
352  }
353 
354  // Empty node. Cannot have 'mixed' as its type since not divided
355  // up and has no items inside it.
357  << "Sample:" << sample << " node:" << nodeI
358  << " with bb:" << nodes_[nodeI].bb_ << nl
359  << "Empty subnode has invalid volume type MIXED."
360  << abort(FatalError);
361 
362  return volumeType::UNKNOWN;
363  }
364 
366  << "Sample:" << sample << " at node:" << nodeI
367  << " octant:" << octant
368  << " with bb:" << nod.bb_.subBbox(octant) << nl
369  << "Node has invalid volume type " << octantType
370  << abort(FatalError);
372  return volumeType::UNKNOWN;
373 }
374 
375 
376 template<class Type>
378 (
379  const vector& outsideNormal,
380  const vector& vec
381 )
382 {
383  if ((outsideNormal&vec) >= 0)
384  {
385  return volumeType::OUTSIDE;
386  }
387  else
388  {
389  return volumeType::INSIDE;
390  }
391 }
392 
393 
394 template<class Type>
395 template<class FindNearestOp>
397 (
398  const label nodeI,
399  const point& sample,
400 
401  scalar& nearestDistSqr,
402  label& nearestShapeI,
403  point& nearestPoint,
404 
405  const FindNearestOp& fnOp
406 ) const
407 {
408  const node& nod = nodes_[nodeI];
409 
410  // Determine order to walk through octants
411  // Go into all suboctants (one containing sample first) and update nearest.
412 
413  for (const direction octant : nod.bb_.searchOrder(sample))
414  {
415  labelBits index = nod.subNodes_[octant];
416 
417  if (isNode(index))
418  {
419  label subNodeI = getNode(index);
420 
421  const treeBoundBox& subBb = nodes_[subNodeI].bb_;
422 
423  if (subBb.overlaps(sample, nearestDistSqr))
424  {
425  findNearest
426  (
427  subNodeI,
428  sample,
429 
430  nearestDistSqr,
431  nearestShapeI,
432  nearestPoint,
433 
434  fnOp
435  );
436  }
437  }
438  else if (isContent(index))
439  {
440  if (nod.bb_.subOverlaps(octant, sample, nearestDistSqr))
441  {
442  fnOp
443  (
444  contents_[getContent(index)],
445  sample,
446 
447  nearestDistSqr,
448  nearestShapeI,
449  nearestPoint
450  );
451  }
452  }
453  }
454 }
455 
456 
457 template<class Type>
458 template<class FindNearestOp>
460 (
461  const label nodeI,
462  const linePointRef& ln,
463 
464  treeBoundBox& tightest,
465  label& nearestShapeI,
466  point& linePoint,
467  point& nearestPoint,
468 
469  const FindNearestOp& fnOp
470 ) const
471 {
472  const node& nod = nodes_[nodeI];
473  const treeBoundBox& nodeBb = nod.bb_;
474 
475  // Determine order to walk through octants
476  // Go into all suboctants (one containing sample first) and update nearest.
477 
478  for (const direction octant : nod.bb_.searchOrder(ln.centre()))
479  {
480  labelBits index = nod.subNodes_[octant];
481 
482  if (isNode(index))
483  {
484  const treeBoundBox& subBb = nodes_[getNode(index)].bb_;
485 
486  if (subBb.overlaps(tightest))
487  {
488  findNearest
489  (
490  getNode(index),
491  ln,
492 
493  tightest,
494  nearestShapeI,
495  linePoint,
496  nearestPoint,
497 
498  fnOp
499  );
500  }
501  }
502  else if (isContent(index))
503  {
504  if (nodeBb.subOverlaps(octant, tightest))
505  {
506  fnOp
507  (
508  contents_[getContent(index)],
509  ln,
510 
511  tightest,
512  nearestShapeI,
513  linePoint,
514  nearestPoint
515  );
516  }
517  }
518  }
519 }
520 
521 
522 template<class Type>
524 (
525  const label parentNodeI,
526  const direction octant
527 ) const
528 {
529  // Get type of node at octant
530  const node& nod = nodes_[parentNodeI];
531  labelBits index = nod.subNodes_[octant];
532 
533  if (isNode(index))
534  {
535  // Use stored bb
536  return nodes_[getNode(index)].bb_;
537  }
538 
539  // Calculate subBb
540  return nod.bb_.subBbox(octant);
541 }
542 
543 
544 template<class Type>
546 (
547  const treeBoundBox& bb,
548  const point& pt,
549  const bool pushInside
550 )
551 {
552  // Get local length scale.
553  const vector perturbVec = perturbTol_*bb.span();
554 
555  point perturbedPt(pt);
556 
557  // Modify all components which are close to any face of the bb to be
558  // well inside/outside them.
559 
560  if (pushInside)
561  {
562  for (direction dir = 0; dir < vector::nComponents; dir++)
563  {
564  if (mag(pt[dir]-bb.min()[dir]) < mag(perturbVec[dir]))
565  {
566  // Close to 'left' side. Push well beyond left side.
567  scalar perturbDist = perturbVec[dir] + ROOTVSMALL;
568  perturbedPt[dir] = bb.min()[dir] + perturbDist;
569  }
570  else if (mag(pt[dir]-bb.max()[dir]) < mag(perturbVec[dir]))
571  {
572  // Close to 'right' side. Push well beyond right side.
573  scalar perturbDist = perturbVec[dir] + ROOTVSMALL;
574  perturbedPt[dir] = bb.max()[dir] - perturbDist;
575  }
576  }
577  }
578  else
579  {
580  for (direction dir = 0; dir < vector::nComponents; dir++)
581  {
582  if (mag(pt[dir]-bb.min()[dir]) < mag(perturbVec[dir]))
583  {
584  scalar perturbDist = perturbVec[dir] + ROOTVSMALL;
585  perturbedPt[dir] = bb.min()[dir] - perturbDist;
586  }
587  else if (mag(pt[dir]-bb.max()[dir]) < mag(perturbVec[dir]))
588  {
589  scalar perturbDist = perturbVec[dir] + ROOTVSMALL;
590  perturbedPt[dir] = bb.max()[dir] + perturbDist;
591  }
592  }
593  }
594 
595  if (debug)
596  {
597  if (pushInside != bb.contains(perturbedPt))
598  {
600  << "pushed point:" << pt
601  << " to:" << perturbedPt
602  << " wanted side:" << pushInside
603  << " obtained side:" << bb.contains(perturbedPt)
604  << " of bb:" << bb << nl;
605 
606  if (debug > 1)
607  {
609  }
610  }
611  }
612 
613  return perturbedPt;
614 }
615 
616 
617 template<class Type>
619 (
620  const treeBoundBox& bb,
621  const direction faceID,
622  const point& pt,
623  const bool pushInside
624 )
625 {
626  // Get local length scale.
627  const vector perturbVec = perturbTol_*bb.span();
628 
629  point perturbedPt(pt);
630 
631  // Modify all components which are close to any face of the bb to be
632  // well outside them.
633 
634  if (faceID == 0)
635  {
637  << abort(FatalError);
638  }
639 
640  {
641  constexpr direction dir(0); // vector::X
642 
643  if (faceID & treeBoundBox::LEFTBIT)
644  {
645  perturbedPt[dir] =
646  (
647  pushInside
648  ? (bb.min()[dir] + (perturbVec[dir] + ROOTVSMALL))
649  : (bb.min()[dir] - (perturbVec[dir] + ROOTVSMALL))
650  );
651  }
652  else if (faceID & treeBoundBox::RIGHTBIT)
653  {
654  perturbedPt[dir] =
655  (
656  pushInside
657  ? (bb.max()[dir] - (perturbVec[dir] + ROOTVSMALL))
658  : (bb.max()[dir] + (perturbVec[dir] + ROOTVSMALL))
659  );
660  }
661  }
662 
663  {
664  constexpr direction dir(1); // vector::Y
665 
666  if (faceID & treeBoundBox::BOTTOMBIT)
667  {
668  perturbedPt[dir] =
669  (
670  pushInside
671  ? (bb.min()[dir] + (perturbVec[dir] + ROOTVSMALL))
672  : (bb.min()[dir] - (perturbVec[dir] + ROOTVSMALL))
673  );
674  }
675  else if (faceID & treeBoundBox::TOPBIT)
676  {
677  perturbedPt[dir] =
678  (
679  pushInside
680  ? (bb.max()[dir] - (perturbVec[dir] + ROOTVSMALL))
681  : (bb.max()[dir] + (perturbVec[dir] + ROOTVSMALL))
682  );
683  }
684  }
685 
686  {
687  constexpr direction dir(2); // vector::Z
688 
689  if (faceID & treeBoundBox::BACKBIT)
690  {
691  perturbedPt[dir] =
692  (
693  pushInside
694  ? (bb.min()[dir] + (perturbVec[dir] + ROOTVSMALL))
695  : (bb.min()[dir] - (perturbVec[dir] + ROOTVSMALL))
696  );
697  }
698  else if (faceID & treeBoundBox::FRONTBIT)
699  {
700  perturbedPt[dir] =
701  (
702  pushInside
703  ? (bb.max()[dir] - (perturbVec[dir] + ROOTVSMALL))
704  : (bb.max()[dir] + (perturbVec[dir] + ROOTVSMALL))
705  );
706  }
707  }
708 
709  if (debug)
710  {
711  if (pushInside != bb.contains(perturbedPt))
712  {
714  << "pushed point:" << pt << " on face:" << faceString(faceID)
715  << " to:" << perturbedPt
716  << " wanted side:" << pushInside
717  << " obtained side:" << bb.contains(perturbedPt)
718  << " of bb:" << bb << nl;
719 
720  if (debug > 1)
721  {
723  }
724  }
725  }
726 
727  return perturbedPt;
728 }
729 
730 
731 template<class Type>
733 (
734  const treeBoundBox& bb,
735  const vector& dir, // end-start
736  const point& pt
737 )
738 {
739  if (debug)
740  {
741  if (bb.posBits(pt) != 0)
742  {
744  << " bb:" << bb << endl
745  << "does not contain point " << pt << nl;
746 
747  if (debug > 1)
748  {
750  }
751  }
752  }
753 
754 
755  // Handle two cases:
756  // - point exactly on multiple faces. Push away from all but one.
757  // - point on a single face. Push away from edges of face.
758 
759  direction ptFaceID = bb.faceBits(pt);
760 
761  direction nFaces = 0;
762  FixedList<direction, 3> faceIndices;
763 
764  if (ptFaceID & treeBoundBox::LEFTBIT)
765  {
766  faceIndices[nFaces++] = treeBoundBox::LEFT;
767  }
768  else if (ptFaceID & treeBoundBox::RIGHTBIT)
769  {
770  faceIndices[nFaces++] = treeBoundBox::RIGHT;
771  }
772 
773  if (ptFaceID & treeBoundBox::BOTTOMBIT)
774  {
775  faceIndices[nFaces++] = treeBoundBox::BOTTOM;
776  }
777  else if (ptFaceID & treeBoundBox::TOPBIT)
778  {
779  faceIndices[nFaces++] = treeBoundBox::TOP;
780  }
781 
782  if (ptFaceID & treeBoundBox::BACKBIT)
783  {
784  faceIndices[nFaces++] = treeBoundBox::BACK;
785  }
786  else if (ptFaceID & treeBoundBox::FRONTBIT)
787  {
788  faceIndices[nFaces++] = treeBoundBox::FRONT;
789  }
790 
791 
792  // Determine the face we want to keep the point on
793 
794  direction keepFaceID;
795 
796  if (nFaces == 0)
797  {
798  // Return original point
799  return pt;
800  }
801  else if (nFaces == 1)
802  {
803  // Point is on a single face
804  keepFaceID = faceIndices[0];
805  }
806  else
807  {
808  // Determine best face out of faceIndices[0 .. nFaces-1].
809  // The best face is the one most perpendicular to the ray direction.
810 
811  keepFaceID = faceIndices[0];
812  scalar maxInproduct = mag(treeBoundBox::faceNormals[keepFaceID] & dir);
813 
814  for (direction i = 1; i < nFaces; i++)
815  {
816  direction face = faceIndices[i];
817  scalar s = mag(treeBoundBox::faceNormals[face] & dir);
818  if (s > maxInproduct)
819  {
820  maxInproduct = s;
821  keepFaceID = face;
822  }
823  }
824  }
825 
826 
827  // 1. Push point into bb, away from all corners
828 
829  point facePoint(pushPoint(bb, pt, true));
830  direction faceID = 0;
831 
832  // 2. Snap it back onto the preferred face
833 
834  if (keepFaceID == treeBoundBox::LEFT)
835  {
836  facePoint.x() = bb.min().x();
837  faceID = treeBoundBox::LEFTBIT;
838  }
839  else if (keepFaceID == treeBoundBox::RIGHT)
840  {
841  facePoint.x() = bb.max().x();
842  faceID = treeBoundBox::RIGHTBIT;
843  }
844  else if (keepFaceID == treeBoundBox::BOTTOM)
845  {
846  facePoint.y() = bb.min().y();
847  faceID = treeBoundBox::BOTTOMBIT;
848  }
849  else if (keepFaceID == treeBoundBox::TOP)
850  {
851  facePoint.y() = bb.max().y();
852  faceID = treeBoundBox::TOPBIT;
853  }
854  else if (keepFaceID == treeBoundBox::BACK)
855  {
856  facePoint.z() = bb.min().z();
857  faceID = treeBoundBox::BACKBIT;
858  }
859  else if (keepFaceID == treeBoundBox::FRONT)
860  {
861  facePoint.z() = bb.max().z();
862  faceID = treeBoundBox::FRONTBIT;
863  }
864 
865 
866  if (debug)
867  {
868  if (faceID != bb.faceBits(facePoint))
869  {
871  << "Pushed point from " << pt
872  << " on face:" << ptFaceID << " of bb:" << bb << nl
873  << "onto " << facePoint
874  << " on face:" << faceID
875  << " which is not consistent with geometric face "
876  << bb.faceBits(facePoint) << nl;
877 
878  if (debug > 1)
879  {
881  }
882  }
883  if (bb.posBits(facePoint) != 0)
884  {
886  << " bb:" << bb << nl
887  << "does not contain perturbed point "
888  << facePoint << nl;
889 
890  if (debug > 1)
891  {
893  }
894  }
895  }
896 
897 
898  return facePoint;
899 }
900 
901 
902 template<class Type>
904 (
905  const label nodeI,
906  const direction octant,
907 
908  label& parentNodeI,
909  label& parentOctant
910 ) const
911 {
912  parentNodeI = nodes_[nodeI].parent_;
913 
914  if (parentNodeI == -1)
915  {
916  // Reached edge of tree
917  return false;
918  }
919 
920  const node& parentNode = nodes_[parentNodeI];
921 
922  // Find octant nodeI is in.
923  parentOctant = 255;
924 
925  for (direction i = 0; i < node::nChildren; ++i)
926  {
927  labelBits index = parentNode.subNodes_[i];
928 
929  if (isNode(index) && getNode(index) == nodeI)
930  {
931  parentOctant = i;
932  break;
933  }
934  }
935 
936  if (parentOctant == 255)
937  {
939  << "Problem: no parent found for octant:" << octant
940  << " node:" << nodeI
941  << abort(FatalError);
942  }
943 
944  return true;
945 }
946 
947 
948 template<class Type>
950 (
951  const point& facePoint,
952  const direction faceID, // face(s) that facePoint is on
953  label& nodeI,
954  direction& octant
955 ) const
956 {
957  // Gets current position as node and octant in this node and walks in the
958  // direction given by the facePointBits (combination of
959  // treeBoundBox::LEFTBIT, TOPBIT etc.) Returns false if edge of tree hit.
960 
961  label oldNodeI = nodeI;
962  direction oldOctant = octant;
963 
964  // Find out how to test for octant. Say if we want to go left we need
965  // to traverse up the tree until we hit a node where our octant is
966  // on the right.
967 
968  // Coordinate direction to test
969  const direction X = treeBoundBox::RIGHTHALF;
970  const direction Y = treeBoundBox::TOPHALF;
971  const direction Z = treeBoundBox::FRONTHALF;
972 
973  direction octantMask = 0;
974  direction wantedValue = 0;
975 
976  if ((faceID & treeBoundBox::LEFTBIT) != 0)
977  {
978  // We want to go left so check if in right octant (i.e. x-bit is set)
979  octantMask |= X;
980  wantedValue |= X;
981  }
982  else if ((faceID & treeBoundBox::RIGHTBIT) != 0)
983  {
984  octantMask |= X; // wantedValue already 0
985  }
986 
987  if ((faceID & treeBoundBox::BOTTOMBIT) != 0)
988  {
989  // Want to go down so check for y-bit set.
990  octantMask |= Y;
991  wantedValue |= Y;
992  }
993  else if ((faceID & treeBoundBox::TOPBIT) != 0)
994  {
995  // Want to go up so check for y-bit not set.
996  octantMask |= Y;
997  }
998 
999  if ((faceID & treeBoundBox::BACKBIT) != 0)
1000  {
1001  octantMask |= Z;
1002  wantedValue |= Z;
1003  }
1004  else if ((faceID & treeBoundBox::FRONTBIT) != 0)
1005  {
1006  octantMask |= Z;
1007  }
1008 
1009  // So now we have the coordinate directions in the octant we need to check
1010  // and the resulting value.
1011 
1012  /*
1013  // +---+---+
1014  // | | |
1015  // | | |
1016  // | | |
1017  // +---+-+-+
1018  // | | | |
1019  // | a+-+-+
1020  // | |\| |
1021  // +---+-+-+
1022  // \
1023  //
1024  // e.g. ray is at (a) in octant 0(or 4) with faceIDs : LEFTBIT+TOPBIT.
1025  // If we would be in octant 1(or 5) we could go to the correct octant
1026  // in the same node by just flipping the x and y bits (exoring).
1027  // But if we are not in octant 1/5 we have to go up until we are.
1028  // In general for leftbit+topbit:
1029  // - we have to check for x and y : octantMask = 011
1030  // - the checked bits have to be : wantedValue = ?01
1031  */
1032 
1033  //Pout<< "For point " << facePoint << endl;
1034 
1035  // Go up until we have chance to cross to the wanted direction
1036  while (wantedValue != (octant & octantMask))
1037  {
1038  // Go up to the parent.
1039 
1040  // Remove the directions that are not on the boundary of the parent.
1041  // See diagram above
1042  if (wantedValue & X) // && octantMask&X
1043  {
1044  // Looking for right octant.
1045  if (octant & X)
1046  {
1047  // My octant is one of the left ones so punch out the x check
1048  octantMask &= ~X;
1049  wantedValue &= ~X;
1050  }
1051  }
1052  else
1053  {
1054  if (!(octant & X))
1055  {
1056  // My octant is right but I want to go left.
1057  octantMask &= ~X;
1058  wantedValue &= ~X;
1059  }
1060  }
1061 
1062  if (wantedValue & Y)
1063  {
1064  if (octant & Y)
1065  {
1066  octantMask &= ~Y;
1067  wantedValue &= ~Y;
1068  }
1069  }
1070  else
1071  {
1072  if (!(octant & Y))
1073  {
1074  octantMask &= ~Y;
1075  wantedValue &= ~Y;
1076  }
1077  }
1078 
1079  if (wantedValue & Z)
1080  {
1081  if (octant & Z)
1082  {
1083  octantMask &= ~Z;
1084  wantedValue &= ~Z;
1085  }
1086  }
1087  else
1088  {
1089  if (!(octant & Z))
1090  {
1091  octantMask &= ~Z;
1092  wantedValue &= ~Z;
1093  }
1094  }
1095 
1096 
1097  label parentNodeI;
1098  label parentOctant;
1099  walkToParent(nodeI, octant, parentNodeI, parentOctant);
1100 
1101  if (parentNodeI == -1)
1102  {
1103  // Reached edge of tree
1104  return false;
1105  }
1106 
1107  //Pout<< " walked from node:" << nodeI << " octant:" << octant
1108  // << " bb:" << nodes_[nodeI].bb_.subBbox(octant) << endl
1109  // << " to:" << parentNodeI << " octant:" << parentOctant
1110  // << " bb:" << nodes_[parentNodeI].bb_.subBbox(parentOctant)
1111  // << endl;
1112  //
1113  //Pout<< " octantMask:" << octantMask
1114  // << " wantedValue:" << wantedValue << endl;
1115 
1116  nodeI = parentNodeI;
1117  octant = parentOctant;
1118  }
1119 
1120  // So now we hit the correct parent node. Switch to the correct octant.
1121  // Is done by jumping to the 'other half' so if e.g. in x direction in
1122  // right half we now jump to the left half.
1123  octant ^= octantMask;
1124 
1125  //Pout<< " to node:" << nodeI << " octant:" << octant
1126  // << " subBb:" <<subBbox(nodeI, octant) << endl;
1127 
1128 
1129  if (debug)
1130  {
1131  const treeBoundBox subBb(subBbox(nodeI, octant));
1132 
1133  if (!subBb.contains(facePoint))
1134  {
1136  << "When searching for " << facePoint
1137  << " ended up in node:" << nodeI
1138  << " octant:" << octant
1139  << " with bb:" << subBb << nl;
1140 
1141  if (debug > 1)
1142  {
1144  }
1145  }
1146  }
1147 
1148 
1149  // See if we need to travel down. Note that we already go into the
1150  // the first level ourselves (instead of having findNode decide)
1151  labelBits index = nodes_[nodeI].subNodes_[octant];
1152 
1153  if (isNode(index))
1154  {
1155  labelBits node = findNode(getNode(index), facePoint);
1156 
1157  nodeI = getNode(node);
1158  octant = getOctant(node);
1159  }
1160 
1161 
1162  if (debug)
1163  {
1164  const treeBoundBox subBb(subBbox(nodeI, octant));
1165 
1166  if (nodeI == oldNodeI && octant == oldOctant)
1167  {
1169  << "Did not go to neighbour when searching for " << facePoint
1170  << nl
1171  << " starting from face:" << faceString(faceID)
1172  << " node:" << nodeI
1173  << " octant:" << octant
1174  << " bb:" << subBb << nl;
1175 
1176  if (debug > 1)
1177  {
1179  }
1180  }
1181 
1182  if (!subBb.contains(facePoint))
1183  {
1185  << "When searching for " << facePoint
1186  << " ended up in node:" << nodeI
1187  << " octant:" << octant
1188  << " bb:" << subBb << nl;
1189 
1190  if (debug > 1)
1191  {
1193  }
1194  }
1195  }
1196 
1197 
1198  return true;
1199 }
1200 
1201 
1202 template<class Type>
1204 (
1205  const direction faceID
1206 )
1207 {
1208  word desc;
1209 
1210  if (faceID == 0)
1211  {
1212  desc = "noFace";
1213  }
1214  if (faceID & treeBoundBox::LEFTBIT)
1215  {
1216  if (!desc.empty()) desc += "+";
1217  desc += "left";
1218  }
1219  if (faceID & treeBoundBox::RIGHTBIT)
1220  {
1221  if (!desc.empty()) desc += "+";
1222  desc += "right";
1223  }
1224  if (faceID & treeBoundBox::BOTTOMBIT)
1225  {
1226  if (!desc.empty()) desc += "+";
1227  desc += "bottom";
1228  }
1229  if (faceID & treeBoundBox::TOPBIT)
1230  {
1231  if (!desc.empty()) desc += "+";
1232  desc += "top";
1233  }
1234  if (faceID & treeBoundBox::BACKBIT)
1235  {
1236  if (!desc.empty()) desc += "+";
1237  desc += "back";
1238  }
1239  if (faceID & treeBoundBox::FRONTBIT)
1240  {
1241  if (!desc.empty()) desc += "+";
1242  desc += "front";
1243  }
1244  return desc;
1245 }
1246 
1247 
1248 template<class Type>
1249 template<class FindIntersectOp>
1251 (
1252  const bool findAny,
1253  const point& treeStart,
1254  const vector& treeVec,
1255 
1256  const point& start,
1257  const point& end,
1258  const label nodeI,
1259  const direction octant,
1260 
1261  pointIndexHit& hitInfo,
1262  direction& hitBits,
1263 
1264  const FindIntersectOp& fiOp
1265 ) const
1266 {
1267  // Traverse a node. If intersects a triangle return first intersection
1268  // point:
1269  // hitInfo.index = index of shape
1270  // hitInfo.point = point on shape
1271  // Else return a miss and the bounding box face hit:
1272  // hitInfo.point = coordinate of intersection of ray with bounding box
1273  // hitBits = posbits of point on bounding box
1274 
1275  if (debug)
1276  {
1277  const treeBoundBox octantBb(subBbox(nodeI, octant));
1278 
1279  if (octantBb.posBits(start) != 0)
1280  {
1282  << "Node:" << nodeI << " octant:" << octant
1283  << " bb:" << octantBb << nl
1284  << "does not contain point " << start << nl;
1285 
1286  if (debug > 1)
1287  {
1289  }
1290  }
1291  }
1292 
1293  const node& nod = nodes_[nodeI];
1294 
1295  labelBits index = nod.subNodes_[octant];
1296 
1297  if (isContent(index))
1298  {
1299  const labelList& indices = contents_[getContent(index)];
1300 
1301  if (indices.size())
1302  {
1303  if (findAny)
1304  {
1305  // Find any intersection
1306 
1307  forAll(indices, elemI)
1308  {
1309  label shapeI = indices[elemI];
1310 
1311  point pt;
1312  bool hit = fiOp(shapeI, start, end, pt);
1313 
1314  // Note that intersection of shape might actually be
1315  // in a neighbouring box. For findAny this is not important.
1316  if (hit)
1317  {
1318  // Hit so pt is nearer than nearestPoint.
1319  // Update hit info
1320  hitInfo.hitPoint(pt);
1321  hitInfo.setIndex(shapeI);
1322  return;
1323  }
1324  }
1325  }
1326  else
1327  {
1328  // Find nearest intersection
1329 
1330  const treeBoundBox octantBb(subBbox(nodeI, octant));
1331 
1332  point nearestPoint(end);
1333 
1334  forAll(indices, elemI)
1335  {
1336  label shapeI = indices[elemI];
1337 
1338  point pt;
1339  bool hit = fiOp(shapeI, start, nearestPoint, pt);
1340 
1341  // Note that intersection of shape might actually be
1342  // in a neighbouring box. Since we need to maintain strict
1343  // (findAny=false) ordering skip such an intersection. It
1344  // will be found when we are doing the next box.
1345 
1346  if (hit && octantBb.contains(pt))
1347  {
1348  // Hit so pt is nearer than nearestPoint.
1349  nearestPoint = pt;
1350  // Update hit info
1351  hitInfo.hitPoint(pt);
1352  hitInfo.setIndex(shapeI);
1353  }
1354  }
1355 
1356  if (hitInfo.hit())
1357  {
1358  // Found intersected shape.
1359  return;
1360  }
1361  }
1362  }
1363  }
1364 
1365  // Nothing intersected in this node
1366  // Traverse to end of node. Do by ray tracing back from end and finding
1367  // intersection with bounding box of node.
1368  // start point is now considered to be inside bounding box.
1369 
1370  const treeBoundBox octantBb(subBbox(nodeI, octant));
1371 
1372  point pt;
1373  bool intersected = octantBb.intersects
1374  (
1375  end, //treeStart,
1376  (start-end), //treeVec,
1377 
1378  end,
1379  start,
1380 
1381  pt,
1382  hitBits
1383  );
1384 
1385 
1386  if (intersected)
1387  {
1388  // Return miss. Set misspoint to face.
1389  hitInfo.setPoint(pt);
1390  }
1391  else
1392  {
1393  // Rare case: did not find intersection of ray with octantBb. Can
1394  // happen if end is on face/edge of octantBb. Do like
1395  // lagrangian and re-track with perturbed vector (slightly pushed out
1396  // of bounding box)
1397 
1398  point perturbedEnd(pushPoint(octantBb, end, false));
1399 
1400  traverseNode
1401  (
1402  findAny,
1403  treeStart,
1404  treeVec,
1405  start,
1406  perturbedEnd,
1407  nodeI,
1408  octant,
1409 
1410  hitInfo,
1411  hitBits,
1412 
1413  fiOp
1414  );
1415  }
1416 }
1417 
1418 
1419 template<class Type>
1420 template<class FindIntersectOp>
1422 (
1423  const bool findAny,
1424  const point& treeStart,
1425  const point& treeEnd,
1426  const label startNodeI,
1427  const direction startOctant,
1428  const FindIntersectOp& fiOp,
1429  const bool verbose
1430 ) const
1431 {
1432  const vector treeVec(treeEnd - treeStart);
1433 
1434  // Current node as parent+octant
1435  label nodeI = startNodeI;
1436  direction octant = startOctant;
1437 
1438  if (verbose)
1439  {
1440  Pout<< "findLine : treeStart:" << treeStart
1441  << " treeEnd:" << treeEnd << endl
1442  << "node:" << nodeI
1443  << " octant:" << octant
1444  << " bb:" << subBbox(nodeI, octant) << endl;
1445  }
1446 
1447  // Current position. Initialize to miss
1448  pointIndexHit hitInfo(false, treeStart, -1);
1449 
1450  //while (true)
1451  label i = 0;
1452  for (; i < 100000; i++)
1453  {
1454  // Ray-trace to end of current node. Updates point (either on triangle
1455  // in case of hit or on node bounding box in case of miss)
1456 
1457  const treeBoundBox octantBb(subBbox(nodeI, octant));
1458 
1459  // Make sure point is away from any edges/corners
1460  point startPoint
1461  (
1462  pushPointIntoFace
1463  (
1464  octantBb,
1465  treeVec,
1466  hitInfo.point()
1467  )
1468  );
1469 
1470  if (verbose)
1471  {
1472  Pout<< "iter:" << i
1473  << " at current:" << hitInfo.point()
1474  << " (perturbed:" << startPoint << ")" << endl
1475  << " node:" << nodeI
1476  << " octant:" << octant
1477  << " bb:" << subBbox(nodeI, octant) << endl;
1478  }
1479 
1480 
1481 
1482 
1483  // Faces of current bounding box current point is on
1484  direction hitFaceID = 0;
1485 
1486  traverseNode
1487  (
1488  findAny,
1489  treeStart,
1490  treeVec,
1491 
1492  startPoint, // Note: pass in copy since hitInfo
1493  // also used in return value.
1494  treeEnd, // pass in overall end so is nicely outside bb
1495  nodeI,
1496  octant,
1497 
1498  hitInfo,
1499  hitFaceID,
1500 
1501  fiOp
1502  );
1503 
1504  // Did we hit a triangle?
1505  if (hitInfo.hit())
1506  {
1507  break;
1508  }
1509 
1510  if (hitFaceID == 0 || hitInfo.point() == treeEnd)
1511  {
1512  // endpoint inside the tree. Return miss.
1513  break;
1514  }
1515 
1516  // Create a point on other side of face.
1517  point perturbedPoint
1518  (
1519  pushPoint
1520  (
1521  octantBb,
1522  hitFaceID,
1523  hitInfo.point(),
1524  false // push outside of octantBb
1525  )
1526  );
1527 
1528  if (verbose)
1529  {
1530  Pout<< " iter:" << i
1531  << " hit face:" << faceString(hitFaceID)
1532  << " at:" << hitInfo.point() << nl
1533  << " node:" << nodeI
1534  << " octant:" << octant
1535  << " bb:" << subBbox(nodeI, octant) << nl
1536  << " walking to neighbour containing:" << perturbedPoint
1537  << endl;
1538  }
1539 
1540 
1541  // Nothing hit so we are on face of bounding box (given as node+octant+
1542  // position bits). Traverse to neighbouring node. Use slightly perturbed
1543  // point.
1544 
1545  bool ok = walkToNeighbour
1546  (
1547  perturbedPoint,
1548  hitFaceID, // face(s) that hitInfo is on
1549 
1550  nodeI,
1551  octant
1552  );
1553 
1554  if (!ok)
1555  {
1556  // Hit the edge of the tree. Return miss.
1557  break;
1558  }
1559 
1560  if (verbose)
1561  {
1562  const treeBoundBox octantBb(subBbox(nodeI, octant));
1563  Pout<< " walked for point:" << hitInfo.point() << endl
1564  << " to neighbour node:" << nodeI
1565  << " octant:" << octant
1566  << " face:" << faceString(octantBb.faceBits(hitInfo.point()))
1567  << " of octantBb:" << octantBb << endl
1568  << endl;
1569  }
1570  }
1571 
1572  if (i == 100000)
1573  {
1574  // Probably in loop.
1575  if (!verbose)
1576  {
1577  // Redo intersection but now with verbose flag switched on.
1578  return findLine
1579  (
1580  findAny,
1581  treeStart,
1582  treeEnd,
1583  startNodeI,
1584  startOctant,
1585  fiOp,
1586  true //verbose,
1587  );
1588  }
1589  if (debug)
1590  {
1592  << "Got stuck in loop raytracing from:" << treeStart
1593  << " to:" << treeEnd << endl
1594  << "inside top box:" << subBbox(startNodeI, startOctant)
1595  << abort(FatalError);
1596  }
1597  else
1598  {
1600  << "Got stuck in loop raytracing from:" << treeStart
1601  << " to:" << treeEnd << endl
1602  << "inside top box:" << subBbox(startNodeI, startOctant)
1603  << endl;
1604  }
1605  }
1606 
1607  return hitInfo;
1608 }
1609 
1610 
1611 template<class Type>
1612 template<class FindIntersectOp>
1614 (
1615  const bool findAny,
1616  const point& start,
1617  const point& end,
1618  const FindIntersectOp& fiOp
1619 ) const
1620 {
1621  pointIndexHit hitInfo;
1622 
1623  if (nodes_.size())
1624  {
1625  const treeBoundBox& treeBb = nodes_[0].bb_;
1626 
1627  // No effort is made to deal with points which are on edge of tree
1628  // bounding box for now.
1629 
1630  direction startBit = treeBb.posBits(start);
1631  direction endBit = treeBb.posBits(end);
1632 
1633  if ((startBit & endBit) != 0)
1634  {
1635  // Both start and end outside domain and in same block.
1636  return pointIndexHit(false, Zero, -1);
1637  }
1638 
1639 
1640  // Trim segment to treeBb
1641 
1642  point trackStart(start);
1643  point trackEnd(end);
1644 
1645  if (startBit != 0)
1646  {
1647  // Track start to inside domain.
1648  if (!treeBb.intersects(start, end, trackStart))
1649  {
1650  return pointIndexHit(false, Zero, -1);
1651  }
1652  }
1653 
1654  if (endBit != 0)
1655  {
1656  // Track end to inside domain.
1657  if (!treeBb.intersects(end, trackStart, trackEnd))
1658  {
1659  return pointIndexHit(false, Zero, -1);
1660  }
1661  }
1662 
1663 
1664  // Find lowest level tree node that start is in.
1665  labelBits index = findNode(0, trackStart);
1666 
1667  label parentNodeI = getNode(index);
1668  direction octant = getOctant(index);
1669 
1670  hitInfo = findLine
1671  (
1672  findAny,
1673  trackStart,
1674  trackEnd,
1675  parentNodeI,
1676  octant,
1677  fiOp
1678  );
1679  }
1680 
1681  return hitInfo;
1682 }
1683 
1684 
1685 template<class Type>
1687 (
1688  const label nodeI,
1689  const treeBoundBox& searchBox,
1690  labelHashSet* elements
1691 ) const
1692 {
1693  const node& nod = nodes_[nodeI];
1694  const treeBoundBox& nodeBb = nod.bb_;
1695 
1696  bool foundAny = false;
1697 
1698  for (direction octant = 0; octant < node::nChildren; ++octant)
1699  {
1700  labelBits index = nod.subNodes_[octant];
1701 
1702  if (isNode(index))
1703  {
1704  const treeBoundBox& subBb = nodes_[getNode(index)].bb_;
1705 
1706  if (subBb.overlaps(searchBox))
1707  {
1708  if (findBox(getNode(index), searchBox, elements))
1709  {
1710  // Early exit if not storing results
1711  if (!elements) return true;
1712 
1713  foundAny = true;
1714  }
1715  }
1716  }
1717  else if (isContent(index))
1718  {
1719  if (nodeBb.subOverlaps(octant, searchBox))
1720  {
1721  const labelList& indices = contents_[getContent(index)];
1722 
1723  for (const label index : indices)
1724  {
1725  if (shapes_.overlaps(index, searchBox))
1726  {
1727  // Early exit if not storing results
1728  if (!elements) return true;
1729 
1730  foundAny = true;
1731  elements->insert(index);
1732  }
1733  }
1734  }
1735  }
1736  }
1737 
1738  return foundAny;
1739 }
1740 
1741 
1742 template<class Type>
1744 (
1745  const label nodeI,
1746  const point& centre,
1747  const scalar radiusSqr,
1748  labelHashSet* elements
1749 ) const
1750 {
1751  const node& nod = nodes_[nodeI];
1752  const treeBoundBox& nodeBb = nod.bb_;
1753 
1754  bool foundAny = false;
1755 
1756  for (direction octant = 0; octant < node::nChildren; ++octant)
1757  {
1758  labelBits index = nod.subNodes_[octant];
1759 
1760  if (isNode(index))
1761  {
1762  const treeBoundBox& subBb = nodes_[getNode(index)].bb_;
1763 
1764  if (subBb.overlaps(centre, radiusSqr))
1765  {
1766  if (findSphere(getNode(index), centre, radiusSqr, elements))
1767  {
1768  // Early exit if not storing results
1769  if (!elements) return true;
1770 
1771  foundAny = true;
1772  }
1773  }
1774  }
1775  else if (isContent(index))
1776  {
1777  if (nodeBb.subOverlaps(octant, centre, radiusSqr))
1778  {
1779  const labelList& indices = contents_[getContent(index)];
1780 
1781  for (const label index : indices)
1782  {
1783  if (shapes_.overlaps(index, centre, radiusSqr))
1784  {
1785  // Early exit if not storing results
1786  if (!elements) return true;
1787 
1788  foundAny = true;
1789  elements->insert(index);
1790  }
1791  }
1792  }
1793  }
1794  }
1795 
1796  return foundAny;
1797 }
1798 
1799 
1800 template<class Type>
1801 template<class CompareOp>
1803 (
1804  const scalar nearDist,
1805  const bool okOrder,
1806  const indexedOctree<Type>& tree1,
1807  const labelBits index1,
1808  const treeBoundBox& bb1,
1809  const indexedOctree<Type>& tree2,
1810  const labelBits index2,
1811  const treeBoundBox& bb2,
1812  CompareOp& cop
1813 )
1814 {
1815  const vector nearSpan(nearDist, nearDist, nearDist);
1816 
1817  if (tree1.isNode(index1))
1818  {
1819  const node& nod1 = tree1.nodes()[tree1.getNode(index1)];
1820  const treeBoundBox searchBox
1821  (
1822  bb1.min()-nearSpan,
1823  bb1.max()+nearSpan
1824  );
1825 
1826  if (tree2.isNode(index2))
1827  {
1828  if (bb2.overlaps(searchBox))
1829  {
1830  const node& nod2 = tree2.nodes()[tree2.getNode(index2)];
1831 
1832  for (direction i2 = 0; i2 < node::nChildren; ++i2)
1833  {
1834  labelBits subIndex2 = nod2.subNodes_[i2];
1835  const treeBoundBox subBb2
1836  (
1837  tree2.isNode(subIndex2)
1838  ? tree2.nodes()[tree2.getNode(subIndex2)].bb_
1839  : bb2.subBbox(i2)
1840  );
1841 
1842  findNear
1843  (
1844  nearDist,
1845  !okOrder,
1846  tree2,
1847  subIndex2,
1848  subBb2,
1849  tree1,
1850  index1,
1851  bb1,
1852  cop
1853  );
1854  }
1855  }
1856  }
1857  else if (tree2.isContent(index2))
1858  {
1859  // index2 is leaf, index1 not yet.
1860  for (direction i1 = 0; i1 < node::nChildren; ++i1)
1861  {
1862  labelBits subIndex1 = nod1.subNodes_[i1];
1863  const treeBoundBox subBb1
1864  (
1865  tree1.isNode(subIndex1)
1866  ? tree1.nodes()[tree1.getNode(subIndex1)].bb_
1867  : bb1.subBbox(i1)
1868  );
1869 
1870  findNear
1871  (
1872  nearDist,
1873  !okOrder,
1874  tree2,
1875  index2,
1876  bb2,
1877  tree1,
1878  subIndex1,
1879  subBb1,
1880  cop
1881  );
1882  }
1883  }
1884  }
1885  else if (tree1.isContent(index1))
1886  {
1887  const treeBoundBox searchBox
1888  (
1889  bb1.min()-nearSpan,
1890  bb1.max()+nearSpan
1891  );
1892 
1893  if (tree2.isNode(index2))
1894  {
1895  const node& nod2 =
1896  tree2.nodes()[tree2.getNode(index2)];
1897 
1898  if (bb2.overlaps(searchBox))
1899  {
1900  for (direction i2 = 0; i2 < node::nChildren; ++i2)
1901  {
1902  labelBits subIndex2 = nod2.subNodes_[i2];
1903  const treeBoundBox subBb2
1904  (
1905  tree2.isNode(subIndex2)
1906  ? tree2.nodes()[tree2.getNode(subIndex2)].bb_
1907  : bb2.subBbox(i2)
1908  );
1909 
1910  findNear
1911  (
1912  nearDist,
1913  !okOrder,
1914  tree2,
1915  subIndex2,
1916  subBb2,
1917  tree1,
1918  index1,
1919  bb1,
1920  cop
1921  );
1922  }
1923  }
1924  }
1925  else if (tree2.isContent(index2))
1926  {
1927  // Both are leaves. Check n^2.
1928 
1929  const labelList& indices1 =
1930  tree1.contents()[tree1.getContent(index1)];
1931  const labelList& indices2 =
1932  tree2.contents()[tree2.getContent(index2)];
1933 
1934  forAll(indices1, i)
1935  {
1936  label shape1 = indices1[i];
1937 
1938  forAll(indices2, j)
1939  {
1940  label shape2 = indices2[j];
1941 
1942  if ((&tree1 != &tree2) || (shape1 != shape2))
1943  {
1944  if (okOrder)
1945  {
1946  cop
1947  (
1948  nearDist,
1949  tree1.shapes(),
1950  shape1,
1951  tree2.shapes(),
1952  shape2
1953  );
1954  }
1955  else
1956  {
1957  cop
1958  (
1959  nearDist,
1960  tree2.shapes(),
1961  shape2,
1962  tree1.shapes(),
1963  shape1
1964  );
1965  }
1966  }
1967  }
1968  }
1969  }
1970  }
1971 }
1972 
1973 
1974 template<class Type>
1975 Foam::label Foam::indexedOctree<Type>::countLeafs(const label nodeI) const
1976 {
1977  label total = 0;
1978 
1979  const node& nod = nodes_[nodeI];
1980 
1981  for (direction octant = 0; octant < node::nChildren; ++octant)
1982  {
1983  labelBits index = nod.subNodes_[octant];
1984 
1985  if (isNode(index))
1986  {
1987  total += countLeafs(getNode(index));
1988  }
1989  else if (isContent(index))
1990  {
1991  ++total;
1992  }
1993  }
1994 
1995  return total;
1996 }
1997 
1998 
1999 template<class Type>
2001 (
2002  const labelBits index
2003 ) const
2004 {
2005  label nElems = 0;
2006 
2007  if (isNode(index))
2008  {
2009  // tree node.
2010  label nodeI = getNode(index);
2011 
2012  const node& nod = nodes_[nodeI];
2013 
2014  for (direction octant = 0; octant < node::nChildren; ++octant)
2015  {
2016  nElems += countElements(nod.subNodes_[octant]);
2017  }
2018  }
2019  else if (isContent(index))
2020  {
2021  nElems += contents_[getContent(index)].size();
2022  }
2023  else
2024  {
2025  // empty node
2026  }
2027 
2028  return nElems;
2029 }
2030 
2031 
2032 template<class Type>
2034 (
2035  const label nodeI,
2036  Ostream& os,
2037  label& vertIndex,
2038  const bool leavesOnly,
2039  const bool writeLinesOnly
2040 ) const
2041 {
2042  const node& nod = nodes_[nodeI];
2043  const treeBoundBox& bb = nod.bb_;
2044 
2045  for (direction octant = 0; octant < node::nChildren; ++octant)
2046  {
2047  const treeBoundBox subBb(bb.subBbox(octant));
2048 
2049  labelBits index = nod.subNodes_[octant];
2050 
2051  if (isNode(index))
2052  {
2053  label subNodeI = getNode(index);
2054 
2055  writeOBJ(subNodeI, os, vertIndex, leavesOnly, writeLinesOnly);
2056  }
2057  else if (isContent(index))
2058  {
2059  indexedOctreeBase::writeOBJ(os, subBb, vertIndex, writeLinesOnly);
2060  }
2061  else if (isEmpty(index) && !leavesOnly)
2062  {
2063  indexedOctreeBase::writeOBJ(os, subBb, vertIndex, writeLinesOnly);
2064  }
2065  }
2066 }
2067 
2068 
2069 template<class Type>
2071 (
2072  const label nodeI,
2073  const direction octant
2074 ) const
2075 {
2076  labelBits index = nodes_[nodeI].subNodes_[octant];
2077 
2078  treeBoundBox subBb;
2079 
2080  if (isNode(index))
2081  {
2082  subBb = nodes_[getNode(index)].bb_;
2083  }
2084  else if (isContent(index) || isEmpty(index))
2085  {
2086  subBb = nodes_[nodeI].bb_.subBbox(octant);
2087  }
2088 
2089  OFstream os
2090  (
2091  "node" + Foam::name(nodeI) + "_octant" + Foam::name(octant) + ".obj"
2092  );
2093 
2094  Pout<< "dumpContentNode : writing node:" << nodeI << " octant:" << octant
2095  << " to " << os.name() << endl;
2096 
2097  bool writeLinesOnly(false);
2098  label vertIndex(0);
2099  indexedOctreeBase::writeOBJ(os, subBb, vertIndex, writeLinesOnly);
2100 }
2101 
2102 
2103 template<class Type>
2105 {
2106  if (!nodes_.empty())
2107  {
2108  label vertIndex(0);
2109  // leavesOnly=true, writeLinesOnly=false
2110  writeOBJ(0, os, vertIndex, true, false);
2111  }
2112 }
2113 
2114 
2115 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
2116 
2117 template<class Type>
2118 Foam::indexedOctree<Type>::indexedOctree(const Type& shapes)
2119 :
2120  shapes_(shapes),
2121  nodes_(),
2122  contents_(),
2123  nodeTypes_()
2124 {}
2125 
2126 
2127 template<class Type>
2129 (
2130  const Type& shapes,
2131  const List<node>& nodes,
2132  const List<labelList>& contents
2133 )
2134 :
2135  shapes_(shapes),
2136  nodes_(nodes),
2137  contents_(contents),
2138  nodeTypes_()
2139 {}
2140 
2141 
2142 template<class Type>
2144 (
2145  const Type& shapes,
2146  const treeBoundBox& bb,
2147  const label maxLevels, // maximum number of levels
2148  const scalar maxLeafRatio,
2149  const scalar maxDuplicity
2150 )
2151 :
2152  shapes_(shapes),
2153  nodes_(),
2154  contents_(),
2155  nodeTypes_()
2156 {
2157  int64_t oldMemSize = 0;
2158  if (debug)
2159  {
2160  Pout<< "indexedOctree::indexedOctree:" << nl
2161  << " shapes:" << shapes.size() << nl
2162  << " bb:" << bb << nl
2163  << endl;
2164  oldMemSize = Foam::memInfo{}.size();
2165  }
2166 
2167  if (shapes.size() == 0)
2168  {
2169  return;
2170  }
2171 
2172  // Start off with one node with all shapes in it.
2173  DynamicList<node> nodes(label(shapes.size() / maxLeafRatio));
2174  DynamicList<labelList> contents(label(shapes.size() / maxLeafRatio));
2175  contents.push_back(identity(shapes.size()));
2176 
2177  // Create topnode.
2178  node topNode(divide(bb, contents, 0));
2179  nodes.push_back(topNode);
2180 
2181 
2182  // Now have all contents at level 1. Create levels by splitting levels
2183  // above.
2184 
2185  label nLevels = 1;
2186 
2187  for (; nLevels < maxLevels; nLevels++)
2188  {
2189  // Count number of references into shapes (i.e. contents)
2190  label nEntries = 0;
2191  label minEntries = shapes.size();
2192  label maxEntries = 0;
2193  for (const auto& subContents : contents)
2194  {
2195  const label num = subContents.size();
2196 
2197  nEntries += num;
2198  minEntries = min(minEntries, num);
2199  maxEntries = max(maxEntries, num);
2200  }
2201 
2202  if (debug)
2203  {
2204  Pout<< "indexedOctree::indexedOctree:" << nl
2205  << " nLevels:" << nLevels << nl
2206  << " nEntries:" << nEntries << nl
2207  << " - min per leaf: " << minEntries << nl
2208  << " - max per leaf: " << maxEntries << nl
2209  << " - avg per leaf: "
2210  << scalar(nEntries)/contents.size() << nl
2211  << " - per shape (duplicity): "
2212  << scalar(nEntries)/shapes_.size() << nl
2213  << endl;
2214  }
2215 
2216  if
2217  (
2218  //nEntries < maxLeafRatio*contents.size()
2219  // ||
2220  nEntries > maxDuplicity*shapes.size()
2221  )
2222  {
2223  break;
2224  }
2225 
2226 
2227  // Split nodes with more than maxLeafRatio elements
2228  label nOldNodes = nodes.size();
2229  splitNodes
2230  (
2231  label(maxLeafRatio),
2232  nodes,
2233  contents
2234  );
2235 
2236  if (nOldNodes == nodes.size())
2237  {
2238  break;
2239  }
2240  }
2241 
2242  // Shrink
2243  nodes.shrink();
2244  contents.shrink();
2245 
2246 
2247  // Compact such that deeper level contents are always after the
2248  // ones for a shallower level. This way we can slice a coarser level
2249  // off the tree.
2250  contents_.setSize(contents.size());
2251  label compactI = 0;
2252 
2253  label level = 0;
2254 
2255  while (true)
2256  {
2257  label nNodes = compactContents
2258  (
2259  nodes,
2260  contents,
2261  level,
2262  0,
2263  0,
2264  contents_,
2265  compactI
2266  );
2267 
2268  if (compactI == 0 && nNodes == 0)
2269  {
2270  // Did not put contents anywhere - are outside bb!
2271  break;
2272  }
2273 
2274  if (compactI == contents_.size())
2275  {
2276  // Transferred all contents to contents_ (in order breadth first)
2277  break;
2278  }
2279 
2280  level++;
2281  }
2282  nodes_.transfer(nodes);
2283  nodes.clear();
2284 
2285  if (debug)
2286  {
2287  label nEntries = 0;
2288  label minEntries = shapes.size();
2289  label maxEntries = 0;
2290  for (const auto& subContents : contents_)
2291  {
2292  const label num = subContents.size();
2293 
2294  nEntries += num;
2295  minEntries = min(minEntries, num);
2296  maxEntries = max(maxEntries, num);
2297  }
2298 
2299  int64_t memSize = Foam::memInfo{}.size();
2300 
2301  Pout<< "indexedOctree::indexedOctree"
2302  << " : finished construction of tree of:" << shapes.typeName
2303  << nl
2304  << " bb:" << this->bb() << nl
2305  << " shapes:" << shapes.size() << nl
2306  << " nLevels:" << nLevels << nl
2307  << " treeNodes:" << nodes_.size() << nl
2308  << " nEntries:" << nEntries << nl
2309  << " - min per leaf:" << minEntries << nl
2310  << " - max per leaf:" << maxEntries << nl
2311  << " - avg per leaf:"
2312  << scalar(nEntries)/contents.size() << nl
2313  << " - per shape (duplicity):"
2314  << scalar(nEntries)/shapes.size() << nl
2315  << " total memory:" << memSize-oldMemSize << nl
2316  << endl;
2317  }
2318 }
2319 
2320 
2321 template<class Type>
2323 (
2324  const Type& shapes,
2325  Istream& is
2326 )
2327 :
2328  shapes_(shapes),
2329  nodes_(is),
2330  contents_(is),
2331  nodeTypes_()
2332 {}
2333 
2334 
2335 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
2336 
2337 template<class Type>
2339 (
2340  const point& sample,
2341  const scalar startDistSqr
2342 ) const
2343 {
2344  return findNearest
2345  (
2346  sample,
2347  startDistSqr,
2348  typename Type::findNearestOp(*this)
2349  );
2350 }
2351 
2352 
2353 template<class Type>
2354 template<class FindNearestOp>
2356 (
2357  const point& sample,
2358  const scalar startDistSqr,
2359 
2360  const FindNearestOp& fnOp
2361 ) const
2362 {
2363  scalar nearestDistSqr = startDistSqr;
2364  label nearestShapeI = -1;
2365  point nearestPoint = Zero;
2366 
2367  if (nodes_.size())
2368  {
2369  findNearest
2370  (
2371  0,
2372  sample,
2373 
2374  nearestDistSqr,
2375  nearestShapeI,
2376  nearestPoint,
2377 
2378  fnOp
2379  );
2380  }
2382  return pointIndexHit(nearestShapeI != -1, nearestPoint, nearestShapeI);
2383 }
2384 
2385 
2386 template<class Type>
2388 (
2389  const linePointRef& ln,
2390  treeBoundBox& tightest,
2391  point& linePoint
2392 ) const
2393 {
2394  return findNearest
2395  (
2396  ln,
2397  tightest,
2398  linePoint,
2399  typename Type::findNearestOp(*this)
2400  );
2401 }
2402 
2403 
2404 template<class Type>
2405 template<class FindNearestOp>
2407 (
2408  const linePointRef& ln,
2409  treeBoundBox& tightest,
2410  point& linePoint,
2411 
2412  const FindNearestOp& fnOp
2413 ) const
2414 {
2415  label nearestShapeI = -1;
2416  point nearestPoint = Zero;
2417 
2418  if (nodes_.size())
2419  {
2420  findNearest
2421  (
2422  0,
2423  ln,
2424 
2425  tightest,
2426  nearestShapeI,
2427  linePoint,
2428  nearestPoint,
2429 
2430  fnOp
2431  );
2432  }
2434  return pointIndexHit(nearestShapeI != -1, nearestPoint, nearestShapeI);
2435 }
2436 
2437 
2438 template<class Type>
2440 (
2441  const point& start,
2442  const point& end
2443 ) const
2444 {
2445  return findLine
2446  (
2447  false,
2448  start,
2449  end,
2450  typename Type::findIntersectOp(*this)
2451  );
2452 }
2453 
2454 
2455 template<class Type>
2457 (
2458  const point& start,
2459  const point& end
2460 ) const
2461 {
2462  return findLine
2463  (
2464  true,
2465  start,
2466  end,
2467  typename Type::findIntersectOp(*this)
2468  );
2469 }
2470 
2471 
2472 template<class Type>
2473 template<class FindIntersectOp>
2475 (
2476  const point& start,
2477  const point& end,
2478  const FindIntersectOp& fiOp
2479 ) const
2480 {
2481  return findLine(false, start, end, fiOp);
2482 }
2483 
2484 
2485 template<class Type>
2486 template<class FindIntersectOp>
2488 (
2489  const point& start,
2490  const point& end,
2491  const FindIntersectOp& fiOp
2492 ) const
2494  return findLine(true, start, end, fiOp);
2495 }
2496 
2497 
2498 template<class Type>
2500 (
2501  const treeBoundBox& searchBox
2502 ) const
2503 {
2504  // start node=0, do not store
2505  return !nodes_.empty() && findBox(0, searchBox, nullptr);
2506 }
2507 
2508 
2509 template<class Type>
2511 (
2512  const treeBoundBox& searchBox,
2513  labelHashSet& elements
2514 ) const
2515 {
2516  elements.clear();
2517 
2518  if (!nodes_.empty())
2519  {
2520  // Some arbitrary minimal size estimate (eg, 1/100 are found)
2521  label estimatedCount(max(128, (shapes_.size() / 100)));
2522  elements.reserve(estimatedCount);
2523 
2524  // start node=0, store results
2525  findBox(0, searchBox, &elements);
2526  }
2528  return elements.size();
2529 }
2530 
2531 
2532 template<class Type>
2534 (
2535  const treeBoundBox& searchBox
2536 ) const
2537 {
2538  if (nodes_.empty())
2539  {
2540  return labelList();
2541  }
2542 
2543  labelHashSet elements(0);
2544 
2545  findBox(searchBox, elements);
2546 
2547  //TBD: return sorted ? elements.sortedToc() : elements.toc();
2548  return elements.toc();
2549 }
2550 
2551 
2552 template<class Type>
2554 (
2555  const point& centre,
2556  const scalar radiusSqr
2557 ) const
2558 {
2559  // start node=0, do not store
2560  return !nodes_.empty() && findSphere(0, centre, radiusSqr, nullptr);
2561 }
2562 
2563 
2564 template<class Type>
2566 (
2567  const point& centre,
2568  const scalar radiusSqr,
2569  labelHashSet& elements
2570 ) const
2571 {
2572  elements.clear();
2573 
2574  if (!nodes_.empty())
2575  {
2576  // Some arbitrary minimal size estimate (eg, 1/100 are found)
2577  label estimatedCount(max(128, (shapes_.size() / 100)));
2578  elements.reserve(estimatedCount);
2579 
2580  // start node=0, store results
2581  findSphere(0, centre, radiusSqr, &elements);
2582  }
2583 
2584  return elements.size();
2585 }
2586 
2587 
2588 
2589 template<class Type>
2591 (
2592  const point& centre,
2593  const scalar radiusSqr
2594 ) const
2595 {
2596  if (nodes_.empty())
2597  {
2598  return labelList();
2599  }
2600 
2601  labelHashSet elements(0);
2602 
2603  findSphere(centre, radiusSqr, elements);
2604 
2605  //TBD: return sorted ? elements.sortedToc() : elements.toc();
2606  return elements.toc();
2607 }
2608 
2609 
2610 template<class Type>
2612 (
2613  const label nodeI,
2614  const point& sample
2615 ) const
2616 {
2617  if (nodes_.empty())
2618  {
2619  // Empty tree. Return what?
2620  return nodePlusOctant(nodeI, 0);
2621  }
2622 
2623  const node& nod = nodes_[nodeI];
2624 
2625  direction octant = nod.bb_.subOctant(sample);
2626 
2627  labelBits index = nod.subNodes_[octant];
2628 
2629  if (isNode(index))
2630  {
2631  // Recurse
2632  return findNode(getNode(index), sample);
2633  }
2634  else if (isContent(index))
2635  {
2636  // Content. Return treenode+octant
2637  return nodePlusOctant(nodeI, octant);
2638  }
2639  else
2640  {
2641  // Empty. Return treenode+octant
2642  return nodePlusOctant(nodeI, octant);
2643  }
2644 }
2645 
2646 
2647 template<class Type>
2648 Foam::label Foam::indexedOctree<Type>::findInside(const point& sample) const
2649 {
2650  if (nodes_.empty())
2651  {
2652  return -1;
2653  }
2654 
2655  labelBits index = findNode(0, sample);
2656 
2657  const node& nod = nodes_[getNode(index)];
2658 
2659  labelBits contentIndex = nod.subNodes_[getOctant(index)];
2660 
2661  // Need to check for the presence of content, in-case the node is empty
2662  if (isContent(contentIndex))
2663  {
2664  labelList indices = contents_[getContent(contentIndex)];
2665 
2666  forAll(indices, elemI)
2667  {
2668  label shapeI = indices[elemI];
2669 
2670  if (shapes_.contains(shapeI, sample))
2671  {
2672  return shapeI;
2673  }
2674  }
2675  }
2677  return -1;
2678 }
2679 
2680 
2681 template<class Type>
2683 (
2684  const point& sample
2685 ) const
2686 {
2687  if (nodes_.empty())
2688  {
2689  return labelList::null();
2690  }
2691 
2692  labelBits index = findNode(0, sample);
2693 
2694  const node& nod = nodes_[getNode(index)];
2695 
2696  labelBits contentIndex = nod.subNodes_[getOctant(index)];
2697 
2698  // Need to check for the presence of content, in-case the node is empty
2699  if (isContent(contentIndex))
2700  {
2701  return contents_[getContent(contentIndex)];
2702  }
2704  return labelList::null();
2705 }
2706 
2707 
2708 template<class Type>
2710 (
2711  const point& sample
2712 ) const
2713 {
2714  if (nodes_.empty())
2715  {
2716  return volumeType::UNKNOWN;
2717  }
2718 
2719  if (nodeTypes_.size() != 8*nodes_.size())
2720  {
2721  // Calculate type for every octant of node.
2722 
2723  nodeTypes_.setSize(8*nodes_.size());
2724  nodeTypes_ = volumeType::UNKNOWN;
2725 
2726  calcVolumeType(0);
2727 
2728  if (debug)
2729  {
2730  label nUNKNOWN = 0;
2731  label nMIXED = 0;
2732  label nINSIDE = 0;
2733  label nOUTSIDE = 0;
2734 
2735  forAll(nodeTypes_, i)
2736  {
2737  volumeType type = volumeType::type(nodeTypes_.get(i));
2738 
2739  if (type == volumeType::UNKNOWN)
2740  {
2741  nUNKNOWN++;
2742  }
2743  else if (type == volumeType::MIXED)
2744  {
2745  nMIXED++;
2746  }
2747  else if (type == volumeType::INSIDE)
2748  {
2749  nINSIDE++;
2750  }
2751  else if (type == volumeType::OUTSIDE)
2752  {
2753  nOUTSIDE++;
2754  }
2755  else
2756  {
2758  }
2759  }
2760 
2761  Pout<< "indexedOctree::getVolumeType : "
2762  << " bb:" << bb()
2763  << " nodes_:" << nodes_.size()
2764  << " nodeTypes_:" << nodeTypes_.size()
2765  << " nUNKNOWN:" << nUNKNOWN
2766  << " nMIXED:" << nMIXED
2767  << " nINSIDE:" << nINSIDE
2768  << " nOUTSIDE:" << nOUTSIDE
2769  << endl;
2770  }
2771  }
2772 
2773  return getVolumeType(0, sample);
2774 }
2775 
2776 
2777 template<class Type>
2778 template<class CompareOp>
2780 (
2781  const scalar nearDist,
2782  const indexedOctree<Type>& tree2,
2783  CompareOp& cop
2784 ) const
2785 {
2786  if (!nodes_.empty())
2787  {
2788  findNear
2789  (
2790  nearDist,
2791  true,
2792  *this,
2793  nodePlusOctant(0, 0),
2794  bb(),
2795  tree2,
2796  nodePlusOctant(0, 0),
2797  tree2.bb(),
2798  cop
2799  );
2800  }
2801 }
2802 
2803 
2804 template<class Type>
2806 (
2807  prefixOSstream& os,
2808  const bool printContents,
2809  const label nodeI
2810 ) const
2811 {
2812  if (nodes_.empty())
2813  {
2814  return;
2815  }
2816 
2817  const node& nod = nodes_[nodeI];
2818  const treeBoundBox& bb = nod.bb_;
2819 
2820  os << "nodeI:" << nodeI << " bb:" << bb << nl
2821  << "parent:" << nod.parent_ << nl
2822  << "n:" << countElements(nodePlusOctant(nodeI, 0)) << nl;
2823 
2824  for (direction octant = 0; octant < node::nChildren; ++octant)
2825  {
2826  const treeBoundBox subBb(bb.subBbox(octant));
2827 
2828  labelBits index = nod.subNodes_[octant];
2829 
2830  if (isNode(index))
2831  {
2832  // tree node.
2833  label subNodeI = getNode(index);
2834 
2835  os << "octant:" << octant
2836  << " node: n:" << countElements(index)
2837  << " bb:" << subBb << endl;
2838 
2839  string oldPrefix = os.prefix();
2840  os.prefix() = " " + oldPrefix;
2841 
2842  print(os, printContents, subNodeI);
2843 
2844  os.prefix() = oldPrefix;
2845  }
2846  else if (isContent(index))
2847  {
2848  const labelList& indices = contents_[getContent(index)];
2849 
2850  if (debug)
2851  {
2852  writeOBJ(nodeI, octant);
2853  }
2854 
2855  os << "octant:" << octant
2856  << " content: n:" << indices.size()
2857  << " bb:" << subBb;
2858 
2859  if (printContents)
2860  {
2861  os << " contents:";
2862  forAll(indices, j)
2863  {
2864  os << ' ' << indices[j];
2865  }
2866  }
2867  os << endl;
2868  }
2869  else
2870  {
2871  if (debug)
2872  {
2873  writeOBJ(nodeI, octant);
2874  }
2875 
2876  os << "octant:" << octant << " empty:" << subBb << endl;
2877  }
2878  }
2879 }
2880 
2881 
2882 template<class Type>
2883 Foam::label Foam::indexedOctree<Type>::nLeafs() const
2884 {
2885  if (nodes_.size() < 2)
2886  {
2887  // If 0 or 1 nodes, treat directly as content nodes
2888  return nodes_.size();
2889  }
2890 
2891  return countLeafs(0);
2892 }
2893 
2894 
2895 template<class Type>
2896 bool Foam::indexedOctree<Type>::write(Ostream& os) const
2897 {
2898  os << *this;
2899 
2900  return os.good();
2901 }
2902 
2903 
2904 template<class Type>
2905 Foam::Ostream& Foam::operator<<(Ostream& os, const indexedOctree<Type>& t)
2906 {
2907  return
2908  os << t.bb() << token::SPACE << t.nodes()
2909  << token::SPACE << t.contents();
2910 }
2911 
2912 
2913 // ************************************************************************* //
type
Types of root.
Definition: Roots.H:52
void divide(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
uint8_t direction
Definition: direction.H:46
A line primitive.
Definition: line.H:52
bool overlaps(const treeBoundBox &bb) const
True if any shapes overlap the bounding box.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Mixed uniform/non-uniform (eg, after reduction)
Definition: ListPolicy.H:108
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
pointIndexHit findLineAny(const point &start, const point &end) const
Find any intersection of line between start and end.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
Ostream & print(Ostream &os, UIntType value, char off='0', char on='1')
Print 0/1 bits in the (unsigned) integral type.
Definition: BitOps.H:323
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
Definition: pointIndexHit.H:58
This class describes the interaction of an object (often a face) and a point. It carries the info of ...
Definition: pointIndexHit.H:44
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:196
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
static volumeType getSide(const vector &outsideNormal, const vector &vec)
Helper function to return the side. Returns outside if.
An enumeration wrapper for classification of a location as being inside/outside of a volume...
Definition: volumeType.H:55
labelBits findNode(const label nodeI, const point &) const
Find deepest node (as parent+octant) containing point. Starts.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
Various functions to operate on Lists.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
direction posBits(const point &pt) const
Position of point relative to bounding box.
Definition: treeBoundBox.C:486
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
label size() const noexcept
The number of elements in table.
Definition: HashTable.H:342
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:255
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
line< point, const point & > linePointRef
A line using referred points.
Definition: line.H:66
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
A class for handling words, derived from Foam::string.
Definition: word.H:63
void clear()
Remove all entries from table.
Definition: HashTable.C:725
Vector< scalar > vector
Definition: vector.H:57
const treeBoundBox & bb() const
Top bounding box.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
errorManip< error > abort(error &err)
Definition: errorManip.H:139
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: POSIX.C:1237
indexedOctree(const Type &shapes)
Construct null.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
void print(prefixOSstream &, const bool printContents, const label) const
Print tree. Either print all indices (printContent = true) or.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
label findInside(const point &) const
Find shape containing point. Only implemented for certain.
label facePoint(const int facei, const block &block, const label i, const label j)
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
OBJstream os(runTime.globalPath()/outputName)
Memory usage information for the current process, and the system memory that is free.
Definition: memInfo.H:58
const labelList & findIndices(const point &) const
Find the shape indices that occupy the result of findNode.
vector point
Point is a vector.
Definition: point.H:37
Non-pointer based hierarchical recursive searching.
#define WarningInFunction
Report a warning using Foam::Warning.
PtrList< volScalarField > & Y
void reserve(label numEntries)
Reserve space for at least the specified number of elements (not the number of buckets) and regenerat...
Definition: HashTable.C:712
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:90
A 29bits (or 61bits) integer label with 3bits direction (eg, octant) packed into single label...
Definition: labelBits.H:48
List< Key > toc() const
The table of contents (the keys) in unsorted order.
Definition: HashTable.C:124
List< label > labelList
A List of labels.
Definition: List.H:62
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
static unsigned getOctant(const point &p)
label nLeafs() const
Return the number of leaf nodes.
bool write(Ostream &os) const
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
int64_t size() const noexcept
Memory size at last update - (VmSize in /proc/PID/status)
Definition: memInfo.H:118
Tree node. Has up pointer and down pointers.
Definition: indexedOctree.H:76
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127