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