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