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