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