distributedTriSurfaceMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "mapDistribute.H"
28 #include "Random.H"
30 #include "triangleFuncs.H"
31 #include "matchPoints.H"
32 #include "globalIndex.H"
33 #include "Time.H"
34 
35 #include "IFstream.H"
36 #include "decompositionMethod.H"
37 #include "geomDecomp.H"
38 #include "vectorList.H"
39 #include "PackedBoolList.H"
40 #include "PatchTools.H"
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46  defineTypeNameAndDebug(distributedTriSurfaceMesh, 0);
48  (
49  searchableSurface,
50  distributedTriSurfaceMesh,
51  dict
52  );
53 
54  template<>
55  const char* Foam::NamedEnum
56  <
58  3
59  >::names[] =
60  {
61  "follow",
62  "independent",
63  "frozen"
64  };
65 }
66 
67 
70 
71 
72 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
73 
74 // Read my additional data from the dictionary
75 bool Foam::distributedTriSurfaceMesh::read()
76 {
77  // Get bb of all domains.
78  procBb_.setSize(Pstream::nProcs());
79 
80  procBb_[Pstream::myProcNo()] = List<treeBoundBox>(dict_.lookup("bounds"));
81  Pstream::gatherList(procBb_);
82  Pstream::scatterList(procBb_);
83 
84  // Distribution type
85  distType_ = distributionTypeNames_.read(dict_.lookup("distributionType"));
86 
87  // Merge distance
88  mergeDist_ = readScalar(dict_.lookup("mergeDistance"));
89 
90  return true;
91 }
92 
93 
94 // Is segment fully local?
95 bool Foam::distributedTriSurfaceMesh::isLocal
96 (
97  const List<treeBoundBox>& myBbs,
98  const point& start,
99  const point& end
100 )
101 {
102  forAll(myBbs, bbI)
103  {
104  if (myBbs[bbI].contains(start) && myBbs[bbI].contains(end))
105  {
106  return true;
107  }
108  }
109  return false;
110 }
111 
112 
113 //void Foam::distributedTriSurfaceMesh::splitSegment
114 //(
115 // const label segmentI,
116 // const point& start,
117 // const point& end,
118 // const treeBoundBox& bb,
119 //
120 // DynamicList<segment>& allSegments,
121 // DynamicList<label>& allSegmentMap,
122 // DynamicList<label> sendMap
123 //) const
124 //{
125 // // Work points
126 // point clipPt0, clipPt1;
127 //
128 // if (bb.contains(start))
129 // {
130 // // start within, trim end to bb
131 // bool clipped = bb.intersects(end, start, clipPt0);
132 //
133 // if (clipped)
134 // {
135 // // segment from start to clippedStart passes
136 // // through proc.
137 // sendMap[procI].append(allSegments.size());
138 // allSegmentMap.append(segmentI);
139 // allSegments.append(segment(start, clipPt0));
140 // }
141 // }
142 // else if (bb.contains(end))
143 // {
144 // // end within, trim start to bb
145 // bool clipped = bb.intersects(start, end, clipPt0);
146 //
147 // if (clipped)
148 // {
149 // sendMap[procI].append(allSegments.size());
150 // allSegmentMap.append(segmentI);
151 // allSegments.append(segment(clipPt0, end));
152 // }
153 // }
154 // else
155 // {
156 // // trim both
157 // bool clippedStart = bb.intersects(start, end, clipPt0);
158 //
159 // if (clippedStart)
160 // {
161 // bool clippedEnd = bb.intersects(end, clipPt0, clipPt1);
162 //
163 // if (clippedEnd)
164 // {
165 // // middle part of segment passes through proc.
166 // sendMap[procI].append(allSegments.size());
167 // allSegmentMap.append(segmentI);
168 // allSegments.append(segment(clipPt0, clipPt1));
169 // }
170 // }
171 // }
172 //}
173 
174 
175 void Foam::distributedTriSurfaceMesh::distributeSegment
176 (
177  const label segmentI,
178  const point& start,
179  const point& end,
180 
181  DynamicList<segment>& allSegments,
182  DynamicList<label>& allSegmentMap,
183  List<DynamicList<label> >& sendMap
184 ) const
185 {
186  // 1. Fully local already handled outside. Note: retest is cheap.
187  if (isLocal(procBb_[Pstream::myProcNo()], start, end))
188  {
189  return;
190  }
191 
192 
193  // 2. If fully inside one other processor, then only need to send
194  // to that one processor even if it intersects another. Rare occurrence
195  // but cheap to test.
196  forAll(procBb_, procI)
197  {
198  if (procI != Pstream::myProcNo())
199  {
200  const List<treeBoundBox>& bbs = procBb_[procI];
201 
202  if (isLocal(bbs, start, end))
203  {
204  sendMap[procI].append(allSegments.size());
205  allSegmentMap.append(segmentI);
206  allSegments.append(segment(start, end));
207  return;
208  }
209  }
210  }
211 
212 
213  // 3. If not contained in single processor send to all intersecting
214  // processors.
215  forAll(procBb_, procI)
216  {
217  const List<treeBoundBox>& bbs = procBb_[procI];
218 
219  forAll(bbs, bbI)
220  {
221  const treeBoundBox& bb = bbs[bbI];
222 
223  // Scheme a: any processor that intersects the segment gets
224  // the segment.
225 
226  // Intersection point
227  point clipPt;
228 
229  if (bb.intersects(start, end, clipPt))
230  {
231  sendMap[procI].append(allSegments.size());
232  allSegmentMap.append(segmentI);
233  allSegments.append(segment(start, end));
234  }
235 
236  // Alternative: any processor only gets clipped bit of
237  // segment. This gives small problems with additional
238  // truncation errors.
239  //splitSegment
240  //(
241  // segmentI,
242  // start,
243  // end,
244  // bb,
245  //
246  // allSegments,
247  // allSegmentMap,
248  // sendMap[procI]
249  //);
250  }
251  }
252 }
253 
254 
256 Foam::distributedTriSurfaceMesh::distributeSegments
257 (
258  const pointField& start,
259  const pointField& end,
260 
261  List<segment>& allSegments,
262  labelList& allSegmentMap
263 ) const
264 {
265  // Determine send map
266  // ~~~~~~~~~~~~~~~~~~
267 
268  labelListList sendMap(Pstream::nProcs());
269 
270  {
271  // Since intersection test is quite expensive compared to memory
272  // allocation we use DynamicList to immediately store the segment
273  // in the correct bin.
274 
275  // Segments to test
276  DynamicList<segment> dynAllSegments(start.size());
277  // Original index of segment
278  DynamicList<label> dynAllSegmentMap(start.size());
279  // Per processor indices into allSegments to send
281 
282  forAll(start, segmentI)
283  {
284  distributeSegment
285  (
286  segmentI,
287  start[segmentI],
288  end[segmentI],
289 
290  dynAllSegments,
291  dynAllSegmentMap,
292  dynSendMap
293  );
294  }
295 
296  // Convert dynamicList to labelList
297  sendMap.setSize(Pstream::nProcs());
298  forAll(sendMap, procI)
299  {
300  dynSendMap[procI].shrink();
301  sendMap[procI].transfer(dynSendMap[procI]);
302  }
303 
304  allSegments.transfer(dynAllSegments.shrink());
305  allSegmentMap.transfer(dynAllSegmentMap.shrink());
306  }
307 
308 
309  // Send over how many I need to receive.
310  labelListList sendSizes(Pstream::nProcs());
311  sendSizes[Pstream::myProcNo()].setSize(Pstream::nProcs());
312  forAll(sendMap, procI)
313  {
314  sendSizes[Pstream::myProcNo()][procI] = sendMap[procI].size();
315  }
316  Pstream::gatherList(sendSizes);
317  Pstream::scatterList(sendSizes);
318 
319 
320  // Determine order of receiving
321  labelListList constructMap(Pstream::nProcs());
322 
323  // My local segments first
324  constructMap[Pstream::myProcNo()] = identity
325  (
326  sendMap[Pstream::myProcNo()].size()
327  );
328 
329  label segmentI = constructMap[Pstream::myProcNo()].size();
330  forAll(constructMap, procI)
331  {
332  if (procI != Pstream::myProcNo())
333  {
334  // What I need to receive is what other processor is sending to me.
335  label nRecv = sendSizes[procI][Pstream::myProcNo()];
336  constructMap[procI].setSize(nRecv);
337 
338  for (label i = 0; i < nRecv; i++)
339  {
340  constructMap[procI][i] = segmentI++;
341  }
342  }
343  }
344 
346  (
347  new mapDistribute
348  (
349  segmentI, // size after construction
350  sendMap.xfer(),
351  constructMap.xfer()
352  )
353  );
354 }
355 
356 
357 void Foam::distributedTriSurfaceMesh::findLine
358 (
359  const bool nearestIntersection,
360  const pointField& start,
361  const pointField& end,
362  List<pointIndexHit>& info
363 ) const
364 {
365  const indexedOctree<treeDataTriSurface>& octree = tree();
366 
367  // Initialise
368  info.setSize(start.size());
369  forAll(info, i)
370  {
371  info[i].setMiss();
372  }
373 
374  if (!Pstream::parRun())
375  {
376  forAll(start, i)
377  {
378  if (nearestIntersection)
379  {
380  info[i] = octree.findLine(start[i], end[i]);
381  }
382  else
383  {
384  info[i] = octree.findLineAny(start[i], end[i]);
385  }
386  }
387  }
388  else
389  {
390  // Important:force synchronised construction of indexing
391  const globalIndex& triIndexer = globalTris();
392 
393 
394  // Do any local queries
395  // ~~~~~~~~~~~~~~~~~~~~
396 
397  label nLocal = 0;
398 
399  forAll(start, i)
400  {
401  if (isLocal(procBb_[Pstream::myProcNo()], start[i], end[i]))
402  {
403  if (nearestIntersection)
404  {
405  info[i] = octree.findLine(start[i], end[i]);
406  }
407  else
408  {
409  info[i] = octree.findLineAny(start[i], end[i]);
410  }
411 
412  if (info[i].hit())
413  {
414  info[i].setIndex(triIndexer.toGlobal(info[i].index()));
415  }
416  nLocal++;
417  }
418  }
419 
420 
421  if
422  (
423  returnReduce(nLocal, sumOp<label>())
424  < returnReduce(start.size(), sumOp<label>())
425  )
426  {
427  // Not all can be resolved locally. Build segments and map,
428  // send over segments, do intersections, send back and merge.
429 
430 
431  // Construct queries (segments)
432  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
433 
434  // Segments to test
435  List<segment> allSegments(start.size());
436  // Original index of segment
437  labelList allSegmentMap(start.size());
438 
439  const autoPtr<mapDistribute> mapPtr
440  (
441  distributeSegments
442  (
443  start,
444  end,
445  allSegments,
446  allSegmentMap
447  )
448  );
449  const mapDistribute& map = mapPtr();
450 
451  label nOldAllSegments = allSegments.size();
452 
453 
454  // Exchange the segments
455  // ~~~~~~~~~~~~~~~~~~~~~
456 
457  map.distribute(allSegments);
458 
459 
460  // Do tests I need to do
461  // ~~~~~~~~~~~~~~~~~~~~~
462 
463  // Intersections
464  List<pointIndexHit> intersections(allSegments.size());
465 
466  forAll(allSegments, i)
467  {
468  if (nearestIntersection)
469  {
470  intersections[i] = octree.findLine
471  (
472  allSegments[i].first(),
473  allSegments[i].second()
474  );
475  }
476  else
477  {
478  intersections[i] = octree.findLineAny
479  (
480  allSegments[i].first(),
481  allSegments[i].second()
482  );
483  }
484 
485  // Convert triangle index to global numbering
486  if (intersections[i].hit())
487  {
488  intersections[i].setIndex
489  (
490  triIndexer.toGlobal(intersections[i].index())
491  );
492  }
493  }
494 
495 
496  // Exchange the intersections (opposite to segments)
497  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
498 
499  map.reverseDistribute(nOldAllSegments, intersections);
500 
501 
502  // Extract the hits
503  // ~~~~~~~~~~~~~~~~
504 
505  forAll(intersections, i)
506  {
507  const pointIndexHit& allInfo = intersections[i];
508  label segmentI = allSegmentMap[i];
509  pointIndexHit& hitInfo = info[segmentI];
510 
511  if (allInfo.hit())
512  {
513  if (!hitInfo.hit())
514  {
515  // No intersection yet so take this one
516  hitInfo = allInfo;
517  }
518  else if (nearestIntersection)
519  {
520  // Nearest intersection
521  if
522  (
523  magSqr(allInfo.hitPoint()-start[segmentI])
524  < magSqr(hitInfo.hitPoint()-start[segmentI])
525  )
526  {
527  hitInfo = allInfo;
528  }
529  }
530  }
531  }
532  }
533  }
534 }
535 
536 
537 // Exchanges indices to the processor they come from.
538 // - calculates exchange map
539 // - uses map to calculate local triangle index
541 Foam::distributedTriSurfaceMesh::calcLocalQueries
542 (
543  const List<pointIndexHit>& info,
544  labelList& triangleIndex
545 ) const
546 {
547  triangleIndex.setSize(info.size());
548 
549  const globalIndex& triIndexer = globalTris();
550 
551 
552  // Determine send map
553  // ~~~~~~~~~~~~~~~~~~
554 
555  // Since determining which processor the query should go to is
556  // cheap we do a multi-pass algorithm to save some memory temporarily.
557 
558  // 1. Count
559  labelList nSend(Pstream::nProcs(), 0);
560 
561  forAll(info, i)
562  {
563  if (info[i].hit())
564  {
565  label procI = triIndexer.whichProcID(info[i].index());
566  nSend[procI]++;
567  }
568  }
569 
570  // 2. Size sendMap
571  labelListList sendMap(Pstream::nProcs());
572  forAll(nSend, procI)
573  {
574  sendMap[procI].setSize(nSend[procI]);
575  nSend[procI] = 0;
576  }
577 
578  // 3. Fill sendMap
579  forAll(info, i)
580  {
581  if (info[i].hit())
582  {
583  label procI = triIndexer.whichProcID(info[i].index());
584  triangleIndex[i] = triIndexer.toLocal(procI, info[i].index());
585  sendMap[procI][nSend[procI]++] = i;
586  }
587  else
588  {
589  triangleIndex[i] = -1;
590  }
591  }
592 
593 
594  // Send over how many I need to receive
595  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
596 
597  labelListList sendSizes(Pstream::nProcs());
598  sendSizes[Pstream::myProcNo()].setSize(Pstream::nProcs());
599  forAll(sendMap, procI)
600  {
601  sendSizes[Pstream::myProcNo()][procI] = sendMap[procI].size();
602  }
603  Pstream::gatherList(sendSizes);
604  Pstream::scatterList(sendSizes);
605 
606 
607  // Determine receive map
608  // ~~~~~~~~~~~~~~~~~~~~~
609 
610  labelListList constructMap(Pstream::nProcs());
611 
612  // My local segments first
613  constructMap[Pstream::myProcNo()] = identity
614  (
615  sendMap[Pstream::myProcNo()].size()
616  );
617 
618  label segmentI = constructMap[Pstream::myProcNo()].size();
619  forAll(constructMap, procI)
620  {
621  if (procI != Pstream::myProcNo())
622  {
623  // What I need to receive is what other processor is sending to me.
624  label nRecv = sendSizes[procI][Pstream::myProcNo()];
625  constructMap[procI].setSize(nRecv);
626 
627  for (label i = 0; i < nRecv; i++)
628  {
629  constructMap[procI][i] = segmentI++;
630  }
631  }
632  }
633 
634 
635  // Pack into distribution map
636  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
637 
639  (
640  new mapDistribute
641  (
642  segmentI, // size after construction
643  sendMap.xfer(),
644  constructMap.xfer()
645  )
646  );
647  const mapDistribute& map = mapPtr();
648 
649 
650  // Send over queries
651  // ~~~~~~~~~~~~~~~~~
652 
653  map.distribute(triangleIndex);
654 
655 
656  return mapPtr;
657 }
658 
659 
660 Foam::label Foam::distributedTriSurfaceMesh::calcOverlappingProcs
661 (
662  const point& centre,
663  const scalar radiusSqr,
664  boolList& overlaps
665 ) const
666 {
667  overlaps = false;
668  label nOverlaps = 0;
669 
670  forAll(procBb_, procI)
671  {
672  const List<treeBoundBox>& bbs = procBb_[procI];
673 
674  forAll(bbs, bbI)
675  {
676  if (bbs[bbI].overlaps(centre, radiusSqr))
677  {
678  overlaps[procI] = true;
679  nOverlaps++;
680  break;
681  }
682  }
683  }
684  return nOverlaps;
685 }
686 
687 
688 // Generate queries for parallel distance calculation
689 // - calculates exchange map
690 // - uses map to exchange points and radius
692 Foam::distributedTriSurfaceMesh::calcLocalQueries
693 (
694  const pointField& centres,
695  const scalarField& radiusSqr,
696 
697  pointField& allCentres,
698  scalarField& allRadiusSqr,
699  labelList& allSegmentMap
700 ) const
701 {
702  // Determine queries
703  // ~~~~~~~~~~~~~~~~~
704 
705  labelListList sendMap(Pstream::nProcs());
706 
707  {
708  // Queries
709  DynamicList<point> dynAllCentres(centres.size());
710  DynamicList<scalar> dynAllRadiusSqr(centres.size());
711  // Original index of segment
712  DynamicList<label> dynAllSegmentMap(centres.size());
713  // Per processor indices into allSegments to send
715 
716  // Work array - whether processor bb overlaps the bounding sphere.
717  boolList procBbOverlaps(Pstream::nProcs());
718 
719  forAll(centres, centreI)
720  {
721  // Find the processor this sample+radius overlaps.
722  calcOverlappingProcs
723  (
724  centres[centreI],
725  radiusSqr[centreI],
726  procBbOverlaps
727  );
728 
729  forAll(procBbOverlaps, procI)
730  {
731  if (procI != Pstream::myProcNo() && procBbOverlaps[procI])
732  {
733  dynSendMap[procI].append(dynAllCentres.size());
734  dynAllSegmentMap.append(centreI);
735  dynAllCentres.append(centres[centreI]);
736  dynAllRadiusSqr.append(radiusSqr[centreI]);
737  }
738  }
739  }
740 
741  // Convert dynamicList to labelList
742  sendMap.setSize(Pstream::nProcs());
743  forAll(sendMap, procI)
744  {
745  dynSendMap[procI].shrink();
746  sendMap[procI].transfer(dynSendMap[procI]);
747  }
748 
749  allCentres.transfer(dynAllCentres.shrink());
750  allRadiusSqr.transfer(dynAllRadiusSqr.shrink());
751  allSegmentMap.transfer(dynAllSegmentMap.shrink());
752  }
753 
754 
755  // Send over how many I need to receive.
756  labelListList sendSizes(Pstream::nProcs());
757  sendSizes[Pstream::myProcNo()].setSize(Pstream::nProcs());
758  forAll(sendMap, procI)
759  {
760  sendSizes[Pstream::myProcNo()][procI] = sendMap[procI].size();
761  }
762  Pstream::gatherList(sendSizes);
763  Pstream::scatterList(sendSizes);
764 
765 
766  // Determine order of receiving
767  labelListList constructMap(Pstream::nProcs());
768 
769  // My local segments first
770  constructMap[Pstream::myProcNo()] = identity
771  (
772  sendMap[Pstream::myProcNo()].size()
773  );
774 
775  label segmentI = constructMap[Pstream::myProcNo()].size();
776  forAll(constructMap, procI)
777  {
778  if (procI != Pstream::myProcNo())
779  {
780  // What I need to receive is what other processor is sending to me.
781  label nRecv = sendSizes[procI][Pstream::myProcNo()];
782  constructMap[procI].setSize(nRecv);
783 
784  for (label i = 0; i < nRecv; i++)
785  {
786  constructMap[procI][i] = segmentI++;
787  }
788  }
789  }
790 
792  (
793  new mapDistribute
794  (
795  segmentI, // size after construction
796  sendMap.xfer(),
797  constructMap.xfer()
798  )
799  );
800  return mapPtr;
801 }
802 
803 
804 // Find bounding boxes that guarantee a more or less uniform distribution
805 // of triangles. Decomposition in here is only used to get the bounding
806 // boxes, actual decomposition is done later on.
807 // Returns a per processor a list of bounding boxes that most accurately
808 // describe the shape. For now just a single bounding box per processor but
809 // optimisation might be to determine a better fitting shape.
811 Foam::distributedTriSurfaceMesh::independentlyDistributedBbs
812 (
813  const triSurface& s
814 )
815 {
816  if (!decomposer_.valid())
817  {
818  // Use current decomposer.
819  // Note: or always use hierarchical?
820  IOdictionary decomposeDict
821  (
822  IOobject
823  (
824  "decomposeParDict",
829  false
830  )
831  );
832  decomposer_ = decompositionMethod::New(decomposeDict);
833 
834  if (!decomposer_().parallelAware())
835  {
837  (
838  "distributedTriSurfaceMesh::independentlyDistributedBbs"
839  "(const triSurface&)"
840  ) << "The decomposition method " << decomposer_().typeName
841  << " does not decompose in parallel."
842  << " Please choose one that does." << exit(FatalError);
843  }
844 
845  if (!isA<geomDecomp>(decomposer_()))
846  {
848  (
849  "distributedTriSurfaceMesh::independentlyDistributedBbs"
850  "(const triSurface&)"
851  ) << "The decomposition method " << decomposer_().typeName
852  << " is not a geometric decomposition method." << endl
853  << "Only geometric decomposition methods are currently"
854  << " supported."
855  << exit(FatalError);
856  }
857  }
858 
859  // Do decomposition according to triangle centre
860  pointField triCentres(s.size());
861  forAll(s, triI)
862  {
863  triCentres[triI] = s[triI].centre(s.points());
864  }
865 
866 
867  geomDecomp& decomposer = refCast<geomDecomp>(decomposer_());
868 
869  // Do the actual decomposition
870  labelList distribution(decomposer.decompose(triCentres));
871 
872  // Find bounding box for all triangles on new distribution.
873 
874  // Initialise to inverted box (VGREAT, -VGREAT)
876  forAll(bbs, procI)
877  {
878  bbs[procI].setSize(1);
879  //bbs[procI][0] = boundBox::invertedBox;
880  bbs[procI][0].min() = point( VGREAT, VGREAT, VGREAT);
881  bbs[procI][0].max() = point(-VGREAT, -VGREAT, -VGREAT);
882  }
883 
884  forAll(s, triI)
885  {
886  point& bbMin = bbs[distribution[triI]][0].min();
887  point& bbMax = bbs[distribution[triI]][0].max();
888 
889  const triSurface::FaceType& f = s[triI];
890  forAll(f, fp)
891  {
892  const point& pt = s.points()[f[fp]];
893  bbMin = ::Foam::min(bbMin, pt);
894  bbMax = ::Foam::max(bbMax, pt);
895  }
896  }
897 
898  // Now combine for all processors and convert to correct format.
899  forAll(bbs, procI)
900  {
901  forAll(bbs[procI], i)
902  {
903  reduce(bbs[procI][i].min(), minOp<point>());
904  reduce(bbs[procI][i].max(), maxOp<point>());
905  }
906  }
907  return bbs;
908 }
909 
910 
911 // Does any part of triangle overlap bb.
912 bool Foam::distributedTriSurfaceMesh::overlaps
913 (
914  const List<treeBoundBox>& bbs,
915  const point& p0,
916  const point& p1,
917  const point& p2
918 )
919 {
920  forAll(bbs, bbI)
921  {
922  const treeBoundBox& bb = bbs[bbI];
923 
924  treeBoundBox triBb(p0, p0);
925  triBb.min() = min(triBb.min(), p1);
926  triBb.min() = min(triBb.min(), p2);
927 
928  triBb.max() = max(triBb.max(), p1);
929  triBb.max() = max(triBb.max(), p2);
930 
931  //- Exact test of triangle intersecting bb
932 
933  // Quick rejection. If whole bounding box of tri is outside cubeBb then
934  // there will be no intersection.
935  if (bb.overlaps(triBb))
936  {
937  // Check if one or more triangle point inside
938  if (bb.contains(p0) || bb.contains(p1) || bb.contains(p2))
939  {
940  // One or more points inside
941  return true;
942  }
943 
944  // Now we have the difficult case: all points are outside but
945  // connecting edges might go through cube. Use fast intersection
946  // of bounding box.
947  bool intersect = triangleFuncs::intersectBb(p0, p1, p2, bb);
948 
949  if (intersect)
950  {
951  return true;
952  }
953  }
954  }
955  return false;
956 }
957 
958 
959 void Foam::distributedTriSurfaceMesh::subsetMeshMap
960 (
961  const triSurface& s,
962  const boolList& include,
963  const label nIncluded,
964  labelList& newToOldPoints,
965  labelList& oldToNewPoints,
966  labelList& newToOldFaces
967 )
968 {
969  newToOldFaces.setSize(nIncluded);
970  newToOldPoints.setSize(s.points().size());
971  oldToNewPoints.setSize(s.points().size());
972  oldToNewPoints = -1;
973  {
974  label faceI = 0;
975  label pointI = 0;
976 
977  forAll(include, oldFacei)
978  {
979  if (include[oldFacei])
980  {
981  // Store new faces compact
982  newToOldFaces[faceI++] = oldFacei;
983 
984  // Renumber labels for face
985  const triSurface::FaceType& f = s[oldFacei];
986 
987  forAll(f, fp)
988  {
989  label oldPointI = f[fp];
990 
991  if (oldToNewPoints[oldPointI] == -1)
992  {
993  oldToNewPoints[oldPointI] = pointI;
994  newToOldPoints[pointI++] = oldPointI;
995  }
996  }
997  }
998  }
999  newToOldPoints.setSize(pointI);
1000  }
1001 }
1002 
1003 
1004 Foam::triSurface Foam::distributedTriSurfaceMesh::subsetMesh
1005 (
1006  const triSurface& s,
1007  const labelList& newToOldPoints,
1008  const labelList& oldToNewPoints,
1009  const labelList& newToOldFaces
1010 )
1011 {
1012  // Extract points
1013  pointField newPoints(newToOldPoints.size());
1014  forAll(newToOldPoints, i)
1015  {
1016  newPoints[i] = s.points()[newToOldPoints[i]];
1017  }
1018  // Extract faces
1019  List<labelledTri> newTriangles(newToOldFaces.size());
1020 
1021  forAll(newToOldFaces, i)
1022  {
1023  // Get old vertex labels
1024  const labelledTri& tri = s[newToOldFaces[i]];
1025 
1026  newTriangles[i][0] = oldToNewPoints[tri[0]];
1027  newTriangles[i][1] = oldToNewPoints[tri[1]];
1028  newTriangles[i][2] = oldToNewPoints[tri[2]];
1029  newTriangles[i].region() = tri.region();
1030  }
1031 
1032  // Reuse storage.
1033  return triSurface(newTriangles, s.patches(), newPoints, true);
1034 }
1035 
1036 
1037 Foam::triSurface Foam::distributedTriSurfaceMesh::subsetMesh
1038 (
1039  const triSurface& s,
1040  const boolList& include,
1041  labelList& newToOldPoints,
1042  labelList& newToOldFaces
1043 )
1044 {
1045  label n = 0;
1046 
1047  forAll(include, i)
1048  {
1049  if (include[i])
1050  {
1051  n++;
1052  }
1053  }
1054 
1055  labelList oldToNewPoints;
1056  subsetMeshMap
1057  (
1058  s,
1059  include,
1060  n,
1061  newToOldPoints,
1062  oldToNewPoints,
1063  newToOldFaces
1064  );
1065 
1066  return subsetMesh
1067  (
1068  s,
1069  newToOldPoints,
1070  oldToNewPoints,
1071  newToOldFaces
1072  );
1073 }
1074 
1075 
1076 Foam::triSurface Foam::distributedTriSurfaceMesh::subsetMesh
1077 (
1078  const triSurface& s,
1079  const labelList& newToOldFaces,
1080  labelList& newToOldPoints
1081 )
1082 {
1083  const boolList include
1084  (
1085  createWithValues<boolList>
1086  (
1087  s.size(),
1088  false,
1089  newToOldFaces,
1090  true
1091  )
1092  );
1093 
1094  newToOldPoints.setSize(s.points().size());
1095  labelList oldToNewPoints(s.points().size(), -1);
1096  {
1097  label pointI = 0;
1098 
1099  forAll(include, oldFacei)
1100  {
1101  if (include[oldFacei])
1102  {
1103  // Renumber labels for face
1104  const triSurface::FaceType& f = s[oldFacei];
1105 
1106  forAll(f, fp)
1107  {
1108  label oldPointI = f[fp];
1109 
1110  if (oldToNewPoints[oldPointI] == -1)
1111  {
1112  oldToNewPoints[oldPointI] = pointI;
1113  newToOldPoints[pointI++] = oldPointI;
1114  }
1115  }
1116  }
1117  }
1118  newToOldPoints.setSize(pointI);
1119  }
1120 
1121  return subsetMesh
1122  (
1123  s,
1124  newToOldPoints,
1125  oldToNewPoints,
1126  newToOldFaces
1127  );
1128 }
1129 
1130 
1131 Foam::label Foam::distributedTriSurfaceMesh::findTriangle
1132 (
1133  const List<labelledTri>& allFaces,
1134  const labelListList& allPointFaces,
1135  const labelledTri& otherF
1136 )
1137 {
1138  // allFaces connected to otherF[0]
1139  const labelList& pFaces = allPointFaces[otherF[0]];
1140 
1141  forAll(pFaces, i)
1142  {
1143  const labelledTri& f = allFaces[pFaces[i]];
1144 
1145  if (f.region() == otherF.region())
1146  {
1147  // Find index of otherF[0]
1148  label fp0 = findIndex(f, otherF[0]);
1149  // Check rest of triangle in same order
1150  label fp1 = f.fcIndex(fp0);
1151  label fp2 = f.fcIndex(fp1);
1152 
1153  if (f[fp1] == otherF[1] && f[fp2] == otherF[2])
1154  {
1155  return pFaces[i];
1156  }
1157  }
1158  }
1159  return -1;
1160 }
1161 
1162 
1163 // Merge into allSurf.
1164 void Foam::distributedTriSurfaceMesh::merge
1165 (
1166  const scalar mergeDist,
1167  const List<labelledTri>& subTris,
1168  const pointField& subPoints,
1169 
1170  List<labelledTri>& allTris,
1171  pointField& allPoints,
1172 
1173  labelList& faceConstructMap,
1174  labelList& pointConstructMap
1175 )
1176 {
1177  labelList subToAll;
1178  matchPoints
1179  (
1180  subPoints,
1181  allPoints,
1182  scalarField(subPoints.size(), mergeDist), // match distance
1183  false, // verbose
1184  pointConstructMap
1185  );
1186 
1187  label nOldAllPoints = allPoints.size();
1188 
1189 
1190  // Add all unmatched points
1191  // ~~~~~~~~~~~~~~~~~~~~~~~~
1192 
1193  label allPointI = nOldAllPoints;
1194  forAll(pointConstructMap, pointI)
1195  {
1196  if (pointConstructMap[pointI] == -1)
1197  {
1198  pointConstructMap[pointI] = allPointI++;
1199  }
1200  }
1201 
1202  if (allPointI > nOldAllPoints)
1203  {
1204  allPoints.setSize(allPointI);
1205 
1206  forAll(pointConstructMap, pointI)
1207  {
1208  if (pointConstructMap[pointI] >= nOldAllPoints)
1209  {
1210  allPoints[pointConstructMap[pointI]] = subPoints[pointI];
1211  }
1212  }
1213  }
1214 
1215 
1216  // To check whether triangles are same we use pointFaces.
1217  labelListList allPointFaces;
1218  invertManyToMany(nOldAllPoints, allTris, allPointFaces);
1219 
1220 
1221  // Add all unmatched triangles
1222  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
1223 
1224  label allTriI = allTris.size();
1225  allTris.setSize(allTriI + subTris.size());
1226 
1227  faceConstructMap.setSize(subTris.size());
1228 
1229  forAll(subTris, triI)
1230  {
1231  const labelledTri& subTri = subTris[triI];
1232 
1233  // Get triangle in new numbering
1234  labelledTri mappedTri
1235  (
1236  pointConstructMap[subTri[0]],
1237  pointConstructMap[subTri[1]],
1238  pointConstructMap[subTri[2]],
1239  subTri.region()
1240  );
1241 
1242 
1243  // Check if all points were already in surface
1244  bool fullMatch = true;
1245 
1246  forAll(mappedTri, fp)
1247  {
1248  if (mappedTri[fp] >= nOldAllPoints)
1249  {
1250  fullMatch = false;
1251  break;
1252  }
1253  }
1254 
1255  if (fullMatch)
1256  {
1257  // All three points are mapped to old points. See if same
1258  // triangle.
1259  label i = findTriangle
1260  (
1261  allTris,
1262  allPointFaces,
1263  mappedTri
1264  );
1265 
1266  if (i == -1)
1267  {
1268  // Add
1269  faceConstructMap[triI] = allTriI;
1270  allTris[allTriI] = mappedTri;
1271  allTriI++;
1272  }
1273  else
1274  {
1275  faceConstructMap[triI] = i;
1276  }
1277  }
1278  else
1279  {
1280  // Add
1281  faceConstructMap[triI] = allTriI;
1282  allTris[allTriI] = mappedTri;
1283  allTriI++;
1284  }
1285  }
1286  allTris.setSize(allTriI);
1287 }
1288 
1289 
1290 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1291 
1292 Foam::distributedTriSurfaceMesh::distributedTriSurfaceMesh
1294  const IOobject& io,
1295  const triSurface& s,
1296  const dictionary& dict
1297 )
1298 :
1299  triSurfaceMesh(io, s),
1300  dict_
1301  (
1302  IOobject
1303  (
1304  searchableSurface::name() + "Dict",
1311  ),
1312  dict
1313  )
1314 {
1315  read();
1316 
1317  reduce(bounds().min(), minOp<point>());
1318  reduce(bounds().max(), maxOp<point>());
1319 
1320  if (debug)
1321  {
1322  Info<< "Constructed from triSurface:" << endl;
1323  writeStats(Info);
1324 
1325  labelList nTris(Pstream::nProcs());
1326  nTris[Pstream::myProcNo()] = triSurface::size();
1327  Pstream::gatherList(nTris);
1328  Pstream::scatterList(nTris);
1329 
1330  Info<< endl<< "\tproc\ttris\tbb" << endl;
1331  forAll(nTris, procI)
1332  {
1333  Info<< '\t' << procI << '\t' << nTris[procI]
1334  << '\t' << procBb_[procI] << endl;
1335  }
1336  Info<< endl;
1337  }
1338 }
1339 
1340 
1341 Foam::distributedTriSurfaceMesh::distributedTriSurfaceMesh(const IOobject& io)
1342 :
1343  //triSurfaceMesh(io),
1345  (
1346  IOobject
1347  (
1348  io.name(),
1349  io.time().findInstance(io.local(), word::null),
1350  io.local(),
1351  io.db(),
1352  io.readOpt(),
1353  io.writeOpt(),
1354  io.registerObject()
1355  )
1356  ),
1357  dict_
1358  (
1359  IOobject
1360  (
1361  searchableSurface::name() + "Dict",
1362  searchableSurface::instance(),
1363  searchableSurface::local(),
1364  searchableSurface::db(),
1365  searchableSurface::readOpt(),
1366  searchableSurface::writeOpt(),
1367  searchableSurface::registerObject()
1368  )
1369  )
1370 {
1371  if
1372  (
1373  Pstream::parRun()
1374  && (
1375  dict_.readOpt() == IOobject::MUST_READ
1377  )
1378  && (
1381  )
1382  )
1383  {
1384  FatalErrorIn("Foam::distributedTriSurfaceMesh::read()")
1385  << " distributedTriSurfaceMesh is being constructed\n"
1386  << " using 'timeStampMaster' or 'inotifyMaster.'\n"
1387  << " Modify the entry fileModificationChecking\n"
1388  << " in the etc/controlDict.\n"
1389  << " Use 'timeStamp' instead."
1390  << exit(FatalError);
1391  }
1392 
1393  read();
1394 
1395  reduce(bounds().min(), minOp<point>());
1396  reduce(bounds().max(), maxOp<point>());
1397 
1398  if (debug)
1399  {
1400  Info<< "Read distributedTriSurface from " << io.objectPath()
1401  << ':' << endl;
1402  writeStats(Info);
1403 
1404  labelList nTris(Pstream::nProcs());
1405  nTris[Pstream::myProcNo()] = triSurface::size();
1406  Pstream::gatherList(nTris);
1407  Pstream::scatterList(nTris);
1408 
1409  Info<< endl<< "\tproc\ttris\tbb" << endl;
1410  forAll(nTris, procI)
1411  {
1412  Info<< '\t' << procI << '\t' << nTris[procI]
1413  << '\t' << procBb_[procI] << endl;
1414  }
1415  Info<< endl;
1416  }
1417 }
1418 
1419 
1420 Foam::distributedTriSurfaceMesh::distributedTriSurfaceMesh
1422  const IOobject& io,
1423  const dictionary& dict
1424 )
1425 :
1426  //triSurfaceMesh(io, dict),
1428  (
1429  IOobject
1430  (
1431  io.name(),
1432  io.time().findInstance(io.local(), word::null),
1433  io.local(),
1434  io.db(),
1435  io.readOpt(),
1436  io.writeOpt(),
1437  io.registerObject()
1438  ),
1439  dict
1440  ),
1441  dict_
1442  (
1443  IOobject
1444  (
1445  searchableSurface::name() + "Dict",
1452  )
1453  )
1454 {
1455  if
1456  (
1457  Pstream::parRun()
1458  && (
1459  dict_.readOpt() == IOobject::MUST_READ
1461  )
1462  && (
1465  )
1466  )
1467  {
1468  FatalErrorIn("Foam::distributedTriSurfaceMesh::read()")
1469  << " distributedTriSurfaceMesh is being constructed\n"
1470  << " using 'timeStampMaster' or 'inotifyMaster.'\n"
1471  << " Modify the entry fileModificationChecking\n"
1472  << " in the etc/controlDict.\n"
1473  << " Use 'timeStamp' instead."
1474  << exit(FatalError);
1475  }
1476 
1477  read();
1478 
1479  reduce(bounds().min(), minOp<point>());
1480  reduce(bounds().max(), maxOp<point>());
1481 
1482  if (debug)
1483  {
1484  Info<< "Read distributedTriSurface from " << io.objectPath()
1485  << " and dictionary:" << endl;
1486  writeStats(Info);
1487 
1488  labelList nTris(Pstream::nProcs());
1489  nTris[Pstream::myProcNo()] = triSurface::size();
1490  Pstream::gatherList(nTris);
1491  Pstream::scatterList(nTris);
1492 
1493  Info<< endl<< "\tproc\ttris\tbb" << endl;
1494  forAll(nTris, procI)
1495  {
1496  Info<< '\t' << procI << '\t' << nTris[procI]
1497  << '\t' << procBb_[procI] << endl;
1498  }
1499  Info<< endl;
1500  }
1501 }
1502 
1503 
1504 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1505 
1507 {
1508  clearOut();
1509 }
1510 
1511 
1513 {
1514  globalTris_.clear();
1516 }
1517 
1518 
1519 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1520 
1522 {
1523  if (!globalTris_.valid())
1524  {
1525  globalTris_.reset(new globalIndex(triSurface::size()));
1526  }
1527  return globalTris_;
1528 }
1529 
1530 
1533  const pointField& samples,
1534  const scalarField& nearestDistSqr,
1535  List<pointIndexHit>& info
1536 ) const
1537 {
1538  const indexedOctree<treeDataTriSurface>& octree = tree();
1539 
1540  // Important:force synchronised construction of indexing
1541  const globalIndex& triIndexer = globalTris();
1542 
1543 
1544  // Initialise
1545  // ~~~~~~~~~~
1546 
1547  info.setSize(samples.size());
1548  forAll(info, i)
1549  {
1550  info[i].setMiss();
1551  }
1552 
1553 
1554 
1555  // Do any local queries
1556  // ~~~~~~~~~~~~~~~~~~~~
1557 
1558  label nLocal = 0;
1559 
1560  {
1561  // Work array - whether processor bb overlaps the bounding sphere.
1562  boolList procBbOverlaps(Pstream::nProcs());
1563 
1564  forAll(samples, i)
1565  {
1566  // Find the processor this sample+radius overlaps.
1567  label nProcs = calcOverlappingProcs
1568  (
1569  samples[i],
1570  nearestDistSqr[i],
1571  procBbOverlaps
1572  );
1573 
1574  // Overlaps local processor?
1575  if (procBbOverlaps[Pstream::myProcNo()])
1576  {
1577  info[i] = octree.findNearest(samples[i], nearestDistSqr[i]);
1578  if (info[i].hit())
1579  {
1580  info[i].setIndex(triIndexer.toGlobal(info[i].index()));
1581  }
1582  if (nProcs == 1)
1583  {
1584  // Fully local
1585  nLocal++;
1586  }
1587  }
1588  }
1589  }
1590 
1591 
1592  if
1593  (
1594  Pstream::parRun()
1595  && (
1596  returnReduce(nLocal, sumOp<label>())
1597  < returnReduce(samples.size(), sumOp<label>())
1598  )
1599  )
1600  {
1601  // Not all can be resolved locally. Build queries and map, send over
1602  // queries, do intersections, send back and merge.
1603 
1604  // Calculate queries and exchange map
1605  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1606 
1607  pointField allCentres;
1608  scalarField allRadiusSqr;
1609  labelList allSegmentMap;
1610  autoPtr<mapDistribute> mapPtr
1611  (
1612  calcLocalQueries
1613  (
1614  samples,
1615  nearestDistSqr,
1616 
1617  allCentres,
1618  allRadiusSqr,
1619  allSegmentMap
1620  )
1621  );
1622  const mapDistribute& map = mapPtr();
1623 
1624 
1625  // swap samples to local processor
1626  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1627 
1628  map.distribute(allCentres);
1629  map.distribute(allRadiusSqr);
1630 
1631 
1632  // Do my tests
1633  // ~~~~~~~~~~~
1634 
1635  List<pointIndexHit> allInfo(allCentres.size());
1636  forAll(allInfo, i)
1637  {
1638  allInfo[i] = octree.findNearest
1639  (
1640  allCentres[i],
1641  allRadiusSqr[i]
1642  );
1643  if (allInfo[i].hit())
1644  {
1645  allInfo[i].setIndex(triIndexer.toGlobal(allInfo[i].index()));
1646  }
1647  }
1648 
1649 
1650  // Send back results
1651  // ~~~~~~~~~~~~~~~~~
1652 
1653  map.reverseDistribute(allSegmentMap.size(), allInfo);
1654 
1655 
1656  // Extract information
1657  // ~~~~~~~~~~~~~~~~~~~
1658 
1659  forAll(allInfo, i)
1660  {
1661  if (allInfo[i].hit())
1662  {
1663  label pointI = allSegmentMap[i];
1664 
1665  if (!info[pointI].hit())
1666  {
1667  // No intersection yet so take this one
1668  info[pointI] = allInfo[i];
1669  }
1670  else
1671  {
1672  // Nearest intersection
1673  if
1674  (
1675  magSqr(allInfo[i].hitPoint()-samples[pointI])
1676  < magSqr(info[pointI].hitPoint()-samples[pointI])
1677  )
1678  {
1679  info[pointI] = allInfo[i];
1680  }
1681  }
1682  }
1683  }
1684  }
1685 }
1686 
1687 
1688 void Foam::distributedTriSurfaceMesh::findLine
1690  const pointField& start,
1691  const pointField& end,
1692  List<pointIndexHit>& info
1693 ) const
1694 {
1695  findLine
1696  (
1697  true, // nearestIntersection
1698  start,
1699  end,
1700  info
1701  );
1702 }
1703 
1704 
1707  const pointField& start,
1708  const pointField& end,
1709  List<pointIndexHit>& info
1710 ) const
1711 {
1712  findLine
1713  (
1714  true, // nearestIntersection
1715  start,
1716  end,
1717  info
1718  );
1719 }
1720 
1721 
1724  const pointField& start,
1725  const pointField& end,
1726  List<List<pointIndexHit> >& info
1727 ) const
1728 {
1729  // Reuse fineLine. We could modify all of findLine to do multiple
1730  // intersections but this would mean a lot of data transferred so
1731  // for now we just find nearest intersection and retest from that
1732  // intersection onwards.
1733 
1734  // Work array.
1735  List<pointIndexHit> hitInfo(start.size());
1736 
1737  findLine
1738  (
1739  true, // nearestIntersection
1740  start,
1741  end,
1742  hitInfo
1743  );
1744 
1745  // Tolerances:
1746  // To find all intersections we add a small vector to the last intersection
1747  // This is chosen such that
1748  // - it is significant (SMALL is smallest representative relative tolerance;
1749  // we need something bigger since we're doing calculations)
1750  // - if the start-end vector is zero we still progress
1751  const vectorField dirVec(end-start);
1752  const scalarField magSqrDirVec(magSqr(dirVec));
1753  const vectorField smallVec
1754  (
1755  ROOTSMALL*dirVec
1756  + vector(ROOTVSMALL,ROOTVSMALL,ROOTVSMALL)
1757  );
1758 
1759  // Copy to input and compact any hits
1760  labelList pointMap(start.size());
1761  pointField e0(start.size());
1762  pointField e1(start.size());
1763  label compactI = 0;
1764 
1765  info.setSize(hitInfo.size());
1766  forAll(hitInfo, pointI)
1767  {
1768  if (hitInfo[pointI].hit())
1769  {
1770  info[pointI].setSize(1);
1771  info[pointI][0] = hitInfo[pointI];
1772 
1773  point pt = hitInfo[pointI].hitPoint() + smallVec[pointI];
1774 
1775  if (((pt-start[pointI])&dirVec[pointI]) <= magSqrDirVec[pointI])
1776  {
1777  e0[compactI] = pt;
1778  e1[compactI] = end[pointI];
1779  pointMap[compactI] = pointI;
1780  compactI++;
1781  }
1782  }
1783  else
1784  {
1785  info[pointI].clear();
1786  }
1787  }
1788 
1789  e0.setSize(compactI);
1790  e1.setSize(compactI);
1791  pointMap.setSize(compactI);
1792 
1793  while (returnReduce(e0.size(), sumOp<label>()) > 0)
1794  {
1795  findLine
1796  (
1797  true, // nearestIntersection
1798  e0,
1799  e1,
1800  hitInfo
1801  );
1802 
1803 
1804  // Extract
1805  label compactI = 0;
1806  forAll(hitInfo, i)
1807  {
1808  if (hitInfo[i].hit())
1809  {
1810  label pointI = pointMap[i];
1811 
1812  label sz = info[pointI].size();
1813  info[pointI].setSize(sz+1);
1814  info[pointI][sz] = hitInfo[i];
1815 
1816  point pt = hitInfo[i].hitPoint() + smallVec[pointI];
1817 
1818  if (((pt-start[pointI])&dirVec[pointI]) <= magSqrDirVec[pointI])
1819  {
1820  e0[compactI] = pt;
1821  e1[compactI] = end[pointI];
1822  pointMap[compactI] = pointI;
1823  compactI++;
1824  }
1825  }
1826  }
1827 
1828  // Trim
1829  e0.setSize(compactI);
1830  e1.setSize(compactI);
1831  pointMap.setSize(compactI);
1832  }
1833 }
1834 
1835 
1838  const List<pointIndexHit>& info,
1839  labelList& region
1840 ) const
1841 {
1842  if (!Pstream::parRun())
1843  {
1844  region.setSize(info.size());
1845  forAll(info, i)
1846  {
1847  if (info[i].hit())
1848  {
1849  region[i] = triSurface::operator[](info[i].index()).region();
1850  }
1851  else
1852  {
1853  region[i] = -1;
1854  }
1855  }
1856  return;
1857  }
1858 
1859 
1860  // Get query data (= local index of triangle)
1861  // ~~~~~~~~~~~~~~
1862 
1863  labelList triangleIndex(info.size());
1864  autoPtr<mapDistribute> mapPtr
1865  (
1866  calcLocalQueries
1867  (
1868  info,
1869  triangleIndex
1870  )
1871  );
1872  const mapDistribute& map = mapPtr();
1873 
1874 
1875  // Do my tests
1876  // ~~~~~~~~~~~
1877 
1878  const triSurface& s = static_cast<const triSurface&>(*this);
1879 
1880  region.setSize(triangleIndex.size());
1881 
1882  forAll(triangleIndex, i)
1883  {
1884  label triI = triangleIndex[i];
1885  region[i] = s[triI].region();
1886  }
1887 
1888 
1889  // Send back results
1890  // ~~~~~~~~~~~~~~~~~
1891 
1892  map.reverseDistribute(info.size(), region);
1893 }
1894 
1895 
1898  const List<pointIndexHit>& info,
1899  vectorField& normal
1900 ) const
1901 {
1902  if (!Pstream::parRun())
1903  {
1904  triSurfaceMesh::getNormal(info, normal);
1905  return;
1906  }
1907 
1908 
1909  // Get query data (= local index of triangle)
1910  // ~~~~~~~~~~~~~~
1911 
1912  labelList triangleIndex(info.size());
1913  autoPtr<mapDistribute> mapPtr
1914  (
1915  calcLocalQueries
1916  (
1917  info,
1918  triangleIndex
1919  )
1920  );
1921  const mapDistribute& map = mapPtr();
1922 
1923 
1924  // Do my tests
1925  // ~~~~~~~~~~~
1926 
1927  const triSurface& s = static_cast<const triSurface&>(*this);
1928 
1929  normal.setSize(triangleIndex.size());
1930 
1931  forAll(triangleIndex, i)
1932  {
1933  label triI = triangleIndex[i];
1934  normal[i] = s[triI].normal(s.points());
1935  normal[i] /= mag(normal[i]) + VSMALL;
1936  }
1937 
1938 
1939  // Send back results
1940  // ~~~~~~~~~~~~~~~~~
1941 
1942  map.reverseDistribute(info.size(), normal);
1943 }
1944 
1945 
1948  const List<pointIndexHit>& info,
1949  labelList& values
1950 ) const
1951 {
1952  if (!Pstream::parRun())
1953  {
1954  triSurfaceMesh::getField(info, values);
1955  return;
1956  }
1957 
1958  if (foundObject<triSurfaceLabelField>("values"))
1959  {
1960  const triSurfaceLabelField& fld = lookupObject<triSurfaceLabelField>
1961  (
1962  "values"
1963  );
1964 
1965 
1966  // Get query data (= local index of triangle)
1967  // ~~~~~~~~~~~~~~
1968 
1969  labelList triangleIndex(info.size());
1970  autoPtr<mapDistribute> mapPtr
1971  (
1972  calcLocalQueries
1973  (
1974  info,
1975  triangleIndex
1976  )
1977  );
1978  const mapDistribute& map = mapPtr();
1979 
1980 
1981  // Do my tests
1982  // ~~~~~~~~~~~
1983 
1984  values.setSize(triangleIndex.size());
1985 
1986  forAll(triangleIndex, i)
1987  {
1988  label triI = triangleIndex[i];
1989  values[i] = fld[triI];
1990  }
1991 
1992 
1993  // Send back results
1994  // ~~~~~~~~~~~~~~~~~
1995 
1996  map.reverseDistribute(info.size(), values);
1997  }
1998 }
1999 
2000 
2003  const pointField& points,
2004  List<volumeType>& volType
2005 ) const
2006 {
2007  FatalErrorIn
2008  (
2009  "distributedTriSurfaceMesh::getVolumeType"
2010  "(const pointField&, List<volumeType>&) const"
2011  ) << "Volume type not supported for distributed surfaces."
2012  << exit(FatalError);
2013 }
2014 
2015 
2016 // Subset the part of surface that is overlapping bb.
2019  const triSurface& s,
2020  const List<treeBoundBox>& bbs,
2021 
2022  labelList& subPointMap,
2023  labelList& subFaceMap
2024 )
2025 {
2026  // Determine what triangles to keep.
2027  boolList includedFace(s.size(), false);
2028 
2029  // Create a slightly larger bounding box.
2030  List<treeBoundBox> bbsX(bbs.size());
2031  const scalar eps = 1.0e-4;
2032  forAll(bbs, i)
2033  {
2034  const point mid = 0.5*(bbs[i].min() + bbs[i].max());
2035  const vector halfSpan = (1.0+eps)*(bbs[i].max() - mid);
2036 
2037  bbsX[i].min() = mid - halfSpan;
2038  bbsX[i].max() = mid + halfSpan;
2039  }
2040 
2041  forAll(s, triI)
2042  {
2043  const labelledTri& f = s[triI];
2044  const point& p0 = s.points()[f[0]];
2045  const point& p1 = s.points()[f[1]];
2046  const point& p2 = s.points()[f[2]];
2047 
2048  if (overlaps(bbsX, p0, p1, p2))
2049  {
2050  includedFace[triI] = true;
2051  }
2052  }
2053 
2054  return subsetMesh(s, includedFace, subPointMap, subFaceMap);
2055 }
2056 
2057 
2060  const List<treeBoundBox>& bbs,
2061  const bool keepNonLocal,
2062  autoPtr<mapDistribute>& faceMap,
2063  autoPtr<mapDistribute>& pointMap
2064 )
2065 {
2066  // Get bbs of all domains
2067  // ~~~~~~~~~~~~~~~~~~~~~~
2068 
2069  {
2071 
2072  switch(distType_)
2073  {
2074  case FOLLOW:
2075  newProcBb[Pstream::myProcNo()].setSize(bbs.size());
2076  forAll(bbs, i)
2077  {
2078  newProcBb[Pstream::myProcNo()][i] = bbs[i];
2079  }
2080  Pstream::gatherList(newProcBb);
2081  Pstream::scatterList(newProcBb);
2082  break;
2083 
2084  case INDEPENDENT:
2085  newProcBb = independentlyDistributedBbs(*this);
2086  break;
2087 
2088  case FROZEN:
2089  return;
2090  break;
2091 
2092  default:
2093  FatalErrorIn("distributedTriSurfaceMesh::distribute(..)")
2094  << "Unsupported distribution type." << exit(FatalError);
2095  break;
2096  }
2097 
2098  //if (debug)
2099  //{
2100  // Info<< "old bb:" << procBb_ << endl << endl;
2101  // Info<< "new bb:" << newProcBb << endl << endl;
2102  // Info<< "Same:" << (newProcBb == procBb_) << endl;
2103  //}
2104 
2105  if (newProcBb == procBb_)
2106  {
2107  return;
2108  }
2109  else
2110  {
2111  procBb_.transfer(newProcBb);
2112  dict_.set("bounds", procBb_[Pstream::myProcNo()]);
2113  }
2114  }
2115 
2116 
2117  // Debug information
2118  if (debug)
2119  {
2120  labelList nTris(Pstream::nProcs());
2121  nTris[Pstream::myProcNo()] = triSurface::size();
2122  Pstream::gatherList(nTris);
2123  Pstream::scatterList(nTris);
2124 
2125  Info<< "distributedTriSurfaceMesh::distribute : before distribution:"
2126  << endl
2127  << "\tproc\ttris" << endl;
2128 
2129  forAll(nTris, procI)
2130  {
2131  Info<< '\t' << procI << '\t' << nTris[procI] << endl;
2132  }
2133  Info<< endl;
2134  }
2135 
2136 
2137  // Use procBbs to determine which faces go where
2138  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2139 
2140  labelListList faceSendMap(Pstream::nProcs());
2141  labelListList pointSendMap(Pstream::nProcs());
2142 
2143  forAll(procBb_, procI)
2144  {
2146  (
2147  *this,
2148  procBb_[procI],
2149  pointSendMap[procI],
2150  faceSendMap[procI]
2151  );
2152 
2153  if (debug)
2154  {
2155  //Pout<< "Overlapping with proc " << procI
2156  // << " faces:" << faceSendMap[procI].size()
2157  // << " points:" << pointSendMap[procI].size() << endl << endl;
2158  }
2159  }
2160 
2161  if (keepNonLocal)
2162  {
2163  // Include in faceSendMap/pointSendMap the triangles that are
2164  // not mapped to any processor so they stay local.
2165 
2166  const triSurface& s = static_cast<const triSurface&>(*this);
2167 
2168  boolList includedFace(s.size(), true);
2169 
2170  forAll(faceSendMap, procI)
2171  {
2172  if (procI != Pstream::myProcNo())
2173  {
2174  forAll(faceSendMap[procI], i)
2175  {
2176  includedFace[faceSendMap[procI][i]] = false;
2177  }
2178  }
2179  }
2180 
2181  // Combine includedFace (all triangles that are not on any neighbour)
2182  // with overlap.
2183 
2184  forAll(faceSendMap[Pstream::myProcNo()], i)
2185  {
2186  includedFace[faceSendMap[Pstream::myProcNo()][i]] = true;
2187  }
2188 
2189  subsetMesh
2190  (
2191  s,
2192  includedFace,
2193  pointSendMap[Pstream::myProcNo()],
2194  faceSendMap[Pstream::myProcNo()]
2195  );
2196  }
2197 
2198 
2199  // Send over how many faces/points I need to receive
2200  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2201 
2202  labelListList faceSendSizes(Pstream::nProcs());
2203  faceSendSizes[Pstream::myProcNo()].setSize(Pstream::nProcs());
2204  forAll(faceSendMap, procI)
2205  {
2206  faceSendSizes[Pstream::myProcNo()][procI] = faceSendMap[procI].size();
2207  }
2208  Pstream::gatherList(faceSendSizes);
2209  Pstream::scatterList(faceSendSizes);
2210 
2211 
2212 
2213  // Exchange surfaces
2214  // ~~~~~~~~~~~~~~~~~
2215 
2216  // Storage for resulting surface
2217  List<labelledTri> allTris;
2219 
2220  labelListList faceConstructMap(Pstream::nProcs());
2221  labelListList pointConstructMap(Pstream::nProcs());
2222 
2223 
2224  // My own surface first
2225  // ~~~~~~~~~~~~~~~~~~~~
2226 
2227  {
2228  labelList pointMap;
2229  triSurface subSurface
2230  (
2231  subsetMesh
2232  (
2233  *this,
2234  faceSendMap[Pstream::myProcNo()],
2235  pointMap
2236  )
2237  );
2238 
2239  allTris = subSurface;
2240  allPoints = subSurface.points();
2241 
2242  faceConstructMap[Pstream::myProcNo()] = identity
2243  (
2244  faceSendMap[Pstream::myProcNo()].size()
2245  );
2246  pointConstructMap[Pstream::myProcNo()] = identity
2247  (
2248  pointSendMap[Pstream::myProcNo()].size()
2249  );
2250  }
2251 
2252 
2253 
2254  // Send all
2255  // ~~~~~~~~
2256 
2257  forAll(faceSendSizes, procI)
2258  {
2259  if (procI != Pstream::myProcNo())
2260  {
2261  if (faceSendSizes[Pstream::myProcNo()][procI] > 0)
2262  {
2263  OPstream str(Pstream::blocking, procI);
2264 
2265  labelList pointMap;
2266  triSurface subSurface
2267  (
2268  subsetMesh
2269  (
2270  *this,
2271  faceSendMap[procI],
2272  pointMap
2273  )
2274  );
2275 
2276  //if (debug)
2277  //{
2278  // Pout<< "Sending to " << procI
2279  // << " faces:" << faceSendMap[procI].size()
2280  // << " points:" << subSurface.points().size() << endl
2281  // << endl;
2282  //}
2283 
2284  str << subSurface;
2285  }
2286  }
2287  }
2288 
2289 
2290  // Receive and merge all
2291  // ~~~~~~~~~~~~~~~~~~~~~
2292 
2293  forAll(faceSendSizes, procI)
2294  {
2295  if (procI != Pstream::myProcNo())
2296  {
2297  if (faceSendSizes[procI][Pstream::myProcNo()] > 0)
2298  {
2299  IPstream str(Pstream::blocking, procI);
2300 
2301  // Receive
2302  triSurface subSurface(str);
2303 
2304  //if (debug)
2305  //{
2306  // Pout<< "Received from " << procI
2307  // << " faces:" << subSurface.size()
2308  // << " points:" << subSurface.points().size() << endl
2309  // << endl;
2310  //}
2311 
2312  // Merge into allSurf
2313  merge
2314  (
2315  mergeDist_,
2316  subSurface,
2317  subSurface.points(),
2318 
2319  allTris,
2320  allPoints,
2321  faceConstructMap[procI],
2322  pointConstructMap[procI]
2323  );
2324 
2325  //if (debug)
2326  //{
2327  // Pout<< "Current merged surface : faces:" << allTris.size()
2328  // << " points:" << allPoints.size() << endl << endl;
2329  //}
2330  }
2331  }
2332  }
2333 
2334 
2335  faceMap.reset
2336  (
2337  new mapDistribute
2338  (
2339  allTris.size(),
2340  faceSendMap.xfer(),
2341  faceConstructMap.xfer()
2342  )
2343  );
2344  pointMap.reset
2345  (
2346  new mapDistribute
2347  (
2348  allPoints.size(),
2349  pointSendMap.xfer(),
2350  pointConstructMap.xfer()
2351  )
2352  );
2353 
2354  // Construct triSurface. Reuse storage.
2355  triSurface::operator=(triSurface(allTris, patches(), allPoints, true));
2356 
2357  clearOut();
2358 
2359  // Set the bounds() value to the boundBox of the undecomposed surface
2361 
2362  reduce(bounds().min(), minOp<point>());
2363  reduce(bounds().max(), maxOp<point>());
2364 
2365  // Regions stays same
2366  // Volume type stays same.
2367 
2368  distributeFields<label>(faceMap());
2369  distributeFields<scalar>(faceMap());
2370  distributeFields<vector>(faceMap());
2371  distributeFields<sphericalTensor>(faceMap());
2372  distributeFields<symmTensor>(faceMap());
2373  distributeFields<tensor>(faceMap());
2374 
2375  if (debug)
2376  {
2377  labelList nTris(Pstream::nProcs());
2378  nTris[Pstream::myProcNo()] = triSurface::size();
2379  Pstream::gatherList(nTris);
2380  Pstream::scatterList(nTris);
2381 
2382  Info<< "distributedTriSurfaceMesh::distribute : after distribution:"
2383  << endl
2384  << "\tproc\ttris" << endl;
2385 
2386  forAll(nTris, procI)
2387  {
2388  Info<< '\t' << procI << '\t' << nTris[procI] << endl;
2389  }
2390  Info<< endl;
2391  }
2392 }
2393 
2394 
2395 //- Write using given format, version and compression
2401 ) const
2402 {
2403  // Make sure dictionary goes to same directory as surface
2404  const_cast<fileName&>(dict_.instance()) = searchableSurface::instance();
2405 
2406  // Copy of triSurfaceMesh::writeObject except for the sorting of
2407  // triangles by region. This is done so we preserve region names,
2408  // even if locally we have zero triangles.
2409  {
2411 
2412  if (!mkDir(fullPath.path()))
2413  {
2414  return false;
2415  }
2416 
2417  // Important: preserve any zero-sized patches
2418  triSurface::write(fullPath, true);
2419 
2420  if (!isFile(fullPath))
2421  {
2422  return false;
2423  }
2424  }
2425 
2426  // Dictionary needs to be written in ascii - binary output not supported.
2427  return dict_.writeObject(IOstream::ASCII, ver, cmp);
2428 }
2429 
2430 
2432 {
2433  boundBox bb;
2434  label nPoints;
2435  PatchTools::calcBounds(static_cast<const triSurface&>(*this), bb, nPoints);
2436  reduce(bb.min(), minOp<point>());
2437  reduce(bb.max(), maxOp<point>());
2438 
2439  os << "Triangles : " << returnReduce(triSurface::size(), sumOp<label>())
2440  << endl
2441  << "Vertices : " << returnReduce(nPoints, sumOp<label>()) << endl
2442  << "Bounding Box : " << bb << endl;
2443 }
2444 
2445 
2446 // ************************************************************************* //
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:380
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:260
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp) const
Write using given format, version and compression.
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
Geometrical domain decomposition.
Definition: geomDecomp.H:47
readOption readOpt() const
Definition: IOobject.H:304
void invertManyToMany(const label len, const UList< InList > &, List< OutList > &)
Invert many-to-many.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
vector point
Point is a vector.
Definition: point.H:41
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
virtual void getField(const List< pointIndexHit > &, labelList &) const
WIP. From a set of hits (points and.
static void distribute(const Pstream::commsTypes commsType, const List< labelPair > &schedule, const label constructSize, const labelListList &subMap, const labelListList &constructMap, List< T > &, const int tag=UPstream::msgType())
Distribute data. Note:schedule only used for Pstream::scheduled.
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:390
Base class of (analytical or triangulated) surface. Encapsulates all the search routines. WIP.
dimensioned< scalar > mag(const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
bool contains(const vector &dir, const point &) const
Contains point (inside or on edge) and moving in direction.
Definition: treeBoundBox.C:402
labelList f(nPoints)
virtual label size() const
Range of local indices that can be returned.
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
virtual bool writeObject(IOstream::streamFormat, IOstream::versionNumber, IOstream::compressionType) const
Write using given format, version and compression.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
void reset(T *=0)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
pointIndexHit findLineAny(const point &start, const point &end) const
Find any intersection of line between start and end.
Xfer< List< T > > xfer()
Transfer contents to the Xfer container.
Definition: ListI.H:90
Class containing processor-to-processor mapping information.
Triangulated surface description with patch information.
Definition: triSurface.H:57
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
const fileName & local() const
Definition: IOobject.H:347
A class for handling words, derived from string.
Definition: word.H:59
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:75
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:310
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
void writeStats(Ostream &os) const
Print some stats. Parallel aware version of.
Output inter-processor communications stream.
Definition: OPstream.H:50
const fileName & instance() const
Definition: IOobject.H:337
triSurface()
Construct null.
Definition: triSurface.C:611
const boundBox & bounds() const
Return const reference to boundBox.
virtual void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit > > &) const
Get all intersections in order from start to end.
void clearOut()
Clear storage.
static const Vector max
Definition: Vector.H:82
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
messageStream Info
Accumulating histogram of values. Specified bin resolution automatic generation of bins...
Definition: distribution.H:51
bool hit() const
Is there a hit.
fileName path() const
Return directory path name (part before last /)
Definition: fileName.C:287
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
bool & registerObject()
Register object created from this IOobject with registry if true.
Definition: IOobject.H:290
label nPoints() const
Return number of points supporting patch faces.
Namespace for OpenFOAM.
const Point & hitPoint() const
Return hit point.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Pair< point > segment
void reverseDistribute(const label constructSize, List< T > &, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Reverse distribute data using default commsType.
label n
T & operator[](const label)
Return element of UList.
Definition: UListI.H:163
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:386
void clear()
Clear the list, i.e. set size to zero.
Definition: List.C:379
Triangle with additional region number.
Definition: labelledTri.H:49
label region() const
Return region label.
Definition: labelledTriI.H:68
void setSize(const label)
Reset size of List.
Definition: List.C:318
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const geometricSurfacePatchList & patches() const
Definition: triSurface.H:306
virtual void distribute(const List< treeBoundBox > &, const bool keepNonLocal, autoPtr< mapDistribute > &faceMap, autoPtr< mapDistribute > &pointMap)
Set bounds of surface. Bounds currently set as list of.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
const indexedOctree< treeDataTriSurface > & tree() const
Demand driven construction of the octree.
streamFormat
Enumeration for the format of data in the stream.
Definition: IOstream.H:86
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define forAll(list, i)
Definition: UList.H:421
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
writeOption writeOpt() const
Definition: IOobject.H:314
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:258
fileName objectPath() const
Return complete path + object name.
Definition: IOobject.H:363
virtual void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Return any intersection on segment from start to end.
Macros for easy insertion into run-time selection tables.
void operator=(const triSurface &)
Definition: triSurface.C:1125
static triSurface overlappingSurface(const triSurface &, const List< treeBoundBox > &, labelList &subPointMap, labelList &subFaceMap)
Subset the part of surface that is overlapping bounds.
virtual tmp< pointField > points() const
Get the points that define the surface.
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
static void calcBounds(const PrimitivePatch< Face, FaceList, PointField, PointType > &p, boundBox &bb, label &nPoints)
static autoPtr< decompositionMethod > New(const dictionary &decompositionDict)
Return a reference to the selected decomposition method.
virtual labelList decompose(const pointField &points, const scalarField &pointWeights)=0
Return for every coordinate the wanted processor number.
const word & name() const
Return name.
Definition: IOobject.H:260
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:404
void setIndex(const label index)
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject( name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE ))
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:55
Non-pointer based hierarchical recursive searching.
Definition: treeDataEdge.H:47
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
compressionType
Enumeration for the format of data in the stream.
Definition: IOstream.H:193
virtual void getField(const List< pointIndexHit > &, labelList &) const
WIP. From a set of hits (points and.
static const Vector min
Definition: Vector.H:83
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:866
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
int system(const std::string &command)
Execute the specified command.
Definition: POSIX.C:1162
Determine correspondence between points. See below.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
error FatalError
Input inter-processor communications stream.
Definition: IPstream.H:50
label toLocal(const label i) const
From global to local on current processor.
Definition: globalIndexI.H:117
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
A class for handling file names.
Definition: fileName.H:69
bool matchPoints(const UList< point > &pts0, const UList< point > &pts1, const UList< scalar > &matchDistance, const bool verbose, labelList &from0To1, const point &origin=point::zero)
Determine correspondence between pointFields. Gets passed.
Definition: matchPoints.C:33
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: FixedListI.H:115
label size() const
Return the number of elements in the UList.
word findInstance(const fileName &dir, const word &name=word::null, const IOobject::readOption rOpt=IOobject::MUST_READ, const word &stopInstance=word::null) const
Return the location of "dir" containing the file "name".
Definition: findInstance.C:38
static fileCheckTypes fileModificationChecking
Definition: regIOobject.H:128
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
const Time & time() const
Return time.
Definition: IOobject.C:251
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
bool read(const char *, int32_t &)
Definition: int32IO.C:87
static const NamedEnum< distributionType, 3 > distributionTypeNames_
static bool intersectBb(const point &p0, const point &p1, const point &p2, const treeBoundBox &cubeBb)
Does triangle intersect bounding box.
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:97
bool isFile(const fileName &, const bool checkGzip=true)
Does the name exist as a FILE in the file system?
Definition: POSIX.C:623
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints,-1);const cellModel &hex=*(cellModeller::lookup("hex"));labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){pointMap[i]=i;}for(label i=0;i< nPoints;i++){if(f[i] > 0.0){hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells]=cellShape(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei){if(edges[ei].mag(points)< SMALL){label start=pointMap[edges[ei].start()];while(start!=pointMap[start]){start=pointMap[start];}label end=pointMap[edges[ei].end()];while(end!=pointMap[end]){end=pointMap[end];}label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;}}cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){cellShape &cs=cellShapes[celli];forAll(cs, i){cs[i]=pointMap[cs[i]];}cs.collapse();}label bcIDs[11]={-1, 0, 2, 4,-1, 5,-1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={&wallPolyPatch::typeName,&wallPolyPatch::typeName,&wallPolyPatch::typeName,&wallPolyPatch::typeName,&symmetryPolyPatch::typeName,&wedgePolyPatch::typeName,&polyPatch::typeName,&polyPatch::typeName,&polyPatch::typeName,&polyPatch::typeName,&symmetryPolyPatch::typeName,&oldCyclicPolyPatch::typeName};enum patchTypeNames{PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={"piston","valve","liner","cylinderHead","axis","wedge","inflow","outflow","presin","presout","symmetryPlane","cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
const globalIndex & globalTris() const
Triangle indexing (demand driven)
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:420
virtual void getRegion(const List< pointIndexHit > &, labelList &region) const
From a set of points and indices get the region.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Version number type.
Definition: IOstream.H:96
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:245
virtual void findNearest(const pointField &sample, const scalarField &nearestDistSqr, List< pointIndexHit > &) const
bool overlaps(const boundBox &) const
Overlaps/touches boundingBox?
Definition: boundBoxI.H:120
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
IOoject and searching on triSurface.
const Field< PointType > & points() const
Return reference to global points.
virtual void getVolumeType(const pointField &, List< volumeType > &) const
Determine type (inside/outside/mixed) for point. unknown if.
label whichProcID(const label i) const
Which processor does global come from? Binary search.
Definition: globalIndexI.H:123
static const word null
An empty word.
Definition: word.H:77
defineTypeNameAndDebug(combustionModel, 0)