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