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