fvMeshDistribute.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "fvMeshDistribute.H"
28 #include "fvMeshAdder.H"
29 #include "faceCoupleInfo.H"
30 #include "processorFvPatchField.H"
31 #include "processorFvsPatchField.H"
34 #include "polyTopoChange.H"
35 #include "removeCells.H"
36 #include "polyModifyFace.H"
37 #include "polyRemovePoint.H"
38 #include "mapDistributePolyMesh.H"
39 #include "surfaceFields.H"
40 #include "pointFields.H"
41 #include "syncTools.H"
42 #include "CompactListList.H"
43 #include "fvMeshTools.H"
44 #include "ListOps.H"
45 #include "globalIndex.H"
46 
47 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51  defineTypeNameAndDebug(fvMeshDistribute, 0);
52 
53 //- Less function class that can be used for sorting processor patches
55 {
56  const labelList& nbrProc_;
57  const labelList& referPatchID_;
58 
59 public:
60 
61  lessProcPatches( const labelList& nbrProc, const labelList& referPatchID)
62  :
63  nbrProc_(nbrProc),
64  referPatchID_(referPatchID)
65  {}
66 
67  bool operator()(const label a, const label b)
68  {
69  if (nbrProc_[a] < nbrProc_[b])
70  {
71  return true;
72  }
73  else if (nbrProc_[a] > nbrProc_[b])
74  {
75  return false;
76  }
77  else
78  {
79  // Equal neighbour processor
80  return referPatchID_[a] < referPatchID_[b];
81  }
82  }
83 };
84 
85 }
86 
87 
88 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
89 
90 void Foam::fvMeshDistribute::inplaceRenumberWithFlip
91 (
92  const labelUList& oldToNew,
93  const bool oldToNewHasFlip,
94  const bool lstHasFlip,
95  labelUList& lst
96 )
97 {
98  if (!lstHasFlip && !oldToNewHasFlip)
99  {
100  Foam::inplaceRenumber(oldToNew, lst);
101  }
102  else
103  {
104  // Either input data or map encodes sign so result encodes sign
105 
106  forAll(lst, elemI)
107  {
108  // Extract old value and sign
109  label val = lst[elemI];
110  label sign = 1;
111  if (lstHasFlip)
112  {
113  if (val > 0)
114  {
115  val = val-1;
116  }
117  else if (val < 0)
118  {
119  val = -val-1;
120  sign = -1;
121  }
122  else
123  {
125  << "Problem : zero value " << val
126  << " at index " << elemI << " out of " << lst.size()
127  << " list with flip bit" << exit(FatalError);
128  }
129  }
130 
131 
132  // Lookup new value and possibly change sign
133  label newVal = oldToNew[val];
134 
135  if (oldToNewHasFlip)
136  {
137  if (newVal > 0)
138  {
139  newVal = newVal-1;
140  }
141  else if (newVal < 0)
142  {
143  newVal = -newVal-1;
144  sign = -sign;
145  }
146  else
147  {
149  << "Problem : zero value " << newVal
150  << " at index " << elemI << " out of "
151  << oldToNew.size()
152  << " list with flip bit" << exit(FatalError);
153  }
154  }
155 
156 
157  // Encode new value and sign
158  lst[elemI] = sign*(newVal+1);
159  }
160  }
161 }
162 
163 
164 Foam::labelList Foam::fvMeshDistribute::select
165 (
166  const bool selectEqual,
167  const labelList& values,
168  const label value
169 )
170 {
171  label n = 0;
172 
173  forAll(values, i)
174  {
175  if (selectEqual == (values[i] == value))
176  {
177  n++;
178  }
179  }
180 
181  labelList indices(n);
182  n = 0;
183 
184  forAll(values, i)
185  {
186  if (selectEqual == (values[i] == value))
187  {
188  indices[n++] = i;
189  }
190  }
191  return indices;
192 }
193 
194 
195 // Check all procs have same names and in exactly same order.
196 void Foam::fvMeshDistribute::checkEqualWordList
197 (
198  const string& msg,
199  const wordList& lst
200 )
201 {
202  List<wordList> allNames(Pstream::nProcs());
203  allNames[Pstream::myProcNo()] = lst;
204  Pstream::gatherList(allNames);
205  Pstream::scatterList(allNames);
206 
207  for (label proci = 1; proci < Pstream::nProcs(); proci++)
208  {
209  if (allNames[proci] != allNames[0])
210  {
212  << "When checking for equal " << msg.c_str() << " :" << endl
213  << "processor0 has:" << allNames[0] << endl
214  << "processor" << proci << " has:" << allNames[proci] << endl
215  << msg.c_str() << " need to be synchronised on all processors."
216  << exit(FatalError);
217  }
218  }
219 }
220 
221 
222 Foam::wordList Foam::fvMeshDistribute::mergeWordList(const wordList& procNames)
223 {
224  List<wordList> allNames(Pstream::nProcs());
225  allNames[Pstream::myProcNo()] = procNames;
226  Pstream::gatherList(allNames);
227  Pstream::scatterList(allNames);
228 
229  HashSet<word> mergedNames;
230  forAll(allNames, proci)
231  {
232  forAll(allNames[proci], i)
233  {
234  mergedNames.insert(allNames[proci][i]);
235  }
236  }
237  return mergedNames.toc();
238 }
239 
240 
241 // Print some info on mesh.
243 {
244  Pout<< "Primitives:" << nl
245  << " points :" << mesh.nPoints() << nl
246  << " bb :" << boundBox(mesh.points(), false) << nl
247  << " internalFaces:" << mesh.nInternalFaces() << nl
248  << " faces :" << mesh.nFaces() << nl
249  << " cells :" << mesh.nCells() << nl;
250 
251  const fvBoundaryMesh& patches = mesh.boundary();
252 
253  Pout<< "Patches:" << endl;
254  forAll(patches, patchi)
255  {
256  const polyPatch& pp = patches[patchi].patch();
257 
258  Pout<< " " << patchi << " name:" << pp.name()
259  << " size:" << pp.size()
260  << " start:" << pp.start()
261  << " type:" << pp.type()
262  << endl;
263  }
264 
265  if (mesh.pointZones().size())
266  {
267  Pout<< "PointZones:" << endl;
268  forAll(mesh.pointZones(), zoneI)
269  {
270  const pointZone& pz = mesh.pointZones()[zoneI];
271  Pout<< " " << zoneI << " name:" << pz.name()
272  << " size:" << pz.size()
273  << endl;
274  }
275  }
276  if (mesh.faceZones().size())
277  {
278  Pout<< "FaceZones:" << endl;
279  forAll(mesh.faceZones(), zoneI)
280  {
281  const faceZone& fz = mesh.faceZones()[zoneI];
282  Pout<< " " << zoneI << " name:" << fz.name()
283  << " size:" << fz.size()
284  << endl;
285  }
286  }
287  if (mesh.cellZones().size())
288  {
289  Pout<< "CellZones:" << endl;
290  forAll(mesh.cellZones(), zoneI)
291  {
292  const cellZone& cz = mesh.cellZones()[zoneI];
293  Pout<< " " << zoneI << " name:" << cz.name()
294  << " size:" << cz.size()
295  << endl;
296  }
297  }
298 }
299 
300 
302 (
303  const primitiveMesh& mesh,
304  const labelList& sourceFace,
305  const labelList& sourceProc,
306  const labelList& sourcePatch,
307  const labelList& sourceNewNbrProc
308 )
309 {
310  Pout<< nl
311  << "Current coupling info:"
312  << endl;
313 
314  forAll(sourceFace, bFacei)
315  {
316  label meshFacei = mesh.nInternalFaces() + bFacei;
317 
318  Pout<< " meshFace:" << meshFacei
319  << " fc:" << mesh.faceCentres()[meshFacei]
320  << " connects to proc:" << sourceProc[bFacei]
321  << "/face:" << sourceFace[bFacei]
322  << " which will move to proc:" << sourceNewNbrProc[bFacei]
323  << endl;
324  }
325 }
326 
327 
328 // Finds (non-empty) patch that exposed internal and proc faces can be put into.
329 Foam::label Foam::fvMeshDistribute::findNonEmptyPatch() const
330 {
331  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
332 
333  label nonEmptyPatchi = -1;
334 
335  forAllReverse(patches, patchi)
336  {
337  const polyPatch& pp = patches[patchi];
338 
339  if (!isA<emptyPolyPatch>(pp) && !pp.coupled())
340  {
341  nonEmptyPatchi = patchi;
342  break;
343  }
344  }
345 
346  if (nonEmptyPatchi == -1)
347  {
349  << "Cannot find a patch which is neither of type empty nor"
350  << " coupled in patches " << patches.names() << endl
351  << "There has to be at least one such patch for"
352  << " distribution to work" << abort(FatalError);
353  }
354 
355  if (debug)
356  {
357  Pout<< "findNonEmptyPatch : using patch " << nonEmptyPatchi
358  << " name:" << patches[nonEmptyPatchi].name()
359  << " type:" << patches[nonEmptyPatchi].type()
360  << " to put exposed faces into." << endl;
361  }
362 
363 
364  // Do additional test for processor patches intermingled with non-proc
365  // patches.
366  label procPatchi = -1;
367 
368  forAll(patches, patchi)
369  {
370  if (isA<processorPolyPatch>(patches[patchi]))
371  {
372  procPatchi = patchi;
373  }
374  else if (procPatchi != -1)
375  {
377  << "Processor patches should be at end of patch list."
378  << endl
379  << "Have processor patch " << procPatchi
380  << " followed by non-processor patch " << patchi
381  << " in patches " << patches.names()
382  << abort(FatalError);
383  }
384  }
385 
386  return nonEmptyPatchi;
387 }
388 
389 
390 // Delete all processor patches. Move any processor faces into the last
391 // non-processor patch.
392 Foam::autoPtr<Foam::mapPolyMesh> Foam::fvMeshDistribute::deleteProcPatches
393 (
394  const label destinationPatch
395 )
396 {
397  // New patchID per boundary faces to be repatched. Is -1 (no change)
398  // or new patchID
399  labelList newPatchID(mesh_.nFaces() - mesh_.nInternalFaces(), -1);
400 
401  label nProcPatches = 0;
402 
403  forAll(mesh_.boundaryMesh(), patchi)
404  {
405  const polyPatch& pp = mesh_.boundaryMesh()[patchi];
406 
407  if (isA<processorPolyPatch>(pp))
408  {
409  if (debug)
410  {
411  Pout<< "Moving all faces of patch " << pp.name()
412  << " into patch " << destinationPatch
413  << endl;
414  }
415 
416  label offset = pp.start() - mesh_.nInternalFaces();
417 
418  forAll(pp, i)
419  {
420  newPatchID[offset+i] = destinationPatch;
421  }
422 
423  nProcPatches++;
424  }
425  }
426 
427  // Note: order of boundary faces been kept the same since the
428  // destinationPatch is at the end and we have visited the patches in
429  // incremental order.
430  labelListList dummyFaceMaps;
431  autoPtr<mapPolyMesh> map = repatch(newPatchID, dummyFaceMaps);
432 
433 
434  // Delete (now empty) processor patches.
435  {
436  labelList oldToNew(identity(mesh_.boundaryMesh().size()));
437  label newI = 0;
438  // Non processor patches first
439  forAll(mesh_.boundaryMesh(), patchi)
440  {
441  if (!isA<processorPolyPatch>(mesh_.boundaryMesh()[patchi]))
442  {
443  oldToNew[patchi] = newI++;
444  }
445  }
446  label nNonProcPatches = newI;
447 
448  // Processor patches as last
449  forAll(mesh_.boundaryMesh(), patchi)
450  {
451  if (isA<processorPolyPatch>(mesh_.boundaryMesh()[patchi]))
452  {
453  oldToNew[patchi] = newI++;
454  }
455  }
456  fvMeshTools::reorderPatches(mesh_, oldToNew, nNonProcPatches, false);
457  }
458 
459  return map;
460 }
461 
462 
463 // Repatch the mesh.
464 Foam::autoPtr<Foam::mapPolyMesh> Foam::fvMeshDistribute::repatch
465 (
466  const labelList& newPatchID, // per boundary face -1 or new patchID
467  labelListList& constructFaceMap
468 )
469 {
470  polyTopoChange meshMod(mesh_);
471 
472  forAll(newPatchID, bFacei)
473  {
474  if (newPatchID[bFacei] != -1)
475  {
476  label facei = mesh_.nInternalFaces() + bFacei;
477 
478  label zoneID = mesh_.faceZones().whichZone(facei);
479  bool zoneFlip = false;
480 
481  if (zoneID >= 0)
482  {
483  const faceZone& fZone = mesh_.faceZones()[zoneID];
484  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
485  }
486 
487  meshMod.setAction
488  (
490  (
491  mesh_.faces()[facei], // modified face
492  facei, // label of face
493  mesh_.faceOwner()[facei], // owner
494  -1, // neighbour
495  false, // face flip
496  newPatchID[bFacei], // patch for face
497  false, // remove from zone
498  zoneID, // zone for face
499  zoneFlip // face flip in zone
500  )
501  );
502  }
503  }
504 
505 
506  // Do mapping of fields from one patchField to the other ourselves since
507  // is currently not supported by updateMesh.
508 
509  // Store boundary fields (we only do this for surfaceFields)
511  saveBoundaryFields<scalar, surfaceMesh>(sFlds);
513  saveBoundaryFields<vector, surfaceMesh>(vFlds);
515  saveBoundaryFields<sphericalTensor, surfaceMesh>(sptFlds);
517  saveBoundaryFields<symmTensor, surfaceMesh>(sytFlds);
519  saveBoundaryFields<tensor, surfaceMesh>(tFlds);
520 
521  // Change the mesh (no inflation). Note: parallel comms allowed.
522  //
523  // NOTE: there is one very particular problem with this ordering.
524  // We first create the processor patches and use these to merge out
525  // shared points (see mergeSharedPoints below). So temporarily points
526  // and edges do not match!
527 
528  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh_, false, true);
529 
530  // Update fields. No inflation, parallel sync.
531  mesh_.updateMesh(map);
532 
533  // Map patch fields using stored boundary fields. Note: assumes order
534  // of fields has not changed in object registry!
535  mapBoundaryFields<scalar, surfaceMesh>(map, sFlds);
536  mapBoundaryFields<vector, surfaceMesh>(map, vFlds);
537  mapBoundaryFields<sphericalTensor, surfaceMesh>(map, sptFlds);
538  mapBoundaryFields<symmTensor, surfaceMesh>(map, sytFlds);
539  mapBoundaryFields<tensor, surfaceMesh>(map, tFlds);
540 
541 
542  // Move mesh (since morphing does not do this)
543  if (map().hasMotionPoints())
544  {
545  mesh_.movePoints(map().preMotionPoints());
546  }
547 
548  // Adapt constructMaps.
549 
550  if (debug)
551  {
552  label index = findIndex(map().reverseFaceMap(), -1);
553 
554  if (index != -1)
555  {
557  << "reverseFaceMap contains -1 at index:"
558  << index << endl
559  << "This means that the repatch operation was not just"
560  << " a shuffle?" << abort(FatalError);
561  }
562  }
563 
564  forAll(constructFaceMap, proci)
565  {
566  inplaceRenumberWithFlip
567  (
568  map().reverseFaceMap(),
569  false,
570  true,
571  constructFaceMap[proci]
572  );
573  }
574 
575 
576  return map;
577 }
578 
579 
580 // Detect shared points. Need processor patches to be present.
581 // Background: when adding bits of mesh one can get points which
582 // share the same position but are only detectable to be topologically
583 // the same point when doing parallel analysis. This routine will
584 // merge those points.
585 Foam::autoPtr<Foam::mapPolyMesh> Foam::fvMeshDistribute::mergeSharedPoints
586 (
587  const labelList& pointToGlobalMaster,
588  labelListList& constructPointMap
589 )
590 {
591  // Find out which sets of points get merged and create a map from
592  // mesh point to unique point.
593 
594  label nShared = 0;
595  forAll(pointToGlobalMaster, pointi)
596  {
597  if (pointToGlobalMaster[pointi] != -1)
598  {
599  nShared++;
600  }
601  }
602 
603  Map<label> globalMasterToLocalMaster(2*nShared);
604  Map<label> pointToMaster(2*nShared);
605 
606  forAll(pointToGlobalMaster, pointi)
607  {
608  label globali = pointToGlobalMaster[pointi];
609  if (globali != -1)
610  {
611  Map<label>::const_iterator iter = globalMasterToLocalMaster.find
612  (
613  globali
614  );
615 
616  if (iter == globalMasterToLocalMaster.end())
617  {
618  // Found first point. Designate as master
619  globalMasterToLocalMaster.insert(globali, pointi);
620  pointToMaster.insert(pointi, pointi);
621  }
622  else
623  {
624  pointToMaster.insert(pointi, iter());
625  }
626  }
627  }
628 
629  if (returnReduce(pointToMaster.size(), sumOp<label>()) == 0)
630  {
631  return autoPtr<mapPolyMesh>(nullptr);
632  }
633 
634 
635  polyTopoChange meshMod(mesh_);
636 
637  fvMeshAdder::mergePoints(mesh_, pointToMaster, meshMod);
638 
639  // Change the mesh (no inflation). Note: parallel comms allowed.
640  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh_, false, true);
641 
642  // Update fields. No inflation, parallel sync.
643  mesh_.updateMesh(map);
644 
645  // Adapt constructMaps for merged points.
646  forAll(constructPointMap, proci)
647  {
648  labelList& constructMap = constructPointMap[proci];
649 
650  forAll(constructMap, i)
651  {
652  label oldPointi = constructMap[i];
653 
654  label newPointi = map().reversePointMap()[oldPointi];
655 
656  if (newPointi < -1)
657  {
658  constructMap[i] = -newPointi-2;
659  }
660  else if (newPointi >= 0)
661  {
662  constructMap[i] = newPointi;
663  }
664  else
665  {
667  << "Problem. oldPointi:" << oldPointi
668  << " newPointi:" << newPointi << abort(FatalError);
669  }
670  }
671  }
672  return map;
673 }
674 
675 
676 void Foam::fvMeshDistribute::getCouplingData
677 (
678  const labelList& distribution,
679  labelList& sourceFace,
680  labelList& sourceProc,
681  labelList& sourcePatch,
682  labelList& sourceNewNbrProc,
683  labelList& sourcePointMaster
684 ) const
685 {
686  // Construct the coupling information for all (boundary) faces and
687  // points
688 
689  label nBnd = mesh_.nFaces() - mesh_.nInternalFaces();
690  sourceFace.setSize(nBnd);
691  sourceProc.setSize(nBnd);
692  sourcePatch.setSize(nBnd);
693  sourceNewNbrProc.setSize(nBnd);
694 
695  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
696 
697  // Get neighbouring meshFace labels and new processor of coupled boundaries.
698  labelList nbrFaces(nBnd, -1);
699  labelList nbrNewNbrProc(nBnd, -1);
700 
701  forAll(patches, patchi)
702  {
703  const polyPatch& pp = patches[patchi];
704 
705  if (pp.coupled())
706  {
707  label offset = pp.start() - mesh_.nInternalFaces();
708 
709  // Mesh labels of faces on this side
710  forAll(pp, i)
711  {
712  label bndI = offset + i;
713  nbrFaces[bndI] = pp.start()+i;
714  }
715 
716  // Which processor they will end up on
717  SubList<label>(nbrNewNbrProc, pp.size(), offset) =
718  UIndirectList<label>(distribution, pp.faceCells())();
719  }
720  }
721 
722 
723  // Exchange the boundary data
724  syncTools::swapBoundaryFaceList(mesh_, nbrFaces);
725  syncTools::swapBoundaryFaceList(mesh_, nbrNewNbrProc);
726 
727 
728  forAll(patches, patchi)
729  {
730  const polyPatch& pp = patches[patchi];
731  label offset = pp.start() - mesh_.nInternalFaces();
732 
733  if (isA<processorPolyPatch>(pp))
734  {
735  const processorPolyPatch& procPatch =
736  refCast<const processorPolyPatch>(pp);
737 
738  // Check which of the two faces we store.
739 
740  if (procPatch.owner())
741  {
742  // Use my local face labels
743  forAll(pp, i)
744  {
745  label bndI = offset + i;
746  sourceFace[bndI] = pp.start()+i;
747  sourceProc[bndI] = Pstream::myProcNo();
748  sourceNewNbrProc[bndI] = nbrNewNbrProc[bndI];
749  }
750  }
751  else
752  {
753  // Use my neighbours face labels
754  forAll(pp, i)
755  {
756  label bndI = offset + i;
757  sourceFace[bndI] = nbrFaces[bndI];
758  sourceProc[bndI] = procPatch.neighbProcNo();
759  sourceNewNbrProc[bndI] = nbrNewNbrProc[bndI];
760  }
761  }
762 
763 
764  label patchi = -1;
765  if (isA<processorCyclicPolyPatch>(pp))
766  {
767  patchi = refCast<const processorCyclicPolyPatch>
768  (
769  pp
770  ).referPatchID();
771  }
772 
773  forAll(pp, i)
774  {
775  label bndI = offset + i;
776  sourcePatch[bndI] = patchi;
777  }
778  }
779  else if (isA<cyclicPolyPatch>(pp))
780  {
781  const cyclicPolyPatch& cpp = refCast<const cyclicPolyPatch>(pp);
782 
783  if (cpp.owner())
784  {
785  forAll(pp, i)
786  {
787  label bndI = offset + i;
788  sourceFace[bndI] = pp.start()+i;
789  sourceProc[bndI] = Pstream::myProcNo();
790  sourcePatch[bndI] = patchi;
791  sourceNewNbrProc[bndI] = nbrNewNbrProc[bndI];
792  }
793  }
794  else
795  {
796  forAll(pp, i)
797  {
798  label bndI = offset + i;
799  sourceFace[bndI] = nbrFaces[bndI];
800  sourceProc[bndI] = Pstream::myProcNo();
801  sourcePatch[bndI] = patchi;
802  sourceNewNbrProc[bndI] = nbrNewNbrProc[bndI];
803  }
804  }
805  }
806  else
807  {
808  // Normal (physical) boundary
809  forAll(pp, i)
810  {
811  label bndI = offset + i;
812  sourceFace[bndI] = -1;
813  sourceProc[bndI] = -1;
814  sourcePatch[bndI] = patchi;
815  sourceNewNbrProc[bndI] = -1;
816  }
817  }
818  }
819 
820 
821  // Collect coupled (collocated) points
822  sourcePointMaster.setSize(mesh_.nPoints());
823  sourcePointMaster = -1;
824  {
825  // Assign global master point
826  const globalIndex globalPoints(mesh_.nPoints());
827 
828  const globalMeshData& gmd = mesh_.globalData();
829  const indirectPrimitivePatch& cpp = gmd.coupledPatch();
830  const labelList& meshPoints = cpp.meshPoints();
831  const mapDistribute& slavesMap = gmd.globalCoPointSlavesMap();
832  const labelListList& slaves = gmd.globalCoPointSlaves();
833 
834  labelList elems(slavesMap.constructSize(), -1);
835  forAll(meshPoints, pointi)
836  {
837  const labelList& slots = slaves[pointi];
838 
839  if (slots.size())
840  {
841  // pointi is a master. Assign a unique label.
842 
843  label globalPointi = globalPoints.toGlobal(meshPoints[pointi]);
844  elems[pointi] = globalPointi;
845  forAll(slots, i)
846  {
847  label sloti = slots[i];
848  if (sloti >= meshPoints.size())
849  {
850  // Filter out local collocated points. We don't want
851  // to merge these
852  elems[slots[i]] = globalPointi;
853  }
854  }
855  }
856  }
857 
858  // Push slave-slot data back to slaves
859  slavesMap.reverseDistribute(elems.size(), elems, false);
860 
861  // Extract back onto mesh
862  forAll(meshPoints, pointi)
863  {
864  sourcePointMaster[meshPoints[pointi]] = elems[pointi];
865  }
866  }
867 }
868 
869 
870 // Subset the neighbourCell/neighbourProc fields
871 void Foam::fvMeshDistribute::subsetCouplingData
872 (
873  const fvMesh& mesh,
874  const labelList& pointMap,
875  const labelList& faceMap,
876  const labelList& cellMap,
877 
878  const labelList& oldDistribution,
879  const labelList& oldFaceOwner,
880  const labelList& oldFaceNeighbour,
881  const label oldInternalFaces,
882 
883  const labelList& sourceFace,
884  const labelList& sourceProc,
885  const labelList& sourcePatch,
886  const labelList& sourceNewNbrProc,
887  const labelList& sourcePointMaster,
888 
889  labelList& subFace,
890  labelList& subProc,
891  labelList& subPatch,
892  labelList& subNewNbrProc,
893  labelList& subPointMaster
894 )
895 {
896  subFace.setSize(mesh.nFaces() - mesh.nInternalFaces());
897  subProc.setSize(mesh.nFaces() - mesh.nInternalFaces());
898  subPatch.setSize(mesh.nFaces() - mesh.nInternalFaces());
899  subNewNbrProc.setSize(mesh.nFaces() - mesh.nInternalFaces());
900 
901  forAll(subFace, newBFacei)
902  {
903  label newFacei = newBFacei + mesh.nInternalFaces();
904 
905  label oldFacei = faceMap[newFacei];
906 
907  // Was oldFacei internal face? If so which side did we get.
908  if (oldFacei < oldInternalFaces)
909  {
910  subFace[newBFacei] = oldFacei;
911  subProc[newBFacei] = Pstream::myProcNo();
912  subPatch[newBFacei] = -1;
913 
914  label oldOwn = oldFaceOwner[oldFacei];
915  label oldNei = oldFaceNeighbour[oldFacei];
916 
917  if (oldOwn == cellMap[mesh.faceOwner()[newFacei]])
918  {
919  // We kept the owner side. Where does the neighbour move to?
920  subNewNbrProc[newBFacei] = oldDistribution[oldNei];
921  }
922  else
923  {
924  // We kept the neighbour side.
925  subNewNbrProc[newBFacei] = oldDistribution[oldOwn];
926  }
927  }
928  else
929  {
930  // Was boundary face. Take over boundary information
931  label oldBFacei = oldFacei - oldInternalFaces;
932 
933  subFace[newBFacei] = sourceFace[oldBFacei];
934  subProc[newBFacei] = sourceProc[oldBFacei];
935  subPatch[newBFacei] = sourcePatch[oldBFacei];
936  subNewNbrProc[newBFacei] = sourceNewNbrProc[oldBFacei];
937  }
938  }
939 
940 
941  subPointMaster = UIndirectList<label>(sourcePointMaster, pointMap);
942 }
943 
944 
945 // Find cells on mesh whose faceID/procID match the neighbour cell/proc of
946 // domainMesh. Store the matching face.
947 void Foam::fvMeshDistribute::findCouples
948 (
949  const primitiveMesh& mesh,
950  const labelList& sourceFace,
951  const labelList& sourceProc,
952  const labelList& sourcePatch,
953 
954  const label domain,
955  const primitiveMesh& domainMesh,
956  const labelList& domainFace,
957  const labelList& domainProc,
958  const labelList& domainPatch,
959 
960  labelList& masterCoupledFaces,
961  labelList& slaveCoupledFaces
962 )
963 {
964  // Store domain neighbour as map so we can easily look for pair
965  // with same face+proc.
967 
968  forAll(domainProc, bFacei)
969  {
970  if (domainProc[bFacei] != -1 && domainPatch[bFacei] == -1)
971  {
972  map.insert
973  (
974  labelPair(domainFace[bFacei], domainProc[bFacei]),
975  bFacei
976  );
977  }
978  }
979 
980 
981  // Try to match mesh data.
982 
983  masterCoupledFaces.setSize(domainFace.size());
984  slaveCoupledFaces.setSize(domainFace.size());
985  label coupledI = 0;
986 
987  forAll(sourceFace, bFacei)
988  {
989  if (sourceProc[bFacei] != -1 && sourcePatch[bFacei] == -1)
990  {
991  labelPair myData(sourceFace[bFacei], sourceProc[bFacei]);
992 
994  iter = map.find(myData);
995 
996  if (iter != map.end())
997  {
998  label nbrBFacei = iter();
999 
1000  masterCoupledFaces[coupledI] = mesh.nInternalFaces() + bFacei;
1001  slaveCoupledFaces[coupledI] =
1002  domainMesh.nInternalFaces()
1003  + nbrBFacei;
1004 
1005  coupledI++;
1006  }
1007  }
1008  }
1009 
1010  masterCoupledFaces.setSize(coupledI);
1011  slaveCoupledFaces.setSize(coupledI);
1012 
1013  if (debug)
1014  {
1015  Pout<< "findCouples : found " << coupledI
1016  << " faces that will be stitched" << nl << endl;
1017  }
1018 }
1019 
1020 
1021 // Map data on boundary faces to new mesh (resulting from adding two meshes)
1022 Foam::labelList Foam::fvMeshDistribute::mapBoundaryData
1023 (
1024  const primitiveMesh& mesh, // mesh after adding
1025  const mapAddedPolyMesh& map,
1026  const labelList& boundaryData0, // on mesh before adding
1027  const label nInternalFaces1,
1028  const labelList& boundaryData1 // on added mesh
1029 )
1030 {
1031  labelList newBoundaryData(mesh.nFaces() - mesh.nInternalFaces());
1032 
1033  forAll(boundaryData0, oldBFacei)
1034  {
1035  label newFacei = map.oldFaceMap()[oldBFacei + map.nOldInternalFaces()];
1036 
1037  // Face still exists (is necessary?) and still boundary face
1038  if (newFacei >= 0 && newFacei >= mesh.nInternalFaces())
1039  {
1040  newBoundaryData[newFacei - mesh.nInternalFaces()] =
1041  boundaryData0[oldBFacei];
1042  }
1043  }
1044 
1045  forAll(boundaryData1, addedBFacei)
1046  {
1047  label newFacei = map.addedFaceMap()[addedBFacei + nInternalFaces1];
1048 
1049  if (newFacei >= 0 && newFacei >= mesh.nInternalFaces())
1050  {
1051  newBoundaryData[newFacei - mesh.nInternalFaces()] =
1052  boundaryData1[addedBFacei];
1053  }
1054  }
1055 
1056  return newBoundaryData;
1057 }
1058 
1059 
1060 Foam::labelList Foam::fvMeshDistribute::mapPointData
1061 (
1062  const primitiveMesh& mesh, // mesh after adding
1063  const mapAddedPolyMesh& map,
1064  const labelList& boundaryData0, // on mesh before adding
1065  const labelList& boundaryData1 // on added mesh
1066 )
1067 {
1068  labelList newBoundaryData(mesh.nPoints());
1069 
1070  forAll(boundaryData0, oldPointi)
1071  {
1072  label newPointi = map.oldPointMap()[oldPointi];
1073 
1074  // Point still exists (is necessary?)
1075  if (newPointi >= 0)
1076  {
1077  newBoundaryData[newPointi] = boundaryData0[oldPointi];
1078  }
1079  }
1080 
1081  forAll(boundaryData1, addedPointi)
1082  {
1083  label newPointi = map.addedPointMap()[addedPointi];
1084 
1085  if (newPointi >= 0)
1086  {
1087  newBoundaryData[newPointi] = boundaryData1[addedPointi];
1088  }
1089  }
1090 
1091  return newBoundaryData;
1092 }
1093 
1094 
1095 // Remove cells. Add all exposed faces to patch oldInternalPatchi
1096 Foam::autoPtr<Foam::mapPolyMesh> Foam::fvMeshDistribute::doRemoveCells
1097 (
1098  const labelList& cellsToRemove,
1099  const label oldInternalPatchi
1100 )
1101 {
1102  // Mesh change engine
1103  polyTopoChange meshMod(mesh_);
1104 
1105  // Cell removal topo engine. Do NOT synchronize parallel since
1106  // we are doing a local cell removal.
1107  removeCells cellRemover(mesh_, false);
1108 
1109  // Get all exposed faces
1110  labelList exposedFaces(cellRemover.getExposedFaces(cellsToRemove));
1111 
1112  // Insert the topo changes
1113  cellRemover.setRefinement
1114  (
1115  cellsToRemove,
1116  exposedFaces,
1117  labelList(exposedFaces.size(), oldInternalPatchi), // patch for exposed
1118  // faces.
1119  meshMod
1120  );
1121 
1122 
1124  // tmp<surfaceScalarField> sfld(generateTestField(mesh_));
1125 
1126  // Save internal fields (note: not as DimensionedFields since would
1127  // get mapped)
1128  PtrList<Field<scalar>> sFlds;
1129  saveInternalFields(sFlds);
1130  PtrList<Field<vector>> vFlds;
1131  saveInternalFields(vFlds);
1133  saveInternalFields(sptFlds);
1134  PtrList<Field<symmTensor>> sytFlds;
1135  saveInternalFields(sytFlds);
1136  PtrList<Field<tensor>> tFlds;
1137  saveInternalFields(tFlds);
1138 
1139  // Change the mesh. No inflation. Note: no parallel comms allowed.
1140  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh_, false, false);
1141 
1142  // Update fields
1143  mesh_.updateMesh(map);
1144 
1145 
1146  // Any exposed faces in a surfaceField will not be mapped. Map the value
1147  // of these separately (until there is support in all PatchFields for
1148  // mapping from internal faces ...)
1149 
1150  mapExposedFaces(map(), sFlds);
1151  mapExposedFaces(map(), vFlds);
1152  mapExposedFaces(map(), sptFlds);
1153  mapExposedFaces(map(), sytFlds);
1154  mapExposedFaces(map(), tFlds);
1155 
1156 
1158  // testField(sfld);
1159 
1160 
1161  // Move mesh (since morphing does not do this)
1162  if (map().hasMotionPoints())
1163  {
1164  mesh_.movePoints(map().preMotionPoints());
1165  }
1166 
1167  return map;
1168 }
1169 
1170 
1171 // Delete and add processor patches. Changes mesh and returns per neighbour proc
1172 // the processor patchID.
1173 void Foam::fvMeshDistribute::addProcPatches
1174 (
1175  const labelList& nbrProc, // processor that neighbour is now on
1176  const labelList& referPatchID, // patchID (or -1) I originated from
1177  List<Map<label>>& procPatchID
1178 )
1179 {
1180  // Now use the neighbourFace/Proc to repatch the mesh. These lists
1181  // contain for all current boundary faces the global patchID (for non-proc
1182  // patch) or the processor.
1183 
1184  // Determine a visit order such that the processor patches get added
1185  // in order of increasing neighbour processor (and for same neighbour
1186  // processor (in case of processor cyclics) in order of increasing
1187  // 'refer' patch)
1188  labelList indices;
1189  sortedOrder(nbrProc, indices, lessProcPatches(nbrProc, referPatchID));
1190 
1191  procPatchID.setSize(Pstream::nProcs());
1192 
1193  forAll(indices, i)
1194  {
1195  label bFacei = indices[i];
1196  label proci = nbrProc[bFacei];
1197 
1198  if (proci != -1 && proci != Pstream::myProcNo())
1199  {
1200  if (!procPatchID[proci].found(referPatchID[bFacei]))
1201  {
1202  // No patch for neighbour yet. Is either a normal processor
1203  // patch or a processorCyclic patch.
1204 
1205  if (referPatchID[bFacei] == -1)
1206  {
1207  // Ordinary processor boundary
1208 
1210  (
1211  0, // size
1212  mesh_.nFaces(),
1213  mesh_.boundaryMesh().size(),
1214  mesh_.boundaryMesh(),
1216  proci
1217  );
1218 
1219  procPatchID[proci].insert
1220  (
1221  referPatchID[bFacei],
1223  (
1224  mesh_,
1225  pp,
1226  dictionary(), // optional per field patchField
1228  false // not parallel sync
1229  )
1230  );
1231  }
1232  else
1233  {
1234  const coupledPolyPatch& pcPatch
1235  = refCast<const coupledPolyPatch>
1236  (
1237  mesh_.boundaryMesh()[referPatchID[bFacei]]
1238  );
1240  (
1241  0, // size
1242  mesh_.nFaces(),
1243  mesh_.boundaryMesh().size(),
1244  mesh_.boundaryMesh(),
1246  proci,
1247  pcPatch.name(),
1248  pcPatch.transform()
1249  );
1250 
1251  procPatchID[proci].insert
1252  (
1253  referPatchID[bFacei],
1255  (
1256  mesh_,
1257  pp,
1258  dictionary(), // optional per field patchField
1260  false // not parallel sync
1261  )
1262  );
1263  }
1264  }
1265  }
1266  }
1267 }
1268 
1269 
1270 // Get boundary faces to be repatched. Is -1 or new patchID
1271 Foam::labelList Foam::fvMeshDistribute::getBoundaryPatch
1272 (
1273  const labelList& nbrProc, // new processor per boundary face
1274  const labelList& referPatchID, // patchID (or -1) I originated from
1275  const List<Map<label>>& procPatchID // per proc the new procPatches
1276 )
1277 {
1278  labelList patchIDs(nbrProc);
1279 
1280  forAll(nbrProc, bFacei)
1281  {
1282  if (nbrProc[bFacei] == Pstream::myProcNo())
1283  {
1284  label origPatchi = referPatchID[bFacei];
1285  patchIDs[bFacei] = origPatchi;
1286  }
1287  else if (nbrProc[bFacei] != -1)
1288  {
1289  label origPatchi = referPatchID[bFacei];
1290  patchIDs[bFacei] = procPatchID[nbrProc[bFacei]][origPatchi];
1291  }
1292  else
1293  {
1294  patchIDs[bFacei] = -1;
1295  }
1296  }
1297  return patchIDs;
1298 }
1299 
1300 
1301 // Send mesh and coupling data.
1302 void Foam::fvMeshDistribute::sendMesh
1303 (
1304  const label domain,
1305  const fvMesh& mesh,
1306 
1307  const wordList& pointZoneNames,
1308  const wordList& faceZoneNames,
1309  const wordList& cellZoneNames,
1310 
1311  const labelList& sourceFace,
1312  const labelList& sourceProc,
1313  const labelList& sourcePatch,
1314  const labelList& sourceNewNbrProc,
1315  const labelList& sourcePointMaster,
1316  Ostream& toDomain
1317 )
1318 {
1319  if (debug)
1320  {
1321  Pout<< "Sending to domain " << domain << nl
1322  << " nPoints:" << mesh.nPoints() << nl
1323  << " nFaces:" << mesh.nFaces() << nl
1324  << " nCells:" << mesh.nCells() << nl
1325  << " nPatches:" << mesh.boundaryMesh().size() << nl
1326  << endl;
1327  }
1328 
1329  // Assume sparse, possibly overlapping point zones. Get contents
1330  // in merged-zone indices.
1331  CompactListList<label> zonePoints;
1332  {
1333  const pointZoneMesh& pointZones = mesh.pointZones();
1334 
1335  labelList rowSizes(pointZoneNames.size(), 0);
1336 
1337  forAll(pointZoneNames, nameI)
1338  {
1339  label myZoneID = pointZones.findZoneID(pointZoneNames[nameI]);
1340 
1341  if (myZoneID != -1)
1342  {
1343  rowSizes[nameI] = pointZones[myZoneID].size();
1344  }
1345  }
1346  zonePoints.setSize(rowSizes);
1347 
1348  forAll(pointZoneNames, nameI)
1349  {
1350  label myZoneID = pointZones.findZoneID(pointZoneNames[nameI]);
1351 
1352  if (myZoneID != -1)
1353  {
1354  zonePoints[nameI].deepCopy(pointZones[myZoneID]);
1355  }
1356  }
1357  }
1358 
1359  // Assume sparse, possibly overlapping face zones
1360  CompactListList<label> zoneFaces;
1361  CompactListList<bool> zoneFaceFlip;
1362  {
1363  const faceZoneMesh& faceZones = mesh.faceZones();
1364 
1365  labelList rowSizes(faceZoneNames.size(), 0);
1366 
1367  forAll(faceZoneNames, nameI)
1368  {
1369  label myZoneID = faceZones.findZoneID(faceZoneNames[nameI]);
1370 
1371  if (myZoneID != -1)
1372  {
1373  rowSizes[nameI] = faceZones[myZoneID].size();
1374  }
1375  }
1376 
1377  zoneFaces.setSize(rowSizes);
1378  zoneFaceFlip.setSize(rowSizes);
1379 
1380  forAll(faceZoneNames, nameI)
1381  {
1382  label myZoneID = faceZones.findZoneID(faceZoneNames[nameI]);
1383 
1384  if (myZoneID != -1)
1385  {
1386  zoneFaces[nameI].deepCopy(faceZones[myZoneID]);
1387  zoneFaceFlip[nameI].deepCopy(faceZones[myZoneID].flipMap());
1388  }
1389  }
1390  }
1391 
1392  // Assume sparse, possibly overlapping cell zones
1393  CompactListList<label> zoneCells;
1394  {
1395  const cellZoneMesh& cellZones = mesh.cellZones();
1396 
1397  labelList rowSizes(cellZoneNames.size(), 0);
1398 
1399  forAll(cellZoneNames, nameI)
1400  {
1401  label myZoneID = cellZones.findZoneID(cellZoneNames[nameI]);
1402 
1403  if (myZoneID != -1)
1404  {
1405  rowSizes[nameI] = cellZones[myZoneID].size();
1406  }
1407  }
1408 
1409  zoneCells.setSize(rowSizes);
1410 
1411  forAll(cellZoneNames, nameI)
1412  {
1413  label myZoneID = cellZones.findZoneID(cellZoneNames[nameI]);
1414 
1415  if (myZoneID != -1)
1416  {
1417  zoneCells[nameI].deepCopy(cellZones[myZoneID]);
1418  }
1419  }
1420  }
1422  // labelList cellZoneID;
1423  // if (hasCellZones)
1424  //{
1425  // cellZoneID.setSize(mesh.nCells());
1426  // cellZoneID = -1;
1427  //
1428  // const cellZoneMesh& cellZones = mesh.cellZones();
1429  //
1430  // forAll(cellZones, zoneI)
1431  // {
1432  // UIndirectList<label>(cellZoneID, cellZones[zoneI]) = zoneI;
1433  // }
1434  //}
1435 
1436  // Send
1437  toDomain
1438  << mesh.points()
1440  << mesh.faceOwner()
1441  << mesh.faceNeighbour()
1442  << mesh.boundaryMesh()
1443 
1444  << zonePoints
1445  << zoneFaces
1446  << zoneFaceFlip
1447  << zoneCells
1448 
1449  << sourceFace
1450  << sourceProc
1451  << sourcePatch
1452  << sourceNewNbrProc
1453  << sourcePointMaster;
1454 
1455 
1456  if (debug)
1457  {
1458  Pout<< "Started sending mesh to domain " << domain
1459  << endl;
1460  }
1461 }
1462 
1463 
1464 // Receive mesh. Opposite of sendMesh
1465 Foam::autoPtr<Foam::fvMesh> Foam::fvMeshDistribute::receiveMesh
1466 (
1467  const label domain,
1468  const wordList& pointZoneNames,
1469  const wordList& faceZoneNames,
1470  const wordList& cellZoneNames,
1471  const Time& runTime,
1472  labelList& domainSourceFace,
1473  labelList& domainSourceProc,
1474  labelList& domainSourcePatch,
1475  labelList& domainSourceNewNbrProc,
1476  labelList& domainSourcePointMaster,
1477  Istream& fromNbr
1478 )
1479 {
1480  pointField domainPoints(fromNbr);
1481  faceList domainFaces = CompactListList<label, face>(fromNbr)();
1482  labelList domainAllOwner(fromNbr);
1483  labelList domainAllNeighbour(fromNbr);
1484  PtrList<entry> patchEntries(fromNbr);
1485 
1486  CompactListList<label> zonePoints(fromNbr);
1487  CompactListList<label> zoneFaces(fromNbr);
1488  CompactListList<bool> zoneFaceFlip(fromNbr);
1489  CompactListList<label> zoneCells(fromNbr);
1490 
1491  fromNbr
1492  >> domainSourceFace
1493  >> domainSourceProc
1494  >> domainSourcePatch
1495  >> domainSourceNewNbrProc
1496  >> domainSourcePointMaster;
1497 
1498  // Construct fvMesh
1499  autoPtr<fvMesh> domainMeshPtr
1500  (
1501  new fvMesh
1502  (
1503  IOobject
1504  (
1506  runTime.timeName(),
1507  runTime,
1509  ),
1510  move(domainPoints),
1511  move(domainFaces),
1512  move(domainAllOwner),
1513  move(domainAllNeighbour),
1514  false // no parallel comms
1515  )
1516  );
1517  fvMesh& domainMesh = domainMeshPtr();
1518 
1519  List<polyPatch*> patches(patchEntries.size());
1520 
1521  forAll(patchEntries, patchi)
1522  {
1524  (
1525  patchEntries[patchi].keyword(),
1526  patchEntries[patchi].dict(),
1527  patchi,
1528  domainMesh.boundaryMesh()
1529  ).ptr();
1530  }
1531  // Add patches; no parallel comms
1532  domainMesh.addFvPatches(patches, false);
1533 
1534  // Construct zones
1535  List<pointZone*> pZonePtrs(pointZoneNames.size());
1536  forAll(pZonePtrs, i)
1537  {
1538  pZonePtrs[i] = new pointZone
1539  (
1540  pointZoneNames[i],
1541  zonePoints[i],
1542  i,
1543  domainMesh.pointZones()
1544  );
1545  }
1546 
1547  List<faceZone*> fZonePtrs(faceZoneNames.size());
1548  forAll(fZonePtrs, i)
1549  {
1550  fZonePtrs[i] = new faceZone
1551  (
1552  faceZoneNames[i],
1553  zoneFaces[i],
1554  zoneFaceFlip[i],
1555  i,
1556  domainMesh.faceZones()
1557  );
1558  }
1559 
1560  List<cellZone*> cZonePtrs(cellZoneNames.size());
1561  forAll(cZonePtrs, i)
1562  {
1563  cZonePtrs[i] = new cellZone
1564  (
1565  cellZoneNames[i],
1566  zoneCells[i],
1567  i,
1568  domainMesh.cellZones()
1569  );
1570  }
1571  domainMesh.addZones(pZonePtrs, fZonePtrs, cZonePtrs);
1572 
1573  return domainMeshPtr;
1574 }
1575 
1576 
1577 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1578 
1579 // Construct from components
1581 :
1582  mesh_(mesh),
1583  mergeTol_(mergeTol)
1584 {}
1585 
1586 
1587 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1588 
1591  const labelList& distribution
1592 )
1593 {
1594  labelList nCells(Pstream::nProcs(), 0);
1595  forAll(distribution, celli)
1596  {
1597  label newProc = distribution[celli];
1598 
1599  if (newProc < 0 || newProc >= Pstream::nProcs())
1600  {
1602  << "Distribution should be in range 0.." << Pstream::nProcs()-1
1603  << endl
1604  << "At index " << celli << " distribution:" << newProc
1605  << abort(FatalError);
1606  }
1607  nCells[newProc]++;
1608  }
1609  return nCells;
1610 }
1611 
1612 
1615  const labelList& distribution
1616 )
1617 {
1618  // Some checks on distribution
1619  if (distribution.size() != mesh_.nCells())
1620  {
1622  << "Size of distribution:"
1623  << distribution.size() << " mesh nCells:" << mesh_.nCells()
1624  << abort(FatalError);
1625  }
1626 
1627 
1628  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
1629 
1630  // Check all processors have same non-proc patches in same order.
1631  if (patches.checkParallelSync(true))
1632  {
1634  << "This application requires all non-processor patches"
1635  << " to be present in the same order on all patches" << nl
1636  << "followed by the processor patches (which of course are unique)."
1637  << nl
1638  << "Local patches:" << mesh_.boundaryMesh().names()
1639  << abort(FatalError);
1640  }
1641 
1642  // Save some data for mapping later on
1643  const label nOldPoints(mesh_.nPoints());
1644  const label nOldFaces(mesh_.nFaces());
1645  const label nOldCells(mesh_.nCells());
1646  labelList oldPatchStarts(patches.size());
1647  labelList oldPatchNMeshPoints(patches.size());
1648  forAll(patches, patchi)
1649  {
1650  oldPatchStarts[patchi] = patches[patchi].start();
1651  oldPatchNMeshPoints[patchi] = patches[patchi].nPoints();
1652  }
1653 
1654 
1655  // Short circuit trivial case.
1656  if (!Pstream::parRun())
1657  {
1658  // Collect all maps and return
1660  (
1662  (
1663  mesh_,
1664 
1665  nOldPoints,
1666  nOldFaces,
1667  nOldCells,
1668  move(oldPatchStarts),
1669  move(oldPatchNMeshPoints),
1670 
1671  labelListList(1, identity(mesh_.nPoints())),
1672  labelListList(1, identity(mesh_.nFaces())),
1673  labelListList(1, identity(mesh_.nCells())),
1674  labelListList(1, identity(patches.size())),
1675 
1676  labelListList(1, identity(mesh_.nPoints())),
1677  labelListList(1, identity(mesh_.nFaces())),
1678  labelListList(1, identity(mesh_.nCells())),
1679  labelListList(1, identity(patches.size()))
1680  )
1681  );
1682  }
1683 
1684 
1685  // Collect any zone names
1686  const wordList pointZoneNames(mergeWordList(mesh_.pointZones().names()));
1687  const wordList faceZoneNames(mergeWordList(mesh_.faceZones().names()));
1688  const wordList cellZoneNames(mergeWordList(mesh_.cellZones().names()));
1689 
1690 
1691  // Local environment of all boundary faces
1692  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1693 
1694  // A face is uniquely defined by
1695  // - proc
1696  // - local face no
1697  //
1698  // To glue the parts of meshes which can get sent from anywhere we
1699  // need to know on boundary faces what the above tuple on both sides is.
1700  // So we need to maintain:
1701  // - original face
1702  // - original processor id (= trivial)
1703  // For coupled boundaries (where the faces are 'duplicate') we take the
1704  // lowest numbered processor as the data to store.
1705  //
1706  // Additionally to create the procboundaries we need to know where the owner
1707  // cell on the other side now is: newNeighbourProc.
1708  //
1709 
1710  // physical boundary:
1711  // sourceProc = -1
1712  // sourceNewNbrProc = -1
1713  // sourceFace = -1
1714  // sourcePatch = patchID
1715  // processor boundary:
1716  // sourceProc = proc (on owner side)
1717  // sourceNewNbrProc = distribution of coupled cell
1718  // sourceFace = face (on owner side)
1719  // sourcePatch = -1
1720  // ?cyclic:
1721  // ? sourceProc = proc
1722  // ? sourceNewNbrProc = distribution of coupled cell
1723  // ? sourceFace = face (on owner side)
1724  // ? sourcePatch = patchID
1725  // processor-cyclic boundary:
1726  // sourceProc = proc (on owner side)
1727  // sourceNewNbrProc = distribution of coupled cell
1728  // sourceFace = face (on owner side)
1729  // sourcePatch = patchID
1730 
1731  labelList sourcePatch;
1732  labelList sourceFace;
1733  labelList sourceProc;
1734  labelList sourceNewNbrProc;
1735  labelList sourcePointMaster;
1736  getCouplingData
1737  (
1738  distribution,
1739  sourceFace,
1740  sourceProc,
1741  sourcePatch,
1742  sourceNewNbrProc,
1743  sourcePointMaster
1744  );
1745 
1746 
1747  // Remove meshPhi. Since this would otherwise disappear anyway
1748  // during topo changes and we have to guarantee that all the fields
1749  // can be sent.
1750  //mesh_.clearOut();
1751  mesh_.resetMotion();
1752 
1753  // Get data to send. Make sure is synchronised
1754  const wordList volScalars(mesh_.names(volScalarField::typeName));
1755  checkEqualWordList("volScalarFields", volScalars);
1756  const wordList volVectors(mesh_.names(volVectorField::typeName));
1757  checkEqualWordList("volVectorFields", volVectors);
1758  const wordList volSphereTensors
1759  (
1761  );
1762  checkEqualWordList("volSphericalTensorFields", volSphereTensors);
1763  const wordList volSymmTensors(mesh_.names(volSymmTensorField::typeName));
1764  checkEqualWordList("volSymmTensorFields", volSymmTensors);
1765  const wordList volTensors(mesh_.names(volTensorField::typeName));
1766  checkEqualWordList("volTensorField", volTensors);
1767 
1768  const wordList surfScalars(mesh_.names(surfaceScalarField::typeName));
1769  checkEqualWordList("surfaceScalarFields", surfScalars);
1770  const wordList surfVectors(mesh_.names(surfaceVectorField::typeName));
1771  checkEqualWordList("surfaceVectorFields", surfVectors);
1772  const wordList surfSphereTensors
1773  (
1775  );
1776  checkEqualWordList("surfaceSphericalTensorFields", surfSphereTensors);
1777  const wordList surfSymmTensors
1778  (
1780  );
1781  checkEqualWordList("surfaceSymmTensorFields", surfSymmTensors);
1782  const wordList surfTensors(mesh_.names(surfaceTensorField::typeName));
1783  checkEqualWordList("surfaceTensorFields", surfTensors);
1784 
1785  const wordList pointScalars(mesh_.names(pointScalarField::typeName));
1786  checkEqualWordList("pointScalarFields", pointScalars);
1787  const wordList pointVectors(mesh_.names(pointVectorField::typeName));
1788  checkEqualWordList("pointVectorFields", pointVectors);
1789  const wordList pointSphereTensors
1790  (
1792  );
1793  checkEqualWordList("pointSphericalTensorFields", pointSphereTensors);
1794  const wordList pointSymmTensors
1795  (
1797  );
1798  checkEqualWordList("pointSymmTensorFields", pointSymmTensors);
1799  const wordList pointTensors(mesh_.names(pointTensorField::typeName));
1800  checkEqualWordList("pointTensorFields", pointTensors);
1801 
1802 
1803 
1804  typedef volScalarField::Internal dimScalType;
1805  const wordList dimScalars(mesh_.names(dimScalType::typeName));
1806  checkEqualWordList("volScalarField::Internal", dimScalars);
1807 
1808  typedef volVectorField::Internal dimVecType;
1809  const wordList dimVectors(mesh_.names(dimVecType::typeName));
1810  checkEqualWordList("volVectorField::Internal", dimVectors);
1811 
1812  typedef volSphericalTensorField::Internal dimSphereType;
1813  const wordList dimSphereTensors(mesh_.names(dimSphereType::typeName));
1814  checkEqualWordList
1815  (
1816  "volSphericalTensorField::Internal",
1817  dimSphereTensors
1818  );
1819 
1820  typedef volSymmTensorField::Internal dimSymmTensorType;
1821  const wordList dimSymmTensors(mesh_.names(dimSymmTensorType::typeName));
1822  checkEqualWordList
1823  (
1824  "volSymmTensorField::Internal",
1825  dimSymmTensors
1826  );
1827 
1828  typedef volTensorField::Internal dimTensorType;
1829  const wordList dimTensors(mesh_.names(dimTensorType::typeName));
1830  checkEqualWordList("volTensorField::Internal", dimTensors);
1831 
1832 
1833 
1834  // Find patch to temporarily put exposed and processor faces into.
1835  label oldInternalPatchi = findNonEmptyPatch();
1836 
1837 
1838 
1839  // Delete processor patches, starting from the back. Move all faces into
1840  // oldInternalPatchi.
1841  labelList repatchFaceMap;
1842  {
1843  autoPtr<mapPolyMesh> repatchMap = deleteProcPatches(oldInternalPatchi);
1844 
1845  // Store face map (only face ordering that changed)
1846  repatchFaceMap = repatchMap().faceMap();
1847 
1848  // Reorder all boundary face data (sourceProc, sourceFace etc.)
1849  labelList bFaceMap
1850  (
1852  (
1853  repatchMap().reverseFaceMap(),
1854  mesh_.nFaces() - mesh_.nInternalFaces(),
1855  mesh_.nInternalFaces()
1856  )
1857  - mesh_.nInternalFaces()
1858  );
1859 
1860  inplaceReorder(bFaceMap, sourceFace);
1861  inplaceReorder(bFaceMap, sourceProc);
1862  inplaceReorder(bFaceMap, sourcePatch);
1863  inplaceReorder(bFaceMap, sourceNewNbrProc);
1864  }
1865 
1866 
1867 
1868  // Print a bit.
1869  if (debug)
1870  {
1871  Pout<< nl << "MESH WITH PROC PATCHES DELETED:" << endl;
1872  printMeshInfo(mesh_);
1873  printFieldInfo<volScalarField>(mesh_);
1874  printFieldInfo<volVectorField>(mesh_);
1875  printFieldInfo<volSphericalTensorField>(mesh_);
1876  printFieldInfo<volSymmTensorField>(mesh_);
1877  printFieldInfo<volTensorField>(mesh_);
1878  printFieldInfo<surfaceScalarField>(mesh_);
1879  printFieldInfo<surfaceVectorField>(mesh_);
1880  printFieldInfo<surfaceSphericalTensorField>(mesh_);
1881  printFieldInfo<surfaceSymmTensorField>(mesh_);
1882  printFieldInfo<surfaceTensorField>(mesh_);
1883  printFieldInfo<pointScalarField>(mesh_);
1884  printFieldInfo<pointVectorField>(mesh_);
1885  printFieldInfo<pointSphericalTensorField>(mesh_);
1886  printFieldInfo<pointSymmTensorField>(mesh_);
1887  printFieldInfo<pointTensorField>(mesh_);
1888  Pout<< nl << endl;
1889  }
1890 
1891 
1892 
1893  // Maps from subsetted mesh (that is sent) back to original maps
1894  labelListList subCellMap(Pstream::nProcs());
1895  labelListList subFaceMap(Pstream::nProcs());
1896  labelListList subPointMap(Pstream::nProcs());
1897  labelListList subPatchMap(Pstream::nProcs());
1898  // Maps from subsetted mesh to reconstructed mesh
1899  labelListList constructCellMap(Pstream::nProcs());
1900  labelListList constructFaceMap(Pstream::nProcs());
1901  labelListList constructPointMap(Pstream::nProcs());
1902  labelListList constructPatchMap(Pstream::nProcs());
1903 
1904 
1905 
1906 
1907  // Find out schedule
1908  // ~~~~~~~~~~~~~~~~~
1909 
1910  labelListList nSendCells(Pstream::nProcs());
1911  nSendCells[Pstream::myProcNo()] = countCells(distribution);
1912  Pstream::gatherList(nSendCells);
1913  Pstream::scatterList(nSendCells);
1914 
1915 
1916  // Allocate buffers
1918 
1919 
1920  // What to send to neighbouring domains
1921  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1922 
1923  bool oldParRun = UPstream::parRun();
1924  UPstream::parRun() = false;
1925 
1926  forAll(nSendCells[Pstream::myProcNo()], recvProc)
1927  {
1928  if
1929  (
1930  recvProc != Pstream::myProcNo()
1931  && nSendCells[Pstream::myProcNo()][recvProc] > 0
1932  )
1933  {
1934  // Send to recvProc
1935 
1936  if (debug)
1937  {
1938  Pout<< nl
1939  << "SUBSETTING FOR DOMAIN " << recvProc
1940  << " cells to send:"
1941  << nSendCells[Pstream::myProcNo()][recvProc]
1942  << nl << endl;
1943  }
1944 
1945  // Pstream for sending mesh and fields
1946  // OPstream str(Pstream::commsTypes::blocking, recvProc);
1947  UOPstream str(recvProc, pBufs);
1948 
1949  // Mesh subsetting engine
1950  fvMeshSubset subsetter(mesh_);
1951 
1952  // Subset the cells of the current domain.
1953  subsetter.setLargeCellSubset
1954  (
1955  distribution,
1956  recvProc,
1957  oldInternalPatchi, // oldInternalFaces patch
1958  false // no parallel sync
1959  );
1960 
1961  subCellMap[recvProc] = subsetter.cellMap();
1962  subFaceMap[recvProc] = subsetter.faceFlipMap();
1963  inplaceRenumberWithFlip
1964  (
1965  repatchFaceMap,
1966  false, // oldToNew has flip
1967  true, // subFaceMap has flip
1968  subFaceMap[recvProc]
1969  );
1970  subPointMap[recvProc] = subsetter.pointMap();
1971  subPatchMap[recvProc] = subsetter.patchMap();
1972 
1973 
1974  // Subset the boundary fields (owner/neighbour/processor)
1975  labelList procSourceFace;
1976  labelList procSourceProc;
1977  labelList procSourcePatch;
1978  labelList procSourceNewNbrProc;
1979  labelList procSourcePointMaster;
1980 
1981  subsetCouplingData
1982  (
1983  subsetter.subMesh(),
1984  subsetter.pointMap(), // from subMesh to mesh
1985  subsetter.faceMap(), // ,, ,,
1986  subsetter.cellMap(), // ,, ,,
1987 
1988  distribution, // old mesh distribution
1989  mesh_.faceOwner(), // old owner
1990  mesh_.faceNeighbour(),
1991  mesh_.nInternalFaces(),
1992 
1993  sourceFace,
1994  sourceProc,
1995  sourcePatch,
1996  sourceNewNbrProc,
1997  sourcePointMaster,
1998 
1999  procSourceFace,
2000  procSourceProc,
2001  procSourcePatch,
2002  procSourceNewNbrProc,
2003  procSourcePointMaster
2004  );
2005 
2006 
2007  // Send to neighbour
2008  sendMesh
2009  (
2010  recvProc,
2011  subsetter.subMesh(),
2012 
2013  pointZoneNames,
2014  faceZoneNames,
2015  cellZoneNames,
2016 
2017  procSourceFace,
2018  procSourceProc,
2019  procSourcePatch,
2020  procSourceNewNbrProc,
2021  procSourcePointMaster,
2022 
2023  str
2024  );
2025 
2026  // volFields
2027  sendFields<volScalarField>(recvProc, volScalars, subsetter, str);
2028  sendFields<volVectorField>(recvProc, volVectors, subsetter, str);
2029  sendFields<volSphericalTensorField>
2030  (
2031  recvProc,
2032  volSphereTensors,
2033  subsetter,
2034  str
2035  );
2036  sendFields<volSymmTensorField>
2037  (
2038  recvProc,
2039  volSymmTensors,
2040  subsetter,
2041  str
2042  );
2043  sendFields<volTensorField>(recvProc, volTensors, subsetter, str);
2044 
2045  // surfaceFields
2046  sendFields<surfaceScalarField>
2047  (
2048  recvProc,
2049  surfScalars,
2050  subsetter,
2051  str
2052  );
2053  sendFields<surfaceVectorField>
2054  (
2055  recvProc,
2056  surfVectors,
2057  subsetter,
2058  str
2059  );
2060  sendFields<surfaceSphericalTensorField>
2061  (
2062  recvProc,
2063  surfSphereTensors,
2064  subsetter,
2065  str
2066  );
2067  sendFields<surfaceSymmTensorField>
2068  (
2069  recvProc,
2070  surfSymmTensors,
2071  subsetter,
2072  str
2073  );
2074  sendFields<surfaceTensorField>
2075  (
2076  recvProc,
2077  surfTensors,
2078  subsetter,
2079  str
2080  );
2081 
2082  // pointFields
2083  sendFields<pointScalarField>
2084  (
2085  recvProc,
2086  pointScalars,
2087  subsetter,
2088  str
2089  );
2090  sendFields<pointVectorField>
2091  (
2092  recvProc,
2093  pointVectors,
2094  subsetter,
2095  str
2096  );
2097  sendFields<pointSphericalTensorField>
2098  (
2099  recvProc,
2100  pointSphereTensors,
2101  subsetter,
2102  str
2103  );
2104  sendFields<pointSymmTensorField>
2105  (
2106  recvProc,
2107  pointSymmTensors,
2108  subsetter,
2109  str
2110  );
2111  sendFields<pointTensorField>
2112  (
2113  recvProc,
2114  pointTensors,
2115  subsetter,
2116  str
2117  );
2118 
2119  // dimensionedFields
2120  sendFields<volScalarField::Internal>
2121  (
2122  recvProc,
2123  dimScalars,
2124  subsetter,
2125  str
2126  );
2127  sendFields<volVectorField::Internal>
2128  (
2129  recvProc,
2130  dimVectors,
2131  subsetter,
2132  str
2133  );
2134  sendFields<volSphericalTensorField::Internal>
2135  (
2136  recvProc,
2137  dimSphereTensors,
2138  subsetter,
2139  str
2140  );
2141  sendFields<volSymmTensorField::Internal>
2142  (
2143  recvProc,
2144  dimSymmTensors,
2145  subsetter,
2146  str
2147  );
2148  sendFields<volTensorField::Internal>
2149  (
2150  recvProc,
2151  dimTensors,
2152  subsetter,
2153  str
2154  );
2155  }
2156  }
2157 
2158 
2159  UPstream::parRun() = oldParRun;
2160 
2161 
2162  // Start sending&receiving from buffers
2163  pBufs.finishedSends();
2164 
2165 
2166  // Subset the part that stays
2167  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
2168 
2169  {
2170  // Save old mesh maps before changing mesh
2171  const labelList oldFaceOwner(mesh_.faceOwner());
2172  const labelList oldFaceNeighbour(mesh_.faceNeighbour());
2173  const label oldInternalFaces = mesh_.nInternalFaces();
2174 
2175  // Remove cells.
2176  autoPtr<mapPolyMesh> subMap
2177  (
2178  doRemoveCells
2179  (
2180  select(false, distribution, Pstream::myProcNo()),
2181  oldInternalPatchi
2182  )
2183  );
2184 
2185  // Addressing from subsetted mesh
2186  subCellMap[Pstream::myProcNo()] = subMap().cellMap();
2187  subFaceMap[Pstream::myProcNo()] = renumber
2188  (
2189  repatchFaceMap,
2190  subMap().faceMap()
2191  );
2192  // Insert the sign bit from face flipping
2193  labelList& faceMap = subFaceMap[Pstream::myProcNo()];
2194  forAll(faceMap, faceI)
2195  {
2196  faceMap[faceI] += 1;
2197  }
2198  const labelHashSet& flip = subMap().flipFaceFlux();
2199  forAllConstIter(labelHashSet, flip, iter)
2200  {
2201  label faceI = iter.key();
2202  faceMap[faceI] = -faceMap[faceI];
2203  }
2204  subPointMap[Pstream::myProcNo()] = subMap().pointMap();
2205  subPatchMap[Pstream::myProcNo()] = identity(patches.size());
2206 
2207  // Initialize all addressing into current mesh
2208  constructCellMap[Pstream::myProcNo()] = identity(mesh_.nCells());
2209  constructFaceMap[Pstream::myProcNo()] = identity(mesh_.nFaces()) + 1;
2210  constructPointMap[Pstream::myProcNo()] = identity(mesh_.nPoints());
2211  constructPatchMap[Pstream::myProcNo()] = identity(patches.size());
2212 
2213  // Subset the mesh data: neighbourCell/neighbourProc
2214  // fields
2215  labelList domainSourceFace;
2216  labelList domainSourceProc;
2217  labelList domainSourcePatch;
2218  labelList domainSourceNewNbrProc;
2219  labelList domainSourcePointMaster;
2220 
2221  subsetCouplingData
2222  (
2223  mesh_, // new mesh
2224  subMap().pointMap(), // from new to original mesh
2225  subMap().faceMap(), // from new to original mesh
2226  subMap().cellMap(),
2227 
2228  distribution, // distribution before subsetting
2229  oldFaceOwner, // owner before subsetting
2230  oldFaceNeighbour, // neighbour ,,
2231  oldInternalFaces, // nInternalFaces ,,
2232 
2233  sourceFace,
2234  sourceProc,
2235  sourcePatch,
2236  sourceNewNbrProc,
2237  sourcePointMaster,
2238 
2239  domainSourceFace,
2240  domainSourceProc,
2241  domainSourcePatch,
2242  domainSourceNewNbrProc,
2243  domainSourcePointMaster
2244  );
2245 
2246  sourceFace.transfer(domainSourceFace);
2247  sourceProc.transfer(domainSourceProc);
2248  sourcePatch.transfer(domainSourcePatch);
2249  sourceNewNbrProc.transfer(domainSourceNewNbrProc);
2250  sourcePointMaster.transfer(domainSourcePointMaster);
2251  }
2252 
2253 
2254  // Print a bit.
2255  if (debug)
2256  {
2257  Pout<< nl << "STARTING MESH:" << endl;
2258  printMeshInfo(mesh_);
2259  printFieldInfo<volScalarField>(mesh_);
2260  printFieldInfo<volVectorField>(mesh_);
2261  printFieldInfo<volSphericalTensorField>(mesh_);
2262  printFieldInfo<volSymmTensorField>(mesh_);
2263  printFieldInfo<volTensorField>(mesh_);
2264  printFieldInfo<surfaceScalarField>(mesh_);
2265  printFieldInfo<surfaceVectorField>(mesh_);
2266  printFieldInfo<surfaceSphericalTensorField>(mesh_);
2267  printFieldInfo<surfaceSymmTensorField>(mesh_);
2268  printFieldInfo<surfaceTensorField>(mesh_);
2269  printFieldInfo<pointScalarField>(mesh_);
2270  printFieldInfo<pointVectorField>(mesh_);
2271  printFieldInfo<pointSphericalTensorField>(mesh_);
2272  printFieldInfo<pointSymmTensorField>(mesh_);
2273  printFieldInfo<pointTensorField>(mesh_);
2274  Pout<< nl << endl;
2275  }
2276 
2277 
2278 
2279  // Receive and add what was sent
2280  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2281 
2282  oldParRun = UPstream::parRun();
2283  UPstream::parRun() = false;
2284 
2285  forAll(nSendCells, sendProc)
2286  {
2287  // Did processor sendProc send anything to me?
2288  if
2289  (
2290  sendProc != Pstream::myProcNo()
2291  && nSendCells[sendProc][Pstream::myProcNo()] > 0
2292  )
2293  {
2294  if (debug)
2295  {
2296  Pout<< nl
2297  << "RECEIVING FROM DOMAIN " << sendProc
2298  << " cells to receive:"
2299  << nSendCells[sendProc][Pstream::myProcNo()]
2300  << nl << endl;
2301  }
2302 
2303 
2304  // Pstream for receiving mesh and fields
2305  UIPstream str(sendProc, pBufs);
2306 
2307 
2308  // Receive from sendProc
2309  labelList domainSourceFace;
2310  labelList domainSourceProc;
2311  labelList domainSourcePatch;
2312  labelList domainSourceNewNbrProc;
2313  labelList domainSourcePointMaster;
2314 
2315  autoPtr<fvMesh> domainMeshPtr;
2316 
2322 
2328 
2334 
2340 
2341 
2342  // Opposite of sendMesh
2343  {
2344  domainMeshPtr = receiveMesh
2345  (
2346  sendProc,
2347  pointZoneNames,
2348  faceZoneNames,
2349  cellZoneNames,
2350 
2351  const_cast<Time&>(mesh_.time()),
2352  domainSourceFace,
2353  domainSourceProc,
2354  domainSourcePatch,
2355  domainSourceNewNbrProc,
2356  domainSourcePointMaster,
2357  str
2358  );
2359  fvMesh& domainMesh = domainMeshPtr();
2360  // Force construction of various on mesh.
2361  //(void)domainMesh.globalData();
2362 
2363 
2364  // Receive fields. Read as single dictionary because
2365  // of problems reading consecutive fields from single stream.
2366  dictionary fieldDicts(str);
2367 
2368  // Vol fields
2369  receiveFields<volScalarField>
2370  (
2371  sendProc,
2372  volScalars,
2373  domainMesh,
2374  vsf,
2375  fieldDicts.subDict(volScalarField::typeName)
2376  );
2377  receiveFields<volVectorField>
2378  (
2379  sendProc,
2380  volVectors,
2381  domainMesh,
2382  vvf,
2383  fieldDicts.subDict(volVectorField::typeName)
2384  );
2385  receiveFields<volSphericalTensorField>
2386  (
2387  sendProc,
2388  volSphereTensors,
2389  domainMesh,
2390  vsptf,
2392  );
2393  receiveFields<volSymmTensorField>
2394  (
2395  sendProc,
2396  volSymmTensors,
2397  domainMesh,
2398  vsytf,
2400  );
2401  receiveFields<volTensorField>
2402  (
2403  sendProc,
2404  volTensors,
2405  domainMesh,
2406  vtf,
2407  fieldDicts.subDict(volTensorField::typeName)
2408  );
2409 
2410  // Surface fields
2411  receiveFields<surfaceScalarField>
2412  (
2413  sendProc,
2414  surfScalars,
2415  domainMesh,
2416  ssf,
2418  );
2419  receiveFields<surfaceVectorField>
2420  (
2421  sendProc,
2422  surfVectors,
2423  domainMesh,
2424  svf,
2426  );
2427  receiveFields<surfaceSphericalTensorField>
2428  (
2429  sendProc,
2430  surfSphereTensors,
2431  domainMesh,
2432  ssptf,
2434  );
2435  receiveFields<surfaceSymmTensorField>
2436  (
2437  sendProc,
2438  surfSymmTensors,
2439  domainMesh,
2440  ssytf,
2442  );
2443  receiveFields<surfaceTensorField>
2444  (
2445  sendProc,
2446  surfTensors,
2447  domainMesh,
2448  stf,
2450  );
2451 
2452  // Point fields
2453  pointMesh& domainPointMesh =
2454  const_cast<pointMesh&>(pointMesh::New(domainMesh));
2455  receiveFields<pointScalarField>
2456  (
2457  sendProc,
2458  pointScalars,
2459  domainPointMesh,
2460  psf,
2462  );
2463  receiveFields<pointVectorField>
2464  (
2465  sendProc,
2466  pointVectors,
2467  domainPointMesh,
2468  pvf,
2470  );
2471  receiveFields<pointSphericalTensorField>
2472  (
2473  sendProc,
2474  pointSphereTensors,
2475  domainPointMesh,
2476  psptf,
2478  );
2479  receiveFields<pointSymmTensorField>
2480  (
2481  sendProc,
2482  pointSymmTensors,
2483  domainPointMesh,
2484  psytf,
2486  );
2487  receiveFields<pointTensorField>
2488  (
2489  sendProc,
2490  pointTensors,
2491  domainPointMesh,
2492  ptf,
2494  );
2495 
2496  // Dimensioned fields
2497  receiveFields<volScalarField::Internal>
2498  (
2499  sendProc,
2500  dimScalars,
2501  domainMesh,
2502  dsf,
2503  fieldDicts.subDict
2504  (
2506  )
2507  );
2508  receiveFields<volVectorField::Internal>
2509  (
2510  sendProc,
2511  dimVectors,
2512  domainMesh,
2513  dvf,
2514  fieldDicts.subDict
2515  (
2517  )
2518  );
2519  receiveFields<volSphericalTensorField::Internal>
2520  (
2521  sendProc,
2522  dimSphereTensors,
2523  domainMesh,
2524  dstf,
2525  fieldDicts.subDict
2526  (
2528  typeName
2529  )
2530  );
2531  receiveFields<volSymmTensorField::Internal>
2532  (
2533  sendProc,
2534  dimSymmTensors,
2535  domainMesh,
2536  dsytf,
2537  fieldDicts.subDict
2538  (
2540  )
2541  );
2542  receiveFields<volTensorField::Internal>
2543  (
2544  sendProc,
2545  dimTensors,
2546  domainMesh,
2547  dtf,
2548  fieldDicts.subDict
2549  (
2551  )
2552  );
2553  }
2554  const fvMesh& domainMesh = domainMeshPtr();
2555 
2556 
2557  constructCellMap[sendProc] = identity(domainMesh.nCells());
2558  constructFaceMap[sendProc] = identity(domainMesh.nFaces()) + 1;
2559  constructPointMap[sendProc] = identity(domainMesh.nPoints());
2560  constructPatchMap[sendProc] =
2561  identity(domainMesh.boundaryMesh().size());
2562 
2563 
2564  // Print a bit.
2565  if (debug)
2566  {
2567  Pout<< nl << "RECEIVED MESH FROM:" << sendProc << endl;
2568  printMeshInfo(domainMesh);
2569  printFieldInfo<volScalarField>(domainMesh);
2570  printFieldInfo<volVectorField>(domainMesh);
2571  printFieldInfo<volSphericalTensorField>(domainMesh);
2572  printFieldInfo<volSymmTensorField>(domainMesh);
2573  printFieldInfo<volTensorField>(domainMesh);
2574  printFieldInfo<surfaceScalarField>(domainMesh);
2575  printFieldInfo<surfaceVectorField>(domainMesh);
2576  printFieldInfo<surfaceSphericalTensorField>(domainMesh);
2577  printFieldInfo<surfaceSymmTensorField>(domainMesh);
2578  printFieldInfo<surfaceTensorField>(domainMesh);
2579  printFieldInfo<pointScalarField>(domainMesh);
2580  printFieldInfo<pointVectorField>(domainMesh);
2581  printFieldInfo<pointSphericalTensorField>(domainMesh);
2582  printFieldInfo<pointSymmTensorField>(domainMesh);
2583  printFieldInfo<pointTensorField>(domainMesh);
2584  }
2585 
2586 
2587  // Now this mesh we received (from sendProc) needs to be merged
2588  // with the current mesh. On the current mesh we have for all
2589  // boundaryfaces the original face and processor. See if we can
2590  // match these up to the received domainSourceFace and
2591  // domainSourceProc.
2592  labelList masterCoupledFaces;
2593  labelList slaveCoupledFaces;
2594  findCouples
2595  (
2596  mesh_,
2597 
2598  sourceFace,
2599  sourceProc,
2600  sourcePatch,
2601 
2602  sendProc,
2603  domainMesh,
2604  domainSourceFace,
2605  domainSourceProc,
2606  domainSourcePatch,
2607 
2608  masterCoupledFaces,
2609  slaveCoupledFaces
2610  );
2611 
2612  // Generate additional coupling info (points, edges) from
2613  // faces-that-match
2614  faceCoupleInfo couples
2615  (
2616  mesh_,
2617  masterCoupledFaces,
2618  domainMesh,
2619  slaveCoupledFaces,
2620  mergeTol_, // merge tolerance
2621  true, // faces align
2622  true, // couples are ordered already
2623  false
2624  );
2625 
2626 
2627  // Add domainMesh to mesh
2628  // ~~~~~~~~~~~~~~~~~~~~~~
2629 
2631  (
2632  mesh_,
2633  domainMesh,
2634  couples,
2635  false // no parallel comms
2636  );
2637 
2638  // Update mesh data: sourceFace,sourceProc for added
2639  // mesh.
2640 
2641  sourceFace = mapBoundaryData
2642  (
2643  mesh_,
2644  map(),
2645  sourceFace,
2646  domainMesh.nInternalFaces(),
2647  domainSourceFace
2648  );
2649  sourceProc = mapBoundaryData
2650  (
2651  mesh_,
2652  map(),
2653  sourceProc,
2654  domainMesh.nInternalFaces(),
2655  domainSourceProc
2656  );
2657  sourcePatch = mapBoundaryData
2658  (
2659  mesh_,
2660  map(),
2661  sourcePatch,
2662  domainMesh.nInternalFaces(),
2663  domainSourcePatch
2664  );
2665  sourceNewNbrProc = mapBoundaryData
2666  (
2667  mesh_,
2668  map(),
2669  sourceNewNbrProc,
2670  domainMesh.nInternalFaces(),
2671  domainSourceNewNbrProc
2672  );
2673  // Update pointMaster data
2674  sourcePointMaster = mapPointData
2675  (
2676  mesh_,
2677  map(),
2678  sourcePointMaster,
2679  domainSourcePointMaster
2680  );
2681 
2682 
2683  // Update all addressing so xxProcAddressing points to correct
2684  // item in masterMesh.
2685  const labelList& oldCellMap = map().oldCellMap();
2686  const labelList& oldFaceMap = map().oldFaceMap();
2687  const labelList& oldPointMap = map().oldPointMap();
2688  const labelList& oldPatchMap = map().oldPatchMap();
2689 
2690  // Note: old mesh faces never flipped!
2691  forAll(constructPatchMap, proci)
2692  {
2693  if (proci != sendProc && constructPatchMap[proci].size())
2694  {
2695  // Processor already in mesh (either myProcNo or received)
2696  inplaceRenumber(oldCellMap, constructCellMap[proci]);
2697  inplaceRenumberWithFlip
2698  (
2699  oldFaceMap,
2700  false,
2701  true,
2702  constructFaceMap[proci]
2703  );
2704  inplaceRenumber(oldPointMap, constructPointMap[proci]);
2705  inplaceRenumber(oldPatchMap, constructPatchMap[proci]);
2706  }
2707  }
2708 
2709 
2710  labelHashSet flippedAddedFaces;
2711  {
2712  // Find out if any faces of domain mesh were flipped (boundary
2713  // faces becoming internal)
2714  label nBnd = domainMesh.nFaces()-domainMesh.nInternalFaces();
2715  flippedAddedFaces.resize(nBnd/4);
2716 
2717  for
2718  (
2719  label domainFaceI = domainMesh.nInternalFaces();
2720  domainFaceI < domainMesh.nFaces();
2721  domainFaceI++
2722  )
2723  {
2724  label newFaceI = map().addedFaceMap()[domainFaceI];
2725  label newCellI = mesh_.faceOwner()[newFaceI];
2726 
2727  label domainCellI = domainMesh.faceOwner()[domainFaceI];
2728 
2729  if (newCellI != map().addedCellMap()[domainCellI])
2730  {
2731  flippedAddedFaces.insert(domainFaceI);
2732  }
2733  }
2734  }
2735 
2736 
2737  // Added processor
2738  inplaceRenumber(map().addedCellMap(), constructCellMap[sendProc]);
2739  // Add flip
2740  forAllConstIter(labelHashSet, flippedAddedFaces, iter)
2741  {
2742  label domainFaceI = iter.key();
2743  label& val = constructFaceMap[sendProc][domainFaceI];
2744  val = -val;
2745  }
2746  inplaceRenumberWithFlip
2747  (
2748  map().addedFaceMap(),
2749  false,
2750  true, // constructFaceMap has flip sign
2751  constructFaceMap[sendProc]
2752  );
2753  inplaceRenumber(map().addedPointMap(), constructPointMap[sendProc]);
2754  inplaceRenumber(map().addedPatchMap(), constructPatchMap[sendProc]);
2755 
2756  if (debug)
2757  {
2758  Pout<< nl << "MERGED MESH FROM:" << sendProc << endl;
2759  printMeshInfo(mesh_);
2760  printFieldInfo<volScalarField>(mesh_);
2761  printFieldInfo<volVectorField>(mesh_);
2762  printFieldInfo<volSphericalTensorField>(mesh_);
2763  printFieldInfo<volSymmTensorField>(mesh_);
2764  printFieldInfo<volTensorField>(mesh_);
2765  printFieldInfo<surfaceScalarField>(mesh_);
2766  printFieldInfo<surfaceVectorField>(mesh_);
2767  printFieldInfo<surfaceSphericalTensorField>(mesh_);
2768  printFieldInfo<surfaceSymmTensorField>(mesh_);
2769  printFieldInfo<surfaceTensorField>(mesh_);
2770  printFieldInfo<pointScalarField>(mesh_);
2771  printFieldInfo<pointVectorField>(mesh_);
2772  printFieldInfo<pointSphericalTensorField>(mesh_);
2773  printFieldInfo<pointSymmTensorField>(mesh_);
2774  printFieldInfo<pointTensorField>(mesh_);
2775  Pout<< nl << endl;
2776  }
2777  }
2778  }
2779 
2780  UPstream::parRun() = oldParRun;
2781 
2782  // Print a bit.
2783  if (debug)
2784  {
2785  Pout<< nl << "REDISTRIBUTED MESH:" << endl;
2786  printMeshInfo(mesh_);
2787  printFieldInfo<volScalarField>(mesh_);
2788  printFieldInfo<volVectorField>(mesh_);
2789  printFieldInfo<volSphericalTensorField>(mesh_);
2790  printFieldInfo<volSymmTensorField>(mesh_);
2791  printFieldInfo<volTensorField>(mesh_);
2792  printFieldInfo<surfaceScalarField>(mesh_);
2793  printFieldInfo<surfaceVectorField>(mesh_);
2794  printFieldInfo<surfaceSphericalTensorField>(mesh_);
2795  printFieldInfo<surfaceSymmTensorField>(mesh_);
2796  printFieldInfo<surfaceTensorField>(mesh_);
2797  printFieldInfo<pointScalarField>(mesh_);
2798  printFieldInfo<pointVectorField>(mesh_);
2799  printFieldInfo<pointSphericalTensorField>(mesh_);
2800  printFieldInfo<pointSymmTensorField>(mesh_);
2801  printFieldInfo<pointTensorField>(mesh_);
2802  Pout<< nl << endl;
2803  }
2804 
2805 
2806  // See if any originally shared points need to be merged. Note: does
2807  // parallel comms. After this points and edges should again be consistent.
2808  mergeSharedPoints(sourcePointMaster, constructPointMap);
2809 
2810 
2811  // Add processorPatches
2812  // ~~~~~~~~~~~~~~~~~~~~
2813 
2814  // Per neighbour processor, per originating patch, the patchID
2815  // For faces resulting from internal faces or normal processor patches
2816  // the originating patch is -1. For cyclics this is the cyclic patchID.
2817  List<Map<label>> procPatchID;
2818 
2819  // Add processor and processorCyclic patches.
2820  addProcPatches(sourceNewNbrProc, sourcePatch, procPatchID);
2821 
2822  // Put faces into correct patch. Note that we now have proper
2823  // processorPolyPatches again so repatching will take care of coupled face
2824  // ordering.
2825 
2826  // Get boundary faces to be repatched. Is -1 or new patchID
2827  labelList newPatchID
2828  (
2829  getBoundaryPatch
2830  (
2831  sourceNewNbrProc,
2832  sourcePatch,
2833  procPatchID
2834  )
2835  );
2836 
2837  // Change patches. Since this might change ordering of coupled faces
2838  // we also need to adapt our constructMaps.
2839  repatch(newPatchID, constructFaceMap);
2840 
2841  // Bit of hack: processorFvPatchField does not get reset since created
2842  // from nothing so explicitly reset.
2843  initPatchFields<volScalarField, processorFvPatchField<scalar>>
2844  (
2845  Zero
2846  );
2847  initPatchFields<volVectorField, processorFvPatchField<vector>>
2848  (
2849  Zero
2850  );
2851  initPatchFields
2852  <
2855  >
2856  (
2857  Zero
2858  );
2859  initPatchFields<volSymmTensorField, processorFvPatchField<symmTensor>>
2860  (
2861  Zero
2862  );
2863  initPatchFields<volTensorField, processorFvPatchField<tensor>>
2864  (
2865  Zero
2866  );
2867 
2868 
2869  mesh_.setInstance(mesh_.time().timeName());
2870 
2871 
2872  // Print a bit
2873  if (debug)
2874  {
2875  Pout<< nl << "FINAL MESH:" << endl;
2876  printMeshInfo(mesh_);
2877  printFieldInfo<volScalarField>(mesh_);
2878  printFieldInfo<volVectorField>(mesh_);
2879  printFieldInfo<volSphericalTensorField>(mesh_);
2880  printFieldInfo<volSymmTensorField>(mesh_);
2881  printFieldInfo<volTensorField>(mesh_);
2882  printFieldInfo<surfaceScalarField>(mesh_);
2883  printFieldInfo<surfaceVectorField>(mesh_);
2884  printFieldInfo<surfaceSphericalTensorField>(mesh_);
2885  printFieldInfo<surfaceSymmTensorField>(mesh_);
2886  printFieldInfo<surfaceTensorField>(mesh_);
2887  printFieldInfo<pointScalarField>(mesh_);
2888  printFieldInfo<pointVectorField>(mesh_);
2889  printFieldInfo<pointSphericalTensorField>(mesh_);
2890  printFieldInfo<pointSymmTensorField>(mesh_);
2891  printFieldInfo<pointTensorField>(mesh_);
2892  Pout<< nl << endl;
2893  }
2894 
2895  // Collect all maps and return
2897  (
2899  (
2900  mesh_,
2901 
2902  nOldPoints,
2903  nOldFaces,
2904  nOldCells,
2905  move(oldPatchStarts),
2906  move(oldPatchNMeshPoints),
2907 
2908  move(subPointMap),
2909  move(subFaceMap),
2910  move(subCellMap),
2911  move(subPatchMap),
2912 
2913  move(constructPointMap),
2914  move(constructFaceMap),
2915  move(constructCellMap),
2916  move(constructPatchMap),
2917 
2918  true, // subFaceMap has flip
2919  true // constructFaceMap has flip
2920  )
2921  );
2922 }
2923 
2924 
2925 // ************************************************************************* //
const labelList & oldPointMap() const
From old mesh point/face/cell to new mesh point/face/cell.
Foam::surfaceFields.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
ListType renumber(const labelUList &oldToNew, const ListType &)
Renumber the values (not the indices) of a list.
dimensionedScalar sign(const dimensionedScalar &ds)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
A HashTable with keys but without contents.
Definition: HashSet.H:59
const labelList & patchMap() const
Return patch map.
dictionary dict
static void reorderPatches(fvMesh &, const labelList &oldToNew, const label nPatches, const bool validBoundary)
Reorder and remove trailing patches. If validBoundary call is parallel.
Definition: fvMeshTools.C:140
void inplaceReorder(const labelUList &oldToNew, ListType &)
Inplace reorder the elements of a list.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
const word & name() const
Return name.
const word & name() const
Return name.
Definition: IOobject.H:295
void sortedOrder(const UList< T > &, labelList &order)
Generate the (stable) sort order for the list.
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
void finishedSends(const bool block=true)
Mark all sends as having been done. This will start receives.
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:476
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Class describing modification of a face.
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112
void setRefinement(const labelList &cellsToRemove, const labelList &facesToExpose, const labelList &patchIDs, polyTopoChange &) const
Play commands into polyTopoChange to remove cells.
Definition: removeCells.C:179
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
void addFvPatches(const List< polyPatch *> &, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:455
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
label nInternalFaces() const
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1175
static const char *const typeName
Definition: Field.H:105
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:74
void resetMotion() const
Reset motion.
Definition: polyMesh.C:1386
label nFaces() const
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:278
const labelList & faceFlipMap() const
Return face map with sign to encode flipped faces.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:248
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
Given list of cells to remove insert all the topology changes.
Definition: removeCells.H:59
label nCells() const
wordList names() const
Return the list of names of the IOobjects.
engineTime & runTime
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:309
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:446
Class containing mesh-to-mesh mapping information after a mesh addition where we add a mesh (&#39;added m...
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
Less function class that can be used for sorting processor patches.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:306
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
patches[0]
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
const labelList & addedPointMap() const
From added mesh point/face/cell to new mesh point/face/cell.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:108
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
GeometricField< sphericalTensor, fvPatchField, volMesh > volSphericalTensorField
Definition: volFieldsFwd.H:57
Combination-Reduction operation for a parallel run. The information from all nodes is collected on th...
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Various functions to operate on Lists.
Container for information needed to couple to meshes. When constructed from two meshes and a geometri...
static const pointMesh & New(const polyMesh &mesh)
Definition: MeshObject.C:44
Neighbour processor patch.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1131
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:699
fvMeshDistribute(fvMesh &mesh, const scalar mergeTol)
Construct from mesh and absolute merge tolerance.
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
Input inter-processor communications stream operating on external buffer.
Definition: UIPstream.H:53
A list of faces which address into the list of points.
virtual bool owner() const
Does this side own the patch ?
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
A List obtained as a section of another List.
Definition: SubList.H:53
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:313
bool operator()(const label a, const label b)
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:315
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: ZoneMesh.C:341
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
static void mergePoints(const polyMesh &, const Map< label > &pointToMaster, polyTopoChange &meshMod)
Helper: Merge points.
const pointZoneMesh & pointZones() const
Return point zone mesh.
Definition: polyMesh.H:470
label nOldInternalFaces() const
Number of old internal faces.
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:482
static label addPatch(fvMesh &mesh, const polyPatch &patch, const dictionary &patchFieldDict, const word &defaultPatchFieldType, const bool validBoundary)
Add patch. Inserts patch before all processor patches.
Definition: fvMeshTools.C:33
Cyclic plane patch.
wordList names() const
Return a list of patch names.
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1169
const labelList & oldFaceMap() const
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
const labelList & faceMap() const
Return face map.
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:97
const word & name() const
Return name.
Definition: zone.H:147
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1156
int neighbProcNo() const
Return neighbour processor number.
An STL-conforming hash table.
Definition: HashTable.H:61
void addZones(const List< pointZone *> &pz, const List< faceZone *> &fz, const List< cellZone *> &cz)
Add mesh zones.
Definition: polyMesh.C:953
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
static autoPtr< mapAddedPolyMesh > add(fvMesh &mesh0, const fvMesh &mesh1, const faceCoupleInfo &coupleInfo, const bool validBoundary=true)
Inplace add mesh to fvMesh. Maps all stored fields. Returns map.
Definition: fvMeshAdder.C:71
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
Foam::polyBoundaryMesh.
A packed storage unstructured matrix of objects of type <T> using an offset table for access...
void setSize(const label mRows)
Reset size of CompactListList.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
Output inter-processor communications stream operating on external buffer.
Definition: UOPstream.H:54
bool checkParallelSync(const bool report=false) const
Check whether all procs have all patches and in same order. Return.
static const char nl
Definition: Ostream.H:265
defineTypeNameAndDebug(combustionModel, 0)
autoPtr< mapDistributePolyMesh > distribute(const labelList &dist)
Send cells to neighbours according to distribution.
Post-processing mesh subset tool. Given the original mesh and the list of selected cells...
Definition: fvMeshSubset.H:73
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
A subset of mesh cells.
Definition: cellZone.H:61
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:32
This boundary condition enables processor communication across cyclic patches.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
void reverseDistribute(const label constructSize, List< T > &, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Reverse distribute data using default commsType.
const vectorField & faceCentres() const
label constructSize() const
Constructed data size.
Calculates points shared by more than two processor patches or cyclic patches.
Definition: globalPoints.H:100
void setSize(const label)
Reset size of List.
Definition: List.C:281
static void printMeshInfo(const fvMesh &)
Print some info on mesh.
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
label patchi
Class containing processor-to-processor mapping information.
const fvMesh & subMesh() const
Return reference to subset mesh.
A subset of mesh points. The labels of points in the zone can be obtained from the addressing() list...
Definition: pointZone.H:62
label newPointi
Definition: readKivaGrid.H:501
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
Foam::fvBoundaryMesh.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
static void printCoupleInfo(const primitiveMesh &, const labelList &, const labelList &, const labelList &, const labelList &)
Print some info on coupling data.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A List with indirect addressing.
Definition: fvMatrix.H:106
Direct mesh changes based on v1.3 polyTopoChange syntax.
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:303
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
lessProcPatches(const labelList &nbrProc, const labelList &referPatchID)
wordList names() const
Return a list of zone names.
Definition: ZoneMesh.C:256
label n
virtual transformType transform() const
Type of transform.
label nPoints() const
void resize(const label newSize)
Resize the hash table for efficiency.
Definition: HashTable.C:432
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return a pointer to a new patch created on freestore from.
Definition: polyPatchNew.C:32
This boundary condition enables processor communication across patches.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
labelList getExposedFaces(const labelList &cellsToRemove) const
Get labels of exposed faces.
Definition: removeCells.C:73
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled boundary face values.
Definition: syncTools.H:430
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:202
autoPtr< mapPolyMesh > changeMesh(polyMesh &mesh, const bool inflate, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
label size() const
Return the number of elements in the UList.
Definition: UListI.H:299
static labelList countCells(const labelList &)
Helper function: count cells per processor in wanted distribution.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
bool found
const labelList & addedFaceMap() const
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
void setLargeCellSubset(const labelList &region, const label currentRegion, const label patchID=-1, const bool syncCouples=true)
Set the subset from all cells with region == currentRegion.
Definition: fvMeshSubset.C:885
const labelList & cellMap() const
Return cell map.
const labelList & pointMap() const
Return point map.
Namespace for OpenFOAM.
virtual bool owner() const
Does the processor own the patch ?
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:540