fvMeshTools.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2012-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "fvMeshTools.H"
30 #include "pointSet.H"
31 #include "faceSet.H"
32 #include "cellSet.H"
33 #include "fileOperation.H"
34 #include "BitOps.H"
35 #include "IOobjectList.H"
36 #include "basicFvGeometryScheme.H"
37 #include "processorPolyPatch.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 // Adds patch if not yet there. Returns patchID.
45 (
46  fvMesh& mesh,
47  const polyPatch& patch,
48  const dictionary& patchFieldDict,
49  const word& defaultPatchFieldType,
50  const bool validBoundary
51 )
52 {
53  auto& polyPatches = const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
54 
55  label patchi = polyPatches.findPatchID(patch.name());
56  if (patchi != -1)
57  {
58  // Already there
59  return patchi;
60  }
61 
62 
63  // Append at end unless there are processor patches
64  label insertPatchi = polyPatches.size();
65  label startFacei = mesh.nFaces();
66 
67  if (!isA<processorPolyPatch>(patch))
68  {
69  forAll(polyPatches, patchi)
70  {
71  const polyPatch& pp = polyPatches[patchi];
72 
73  if (isA<processorPolyPatch>(pp))
74  {
75  insertPatchi = patchi;
76  startFacei = pp.start();
77  break;
78  }
79  }
80  }
81 
82 
83  // Below is all quite a hack. Feel free to change once there is a better
84  // mechanism to insert and reorder patches.
85 
86  // Clear local fields and e.g. polyMesh parallelInfo.
87  mesh.clearOut();
88 
89  label sz = polyPatches.size();
90 
91  auto& fvPatches = const_cast<fvBoundaryMesh&>(mesh.boundary());
92 
93  // Add polyPatch at the end
94  polyPatches.resize(sz+1);
95  polyPatches.set
96  (
97  sz,
98  patch.clone
99  (
100  polyPatches,
101  insertPatchi, //index
102  0, //size
103  startFacei //start
104  )
105  );
106  fvPatches.resize(sz+1);
107  fvPatches.set
108  (
109  sz,
111  (
112  polyPatches[sz], // point to newly added polyPatch
113  mesh.boundary()
114  )
115  );
116 
117 
118  do
119  {
120  #undef doLocalCode
121  #define doLocalCode(FieldType) \
122  { \
123  addPatchFields<FieldType> \
124  ( \
125  mesh, patchFieldDict, defaultPatchFieldType, Zero \
126  ); \
127  }
128 
129  // Volume fields
135 
136  // Surface fields
142 
143  #undef doLocalCode
144  }
145  while (false);
146 
147 
148  // Create reordering list
149  // - patches before insert position stay as is
150  // - patches after insert position move one up
151  // - appended patch gets moved to insert position
152 
153  labelList oldToNew(identity(sz+1));
154  for (label i = insertPatchi; i < sz; ++i)
155  {
156  oldToNew[i] = i+1;
157  }
158  oldToNew[sz] = insertPatchi;
159 
160 
161  // Shuffle into place
162  polyPatches.reorder(oldToNew, validBoundary);
163  fvPatches.reorder(oldToNew);
164 
165  do
166  {
167  #undef doLocalCode
168  #define doLocalCode(FieldType) \
169  { \
170  reorderPatchFields<FieldType>(mesh, oldToNew); \
171  }
172 
173  // Volume fields
179 
180  // Surface fields
186 
187  #undef doLocalCode
188  }
189  while (false);
190 
191  return insertPatchi;
192 }
193 
194 
195 void Foam::fvMeshTools::setPatchFields
196 (
197  fvMesh& mesh,
198  const label patchi,
199  const dictionary& patchFieldDict
200 )
201 {
202  do
203  {
204  #undef doLocalCode
205  #define doLocalCode(FieldType) \
206  { \
207  setPatchFields<FieldType>(mesh, patchi, patchFieldDict); \
208  }
209 
210  // Volume fields
216 
217  // Surface fields
224  #undef doLocalCode
225  }
226  while (false);
227 }
228 
229 
230 void Foam::fvMeshTools::zeroPatchFields(fvMesh& mesh, const label patchi)
231 {
232  do
233  {
234  #undef doLocalCode
235  #define doLocalCode(FieldType) \
236  { \
237  setPatchFields<FieldType>(mesh, patchi, Zero); \
238  }
239 
240  // Volume fields
246 
247  // Surface fields
253 
254  #undef doLocalCode
255  }
256  while (false);
257 }
258 
259 
260 // Deletes last patch
261 void Foam::fvMeshTools::trimPatches(fvMesh& mesh, const label nPatches)
262 {
263  // Clear local fields and e.g. polyMesh globalMeshData.
264  mesh.clearOut();
265 
266  auto& polyPatches = const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
267  auto& fvPatches = const_cast<fvBoundaryMesh&>(mesh.boundary());
268 
269  if (polyPatches.empty())
270  {
272  << "No patches in mesh"
273  << abort(FatalError);
274  }
275 
276  label nFaces = 0;
277  for (label patchi = nPatches; patchi < polyPatches.size(); patchi++)
278  {
279  nFaces += polyPatches[patchi].size();
280  }
281  reduce(nFaces, sumOp<label>());
282 
283  if (nFaces)
284  {
286  << "There are still " << nFaces
287  << " faces in " << polyPatches.size()-nPatches
288  << " patches to be deleted" << abort(FatalError);
289  }
290 
291  // Remove actual patches
292  polyPatches.resize(nPatches);
293  fvPatches.resize(nPatches);
294 
295  do
296  {
297  #undef doLocalCode
298  #define doLocalCode(FieldType) \
299  { \
300  trimPatchFields<FieldType>(mesh, nPatches); \
301  }
302 
303  // Volume fields
309 
310  // Surface fields
316 
317  #undef doLocalCode
318  }
319  while (false);
320 }
321 
322 
324 (
325  fvMesh& mesh,
326  const labelList& oldToNew,
327  const label nNewPatches,
328  const bool validBoundary
329 )
330 {
331  auto& polyPatches = const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
332  auto& fvPatches = const_cast<fvBoundaryMesh&>(mesh.boundary());
333 
334  // Shuffle into place
335  polyPatches.reorder(oldToNew, validBoundary);
336  fvPatches.reorder(oldToNew);
337 
338  do
339  {
340  #undef doLocalCode
341  #define doLocalCode(FieldType) \
342  { \
343  reorderPatchFields<FieldType>(mesh, oldToNew); \
344  }
345 
346  // Volume fields
352 
353  // Surface fields
359 
360  #undef doLocalCode
361  }
362  while (false);
364  // Remove last.
365  trimPatches(mesh, nNewPatches);
366 }
367 
368 
370 (
371  fvMesh& mesh,
372  const bool validBoundary
373 )
374 {
376 
377  labelList newToOld(pbm.size());
378  labelList oldToNew(pbm.size(), -1);
379  label newI = 0;
380 
381 
382  // Assumes all non-coupled boundaries are on all processors!
383  forAll(pbm, patchI)
384  {
385  const polyPatch& pp = pbm[patchI];
386 
387  if (!isA<processorPolyPatch>(pp))
388  {
389  label nFaces = pp.size();
390  if (validBoundary)
391  {
392  reduce(nFaces, sumOp<label>());
393  }
394 
395  if (nFaces > 0)
396  {
397  newToOld[newI] = patchI;
398  oldToNew[patchI] = newI++;
399  }
400  }
401  }
402 
403  // Same for processor patches (but need no reduction)
404  forAll(pbm, patchI)
405  {
406  const polyPatch& pp = pbm[patchI];
407 
408  if (isA<processorPolyPatch>(pp) && pp.size())
409  {
410  newToOld[newI] = patchI;
411  oldToNew[patchI] = newI++;
412  }
413  }
414 
415  newToOld.resize(newI);
416 
417  // Move all deletable patches to the end
418  forAll(oldToNew, patchI)
419  {
420  if (oldToNew[patchI] == -1)
421  {
422  oldToNew[patchI] = newI++;
423  }
424  }
426  reorderPatches(mesh, oldToNew, newToOld.size(), validBoundary);
427 
428  return newToOld;
429 }
430 
431 
433 {
434  // Set the fvGeometryScheme to basic since it does not require
435  // any parallel communication to construct the geometry. During
436  // redistributePar one might temporarily end up with processors
437  // with zero procBoundaries. Normally procBoundaries trigger geometry
438  // calculation (e.g. send over cellCentres) so on the processors without
439  // procBoundaries this will not happen. The call to the geometry calculation
440  // is not synchronised and this might lead to a hang for geometry schemes
441  // that do require synchronisation
442 
443  tmp<fvGeometryScheme> basicGeometry
444  (
446  );
447  mesh.geometry(basicGeometry);
448 }
449 
450 
453 (
454  const IOobject& io,
455  const bool masterOnlyReading,
456  const bool verbose
457 )
458 {
459  // Region name
460  // ~~~~~~~~~~~
461 
462  const fileName meshSubDir
463  (
465  );
466 
467 
468  fileName facesInstance;
469  fileName pointsInstance;
470 
471  // Patch types
472  // ~~~~~~~~~~~
473  // Read and scatter master patches (without reading master mesh!)
474 
475  PtrList<entry> patchEntries;
476  if (UPstream::master())
477  {
478  const bool oldParRun = Pstream::parRun(false);
479 
480  facesInstance = io.time().findInstance
481  (
482  meshSubDir,
483  "faces",
485  );
486  pointsInstance = io.time().findInstance
487  (
488  meshSubDir,
489  "points",
491  );
492 
493  patchEntries = polyBoundaryMeshEntries
494  (
495  IOobject
496  (
497  "boundary",
498  facesInstance,
499  meshSubDir,
500  io.db(),
504  )
505  );
506 
507  UPstream::parRun(oldParRun);
508  }
509 
510  // Broadcast information to all
512  (
514  patchEntries,
515  facesInstance,
516  pointsInstance
517  );
518 
519 
520  // Dummy meshes
521  // ~~~~~~~~~~~~
522 
523  // Set up to read-if-present.
524  // Note: does not search for mesh so set instance explicitly
525 
526  IOobject meshIO(io);
527  meshIO.instance() = facesInstance;
528  meshIO.readOpt(IOobject::READ_IF_PRESENT);
529 
530  // For mesh components (points, faces, ...)
531  IOobject cmptIO(meshIO, "points", meshSubDir);
532  cmptIO.readOpt(IOobject::MUST_READ);
533  cmptIO.writeOpt(IOobject::NO_WRITE);
534  cmptIO.registerObject(false);
535 
536 
537  // Check who has a mesh
538  const fileName meshDir = io.time().path()/facesInstance/meshSubDir;
539  bool haveMesh = isDir(meshDir);
540  if (masterOnlyReading && !Pstream::master())
541  {
542  haveMesh = false;
543  meshIO.readOpt(IOobject::NO_READ);
544  }
545 
546  if (!haveMesh)
547  {
548  cmptIO.readOpt(IOobject::NO_READ);
549  }
550 
551 
552  // Read mesh
553  // ~~~~~~~~~
554  // Now all processors use supplied points,faces etc
555  // Note: solution, schemes are also using the supplied IOobject so
556  // on slave will be NO_READ, on master READ_IF_PRESENT. This will
557  // conflict with e.g. timeStampMaster reading so switch off.
558  // Note: v2006 used the READ_IF_PRESENT flag in the meshIO.readOpt(). v2012
559  // (correctly) does no longer so below code explicitly addFvPatches
560  // using the separately read boundary file.
561 
562  const auto oldCheckType = IOobject::fileModificationChecking;
564 
565 
566  // Points
567  cmptIO.instance() = pointsInstance;
568  cmptIO.rename("points");
569  pointIOField points(cmptIO);
570 
571  // All other mesh components use facesInstance ...
572  cmptIO.instance() = facesInstance;
573 
574  // Faces
575  cmptIO.rename("faces");
576  faceCompactIOList faces(cmptIO);
577 
578  // Face owner
579  cmptIO.rename("owner");
580  labelIOList owner(cmptIO);
581 
582  // Face neighbour
583  cmptIO.rename("neighbour");
584  labelIOList neighbour(cmptIO);
585 
586 
588  (
589  meshIO,
590  std::move(points),
591  std::move(faces),
592  std::move(owner),
593  std::move(neighbour)
594  );
595  fvMesh& mesh = *meshPtr;
596 
597  IOobject::fileModificationChecking = oldCheckType;
598 
599 
600  // Some processors without patches? - add patches
601 
603  {
604  DynamicList<polyPatch*> newPatches(patchEntries.size());
605 
606  if (mesh.boundary().size() == patchEntries.size())
607  {
608  // Assumably we have correctly read the boundary and can clone.
609  // Note: for
610  // v2012 onwards this is probably never the case and this whole
611  // section can be removed.
612  forAll(mesh.boundary(), patchi)
613  {
614  newPatches.append
615  (
616  mesh.boundaryMesh()[patchi].clone(mesh.boundaryMesh()).ptr()
617  );
618  }
619  }
620  else
621  {
622  // Use patchEntries, which were read on master and broadcast
623  label nPatches = 0;
624 
625  const bool isEmptyMesh = (mesh.nInternalFaces() == 0);
626 
627  forAll(patchEntries, patchi)
628  {
629  const entry& e = patchEntries[patchi];
630  const word type(e.dict().get<word>("type"));
631  const word& name = e.keyword();
632 
633  if
634  (
635  type == processorPolyPatch::typeName
636  || type == processorCyclicPolyPatch::typeName
637  )
638  {
639  // Unlikely to work with inter-mixed proc-patches anyhow
640  // - could break out of loop here
641  }
642  else
643  {
644  dictionary patchDict(e.dict());
645 
646  if (isEmptyMesh)
647  {
648  patchDict.set("nFaces", 0);
649  patchDict.set("startFace", 0);
650  }
651 
652  newPatches.append
653  (
655  (
656  name,
657  patchDict,
658  nPatches++,
660  ).ptr()
661  );
662  }
663  }
664  }
665 
667  mesh.addFvPatches(newPatches);
668  }
669 
670 
671  // Determine zones
672  // ~~~~~~~~~~~~~~~
673 
674  wordList pointZoneNames(mesh.pointZones().names());
675  wordList faceZoneNames(mesh.faceZones().names());
676  wordList cellZoneNames(mesh.cellZones().names());
677 
679  (
681  pointZoneNames,
682  faceZoneNames,
683  cellZoneNames
684  );
685 
686  if (!haveMesh)
687  {
688  // Add the zones. Make sure to remove the old dummy ones first
689  mesh.pointZones().clear();
690  mesh.faceZones().clear();
691  mesh.cellZones().clear();
692 
693  List<pointZone*> pz(pointZoneNames.size());
694  forAll(pointZoneNames, i)
695  {
696  pz[i] = new pointZone
697  (
698  pointZoneNames[i],
699  i,
700  mesh.pointZones()
701  );
702  }
703  List<faceZone*> fz(faceZoneNames.size());
704  forAll(faceZoneNames, i)
705  {
706  fz[i] = new faceZone
707  (
708  faceZoneNames[i],
709  i,
710  mesh.faceZones()
711  );
712  }
713  List<cellZone*> cz(cellZoneNames.size());
714  forAll(cellZoneNames, i)
715  {
716  cz[i] = new cellZone
717  (
718  cellZoneNames[i],
719  i,
720  mesh.cellZones()
721  );
722  }
723 
724  if (pz.size() || fz.size() || cz.size())
725  {
726  mesh.addZones(pz, fz, cz);
727  }
728  }
729 
730  return meshPtr;
731 }
732 
733 
735 Foam::fvMeshTools::loadOrCreateMeshImpl
736 (
737  const IOobject& io,
738  refPtr<fileOperation>* readHandlerPtr, // Can be nullptr
739  const bool decompose,
740  const bool verbose
741 )
742 {
743  // Region name
744  // ~~~~~~~~~~~
745 
746  const fileName meshSubDir
747  (
749  );
750 
751 
752  // Patch types
753  // ~~~~~~~~~~~
754  // Read and scatter master patches (without reading master mesh!)
755 
756  PtrList<entry> patchEntries;
757  if (UPstream::master())
758  {
759  const bool oldParRun = UPstream::parRun(false);
760  const label oldNumProcs = fileHandler().nProcs();
761  const int oldCache = fileOperation::cacheLevel(0);
762 
763  const fileName facesInstance = io.time().findInstance
764  (
765  meshSubDir,
766  "faces",
768  );
769 
770  patchEntries = polyBoundaryMeshEntries
771  (
772  IOobject
773  (
774  "boundary",
775  facesInstance,
776  meshSubDir,
777  io.db(),
781  )
782  );
783 
784  fileOperation::cacheLevel(oldCache);
785  if (oldParRun)
786  {
787  const_cast<fileOperation&>(fileHandler()).nProcs(oldNumProcs);
788  }
789  UPstream::parRun(oldParRun);
790  }
791 
792  // Broadcast: send patches to all
794 
795 
796  // Check who has or needs a mesh.
797  bool haveLocalMesh = false;
798 
799  if (readHandlerPtr)
800  {
801  // Non-null reference when a mesh exists on given processor
802  haveLocalMesh = (*readHandlerPtr).good();
803  }
804 
805  // Globally consistent information about who has a mesh
806  boolList haveMesh
807  (
808  UPstream::allGatherValues<bool>(haveLocalMesh)
809  );
810 
811  // Make sure all have the same number of processors
812  label masterNProcs = fileHandler().nProcs();
813  Pstream::broadcast(masterNProcs);
814  const_cast<fileOperation&>(fileHandler()).nProcs(masterNProcs);
815 
816 
817  autoPtr<fvMesh> meshPtr;
818 
819  if (!haveLocalMesh)
820  {
821  // No local mesh - need to synthesize one
822 
823  const bool oldParRun = UPstream::parRun(false);
824  const label oldNumProcs = fileHandler().nProcs();
825  const int oldCache = fileOperation::cacheLevel(0);
826 
827  // Create dummy mesh - on procs that don't already have a mesh
828  meshPtr.reset
829  (
830  new fvMesh
831  (
833  Foam::zero{},
834  false
835  )
836  );
837  fvMesh& mesh = *meshPtr;
838 
839  // Add patches
840  polyPatchList patches(patchEntries.size());
841  label nPatches = 0;
842 
843  forAll(patchEntries, patchi)
844  {
845  const entry& e = patchEntries[patchi];
846  const word type(e.dict().get<word>("type"));
847  const word& name = e.keyword();
848 
849  if
850  (
851  type == processorPolyPatch::typeName
852  || type == processorCyclicPolyPatch::typeName
853  )
854  {
855  // Stop at the first processor patch.
856  // - logic fails with inter-mixed proc-patches anyhow
857  break;
858  }
859  else
860  {
861  dictionary patchDict(e.dict());
862  patchDict.set("nFaces", 0);
863  patchDict.set("startFace", 0);
864 
865  patches.set
866  (
867  patchi,
869  (
870  name,
871  patchDict,
872  nPatches++,
874  )
875  );
876  }
877  }
879  mesh.addFvPatches(patches, false); // No parallel comms
880 
881  if (!readHandlerPtr)
882  {
883  // The 'old' way of doing things.
884  // Write the dummy mesh to disk for subsequent re-reading.
885  //
886  // This is not particularly elegant.
887  //
888  // Note: add some dummy zones so upon reading it does not read them
889  // from the undecomposed case. Should be done as extra argument to
890  // regIOobject::readStream?
891 
892  List<pointZone*> pz
893  (
894  1,
895  new pointZone("dummyZone", 0, mesh.pointZones())
896  );
897  List<faceZone*> fz
898  (
899  1,
900  new faceZone("dummyZone", 0, mesh.faceZones())
901  );
902  List<cellZone*> cz
903  (
904  1,
905  new cellZone("dummyZone", 0, mesh.cellZones())
906  );
907  mesh.addZones(pz, fz, cz);
908  mesh.pointZones().clear();
909  mesh.faceZones().clear();
910  mesh.cellZones().clear();
911  //Pout<< "Writing dummy mesh: " << mesh.polyMesh::objectPath()
912  // << endl;
913  mesh.write();
914 
915  // Discard - it will be re-read later
916  meshPtr.reset(nullptr);
917  }
918 
919  fileOperation::cacheLevel(oldCache);
920  if (oldParRun)
921  {
922  const_cast<fileOperation&>(fileHandler()).nProcs(oldNumProcs);
923  }
924  UPstream::parRun(oldParRun); // Restore parallel state
925  }
926  else if (readHandlerPtr && haveLocalMesh)
927  {
928  const labelList meshProcIds(BitOps::sortedToc(haveMesh));
929 
930  UPstream::communicator newCommunicator;
931  const label oldWorldComm = UPstream::commWorld();
932 
933  auto& readHandler = *readHandlerPtr;
934  auto oldHandler = fileOperation::fileHandler(readHandler);
935 
936  // With IO ranks the communicator of the fileOperation will
937  // only include the ranks for the current IO rank.
938  // Instead allocate a new communicator for everyone with a mesh
939 
940  const auto& handlerProcIds = UPstream::procID(fileHandler().comm());
941 
942  // Comparing global ranks in the communicator.
943  // Use std::equal for the List<label> vs List<int> comparison
944 
945  if
946  (
947  meshProcIds.size() == handlerProcIds.size()
948  && std::equal
949  (
950  meshProcIds.cbegin(),
951  meshProcIds.cend(),
952  handlerProcIds.cbegin()
953  )
954  )
955  {
956  const_cast<fileOperation&>(fileHandler()).nProcs(UPstream::nProcs(oldWorldComm));
957  // Can use the handler communicator as is.
959  }
960  else if
961  (
962  UPstream::nProcs(fileHandler().comm())
964  )
965  {
966  // Need a new communicator for the fileHandler.
967 
968  // Warning: MS-MPI currently uses MPI_Comm_create() instead of
969  // MPI_Comm_create_group() so it will block here!
970 
971  newCommunicator.reset(UPstream::worldComm, meshProcIds);
972  UPstream::commWorld(newCommunicator.comm());
973  }
974 
975  // Load but do not initialise
976  meshPtr = autoPtr<fvMesh>::New(io, false);
977 
978  readHandler = fileOperation::fileHandler(oldHandler);
979  UPstream::commWorld(oldWorldComm);
980 
981  // Reset mesh communicator to the real world comm
982  meshPtr().polyMesh::comm() = UPstream::commWorld();
983  }
984 
985 
986  if (!meshPtr)
987  {
988  // Using the 'old' way of doing things (writing to disk and re-reading).
989 
990  // Read mesh from disk
991  //
992  // Now all processors have a (possibly zero size) mesh so can
993  // read in parallel
994 
995  //Pout<< "Reading mesh from " << io.objectRelPath() << endl;
996  // Load but do not initialise
997  meshPtr = autoPtr<fvMesh>::New(io, false);
998  }
999 
1000  fvMesh& mesh = meshPtr();
1001 
1002 
1003  // Check patches
1004  // ~~~~~~~~~~~~~
1005 
1006  #if 0
1007  if (!UPstream::master() && haveLocalMesh)
1008  {
1009  // Check master patch entries names against local ones
1010 
1011  const polyBoundaryMesh& patches = mesh.boundaryMesh();
1012 
1013  forAll(patchEntries, patchi)
1014  {
1015  const entry& e = patchEntries[patchi];
1016  const word type(e.dict().get<word>("type"));
1017  const word& name = e.keyword();
1018 
1019  if
1020  (
1021  type == processorPolyPatch::typeName
1022  || type == processorCyclicPolyPatch::typeName
1023  )
1024  {
1025  break;
1026  }
1027 
1028  if (patchi >= patches.size())
1029  {
1031  << "Non-processor patches not synchronised." << endl
1032  << "Processor " << UPstream::myProcNo()
1033  << " has only " << patches.size()
1034  << " patches, master has "
1035  << patchi << endl
1036  << exit(FatalError);
1037  }
1038 
1039  if
1040  (
1041  type != patches[patchi].type()
1042  || name != patches[patchi].name()
1043  )
1044  {
1046  << "Non-processor patches not synchronised." << endl
1047  << "Master patch " << patchi
1048  << " name:" << type
1049  << " type:" << type << endl
1050  << "Processor " << UPstream::myProcNo()
1051  << " patch " << patchi
1052  << " has name:" << patches[patchi].name()
1053  << " type:" << patches[patchi].type()
1054  << exit(FatalError);
1055  }
1056  }
1057  }
1058  #endif
1059 
1060 
1061  // Synchronize zones
1062  // ~~~~~~~~~~~~~~~~~
1063 
1064  {
1065  wordList pointZoneNames(mesh.pointZones().names());
1066  wordList faceZoneNames(mesh.faceZones().names());
1067  wordList cellZoneNames(mesh.cellZones().names());
1069  (
1071  pointZoneNames,
1072  faceZoneNames,
1073  cellZoneNames
1074  );
1075 
1076 
1077  if (!haveLocalMesh)
1078  {
1079  // Add the zones. Make sure to remove the old dummy ones first
1080  mesh.pointZones().clear();
1081  mesh.faceZones().clear();
1082  mesh.cellZones().clear();
1083 
1084  PtrList<pointZone> pz(pointZoneNames.size());
1085  forAll(pointZoneNames, i)
1086  {
1087  pz.emplace_set(i, pointZoneNames[i], i, mesh.pointZones());
1088  }
1089 
1090  PtrList<faceZone> fz(faceZoneNames.size());
1091  forAll(faceZoneNames, i)
1092  {
1093  fz.emplace_set(i, faceZoneNames[i], i, mesh.faceZones());
1094  }
1095 
1096  PtrList<cellZone> cz(cellZoneNames.size());
1097  forAll(cellZoneNames, i)
1098  {
1099  cz.emplace_set(i, cellZoneNames[i], i, mesh.cellZones());
1100  }
1101  mesh.addZones(std::move(pz), std::move(fz), std::move(cz));
1102  }
1103  }
1104 
1105 
1106  // Synchronize sets (on disk)
1107  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
1108 
1109  if (!readHandlerPtr)
1110  {
1111  wordList pointSetNames;
1112  wordList faceSetNames;
1113  wordList cellSetNames;
1114  if (UPstream::master())
1115  {
1116  // Read sets
1117  const bool oldParRun = UPstream::parRun(false);
1118 
1119  IOobjectList objects(mesh, mesh.facesInstance(), "polyMesh/sets");
1120  UPstream::parRun(oldParRun);
1121 
1122  pointSetNames = objects.sortedNames<pointSet>();
1123  faceSetNames = objects.sortedNames<faceSet>();
1124  cellSetNames = objects.sortedNames<cellSet>();
1125  }
1127  (
1129  pointSetNames,
1130  faceSetNames,
1131  cellSetNames
1132  );
1133 
1134  if (!haveLocalMesh)
1135  {
1136  for (const word& setName : pointSetNames)
1137  {
1138  pointSet(mesh, setName, 0).write();
1139  }
1140  for (const word& setName : faceSetNames)
1141  {
1142  faceSet(mesh, setName, 0).write();
1143  }
1144  for (const word& setName : cellSetNames)
1145  {
1146  cellSet(mesh, setName, 0).write();
1147  }
1148  }
1149  }
1150 
1151 
1152  // Meshes have been done without init so now can do full initialisation
1153 
1155  mesh.init(true);
1156 
1157  // Force early recreation of globalMeshData etc
1158  mesh.globalData();
1159  mesh.tetBasePtIs();
1160  mesh.geometricD();
1161 
1162 
1163  // Do some checks.
1164 
1165  // Check if the boundary definition is unique
1166  // and processor patches are correct
1167  mesh.boundaryMesh().checkDefinition(verbose);
1168  mesh.boundaryMesh().checkParallelSync(verbose);
1169 
1170  // Check names of zones are equal
1171  mesh.cellZones().checkDefinition(verbose);
1172  mesh.cellZones().checkParallelSync(verbose);
1173  mesh.faceZones().checkDefinition(verbose);
1174  mesh.faceZones().checkParallelSync(verbose);
1175  mesh.pointZones().checkDefinition(verbose);
1176  mesh.pointZones().checkParallelSync(verbose);
1178  return meshPtr;
1179 }
1180 
1181 
1184 (
1185  const IOobject& io,
1186  const bool decompose,
1187  const bool verbose
1188 )
1189 {
1190  return fvMeshTools::loadOrCreateMeshImpl
1191  (
1192  io,
1193  nullptr, // fileOperation (ignore)
1194  decompose,
1195  verbose
1196  );
1197 }
1198 
1199 
1202 (
1203  const IOobject& io,
1204  refPtr<fileOperation>& readHandler,
1205  const bool verbose
1206 )
1207 {
1208  return fvMeshTools::loadOrCreateMeshImpl
1209  (
1210  io,
1211  &readHandler,
1212  false, // decompose (ignored)
1213  verbose
1214  );
1215 }
1216 
1217 
1219 (
1220  const objectRegistry& mesh,
1221  const word& regionName,
1222  const bool verbose
1223 )
1224 {
1225  // Create dummy system/fv*
1226  {
1227  IOobject io
1228  (
1229  "fvSchemes",
1230  mesh.time().system(),
1231  regionName,
1232  mesh.thisDb(),
1236  );
1237 
1238  if (!io.typeHeaderOk<IOdictionary>(false))
1239  {
1240  if (verbose)
1241  {
1242  Info<< "Writing dummy " << regionName/io.name() << endl;
1243  }
1244  dictionary dict;
1245  dict.add("divSchemes", dictionary());
1246  dict.add("gradSchemes", dictionary());
1247  dict.add("laplacianSchemes", dictionary());
1248 
1249  IOdictionary(io, dict).regIOobject::write();
1250  }
1251  }
1252  {
1253  IOobject io
1254  (
1255  "fvSolution",
1256  mesh.time().system(),
1257  regionName,
1258  mesh.thisDb(),
1262  );
1263 
1264  if (!io.typeHeaderOk<IOdictionary>(false))
1265  {
1266  if (verbose)
1267  {
1268  Info<< "Writing dummy " << regionName/io.name() << endl;
1269  }
1270  dictionary dict;
1271  IOdictionary(io, dict).regIOobject::write();
1272  }
1273  }
1274 }
1275 
1276 
1277 // ************************************************************************* //
label nPatches
Definition: readKivaGrid.H:396
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
bool checkParallelSync(const bool report=false) const
Check whether all procs have all zones and in same order.
Definition: ZoneMesh.C:869
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:85
const polyBoundaryMesh & pbm
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:227
dictionary dict
static void reorderPatches(fvMesh &, const labelList &oldToNew, const label nPatches, const bool validBoundary)
Reorder and remove trailing patches.
Definition: fvMeshTools.C:317
fileName path() const
Return path = rootPath/caseName. Same as TimePaths::path()
Definition: Time.H:503
void removeFvBoundary()
Remove boundary patches. Warning: fvPatchFields hold ref to these fvPatches.
Definition: fvMesh.C:692
A class for handling file names.
Definition: fileName.H:72
void clear()
Clear the zones.
Definition: ZoneMesh.C:841
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:859
word findInstance(const fileName &dir, const word &name=word::null, IOobjectOption::readOption rOpt=IOobjectOption::MUST_READ, const word &stopInstance=word::null) const
Return time instance (location) of dir that contains the file name (eg, used in reading mesh data)...
Definition: Time.C:725
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:899
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
vectorIOField pointIOField
pointIOField is a vectorIOField.
Definition: pointIOField.H:38
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:86
PtrList< T > clone(Args &&... args) const
Make a copy by cloning each of the list elements.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
GeometricField< sphericalTensor, fvsPatchField, surfaceMesh > surfaceSphericalTensorField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
refPtr< fileOperation > fileHandler(std::nullptr_t)
Delete current file handler - forwards to fileOperation::handler()
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
Definition: dynamicFvMesh.C:84
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:625
Ignore writing from objectRegistry::writeObject()
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
static label worldComm
Communicator for all ranks. May differ from commGlobal() if local worlds are in use.
Definition: UPstream.H:409
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
GeometricField< sphericalTensor, fvPatchField, volMesh > volSphericalTensorField
Definition: volFieldsFwd.H:84
A class for managing references or pointers (no reference counting)
Definition: HashPtrTable.H:49
bool equal(const Matrix< Form1, Type > &A, const Matrix< Form2, Type > &B, const bool verbose=false, const label maxDiffs=10, const scalar relTol=1e-5, const scalar absTol=1e-8)
Compare matrix elements for absolute or relative equality.
Definition: MatrixTools.C:27
static void broadcast(Type &value, const label comm=UPstream::worldComm)
Broadcast content (contiguous or non-contiguous) to all communicator ranks. Does nothing in non-paral...
Default geometry calculation scheme. Slight stabilisation for bad meshes.
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:376
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:865
static void zeroPatchFields(fvMesh &mesh, const label patchi)
Change patchField to zero on registered fields.
Definition: fvMeshTools.C:223
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
bool isDir(const fileName &name, const bool followLink=true)
Does the name exist as a DIRECTORY in the file system?
Definition: POSIX.C:860
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition: polyPatch.H:56
label size() const noexcept
The number of elements in table.
Definition: HashTable.H:342
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
CompactIOList< face, label > faceCompactIOList
Compact IO for a List of face.
Definition: faceIOList.H:35
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run. ...
Definition: UPstream.H:1065
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
IOList< label > labelIOList
IO for a List of label.
Definition: labelIOList.H:32
dynamicFvMesh & mesh
GeometricField< symmTensor, fvsPatchField, surfaceMesh > surfaceSymmTensorField
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const pointField & points
Foam::autoPtr< Foam::dynamicFvMesh > meshPtr
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
A class for handling words, derived from Foam::string.
Definition: word.H:63
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:38
bool checkDefinition(const bool report=false) const
Check boundary definition.
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
const word & system() const noexcept
Return system name.
Definition: TimePathsI.H:118
void addZones(PtrList< pointZone > &&pz, PtrList< faceZone > &&fz, PtrList< cellZone > &&cz)
Add mesh zones.
Definition: polyMesh.C:1009
GeometricField< tensor, fvsPatchField, surfaceMesh > surfaceTensorField
Reading is optional [identical to LAZY_READ].
static autoPtr< fvMesh > newMesh(const IOobject &io, const bool masterOnlyReading, const bool verbose=false)
Read mesh or create dummy mesh (0 cells, >0 patches).
Definition: fvMeshTools.C:446
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:450
const globalMeshData & globalData() const
Return parallel info (demand-driven)
Definition: polyMesh.C:1311
label nInternalFaces() const noexcept
Number of internal faces.
static void createDummyFvMeshFiles(const objectRegistry &parent, const word &regionName, const bool verbose=false)
Create additional fvSchemes/fvSolution files.
Definition: fvMeshTools.C:1212
static autoPtr< fvMesh > loadOrCreateMesh(const IOobject &io, const bool decompose, const bool verbose=false)
Read mesh if available, or create empty mesh with non-proc as per proc0 mesh.
Definition: fvMeshTools.C:1177
errorManip< error > abort(error &err)
Definition: errorManip.H:139
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: PtrList.H:159
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1102
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (respects is_globalIOobject trait) and check its info.
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition: IOobject.C:456
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
bool checkParallelSync(const bool report=false) const
Check whether all procs have all patches and in same order.
static int cacheLevel() noexcept
Return cache level.
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:670
static fileCheckTypes fileModificationChecking
Type of file modification checking.
Definition: IOobject.H:343
void addFvPatches(polyPatchList &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:662
static label commWorld() noexcept
Communicator for all ranks (respecting any local worlds)
Definition: UPstream.H:429
void resize(const label newLen)
Adjust size of PtrList.
Definition: PtrList.C:95
static void broadcasts(const label comm, Type &arg1, Args &&... args)
Broadcast multiple items to all communicator ranks. Does nothing in non-parallel. ...
bool empty() const noexcept
True if the list is empty (ie, size() is zero)
Definition: UPtrListI.H:99
const pointZoneMesh & pointZones() const noexcept
Return point zone mesh.
Definition: polyMesh.H:662
List< word > wordList
List of word.
Definition: fileName.H:59
static List< int > & procID(const label communicator)
The list of ranks within a given communicator.
Definition: UPstream.H:1130
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
Foam::fvBoundaryMesh.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: polyMesh.C:841
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
const polyBoundaryMesh & patches
Nothing to be read.
Automatically write from objectRegistry::writeObject()
const std::string patch
OpenFOAM patch number as a std::string.
virtual const fvGeometryScheme & geometry() const
Return reference to geometry calculation scheme.
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:678
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
bool checkDefinition(const bool report=false) const
Check zone definition. Return true if in error.
Definition: ZoneMesh.C:850
List< label > sortedToc(const UList< bool > &bools)
Return the (sorted) values corresponding to &#39;true&#39; entries.
Definition: BitOps.C:195
wordList names() const
A list of the zone names.
Definition: ZoneMesh.C:362
messageStream Info
Information stream (stdout output on master, null elsewhere)
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
#define doLocalCode(FieldType)
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return pointer to a new patch created on freestore from components.
Definition: polyPatchNew.C:28
static autoPtr< fvPatch > New(const polyPatch &, const fvBoundaryMesh &)
Return a pointer to a new patch created on freestore from polyPatch.
Definition: fvPatchNew.C:28
static const fileOperation & fileHandler()
Return the current file handler. Will create the default file handler if necessary.
List< label > labelList
A List of labels.
Definition: List.H:62
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Registry of regIOobjects.
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
Read and store dictionary entries for boundary patches The object is *never* registered to avoid regi...
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
List< bool > boolList
A List of bools.
Definition: List.H:60
static void setBasicGeometry(fvMesh &mesh)
Set the fvGeometryScheme to basic (to avoid parallel communication)
Definition: fvMeshTools.C:425
Do not request registration (bool: false)
static labelList removeEmptyPatches(fvMesh &, const bool validBoundary)
Remove zero sized patches. All but processor patches are.
Definition: fvMeshTools.C:363
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())