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