loadOrCreateMesh.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-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2021 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 "loadOrCreateMesh.H"
30 #include "processorPolyPatch.H"
32 #include "Time.H"
34 #include "IOobjectList.H"
35 #include "pointSet.H"
36 #include "faceSet.H"
37 #include "cellSet.H"
38 #include "basicFvGeometryScheme.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 // Read mesh if available. Otherwise create empty mesh with same non-proc
43 // patches as proc0 mesh. Requires:
44 // - all processors to have all patches (and in same order).
45 // - io.instance() set to facesInstance
47 (
48  const bool decompose,
49  const IOobject& io
50 )
51 {
52  // Region name
53  // ~~~~~~~~~~~
54 
55  fileName meshSubDir;
56 
57  if (io.name() == polyMesh::defaultRegion)
58  {
59  meshSubDir = polyMesh::meshSubDir;
60  }
61  else
62  {
63  meshSubDir = io.name()/polyMesh::meshSubDir;
64  }
65 
66 
67  // Patch types
68  // ~~~~~~~~~~~
69  // Read and scatter master patches (without reading master mesh!)
70 
71  PtrList<entry> patchEntries;
72  if (Pstream::master())
73  {
74  const bool oldParRun = Pstream::parRun(false);
75 
76  patchEntries = polyBoundaryMeshEntries
77  (
78  IOobject
79  (
80  "boundary",
81  io.instance(), //facesInstance,
82  meshSubDir,
83  io.db(),
84  IOobject::MUST_READ,
85  IOobject::NO_WRITE,
86  false
87  )
88  );
89 
90  Pstream::parRun(oldParRun);
91 
92  // Send patches
93  for (const int slave : Pstream::subProcs())
94  {
95  OPstream toSlave(Pstream::commsTypes::scheduled, slave);
96  toSlave << patchEntries;
97  }
98  }
99  else
100  {
101  // Receive patches
102  IPstream fromMaster
103  (
104  Pstream::commsTypes::scheduled,
105  Pstream::masterNo()
106  );
107  fromMaster >> patchEntries;
108  }
109 
110 
111 
112  // Dummy meshes
113  // ~~~~~~~~~~~~
114 
115  // Check who has a mesh
116  bool haveMesh;
117  if (decompose)
118  {
119  // Mesh needs to be present on the master only
120  haveMesh = Pstream::master();
121  }
122  else
123  {
124  const fileName facesFile
125  (
126  io.time().path()
127  /io.instance() //facesInstance
128  /meshSubDir
129  /"faces"
130  );
131 
132  // Check presence of the searched-for faces file
133  haveMesh = fileHandler().isFile(fileHandler().filePath(facesFile));
134  }
135 
136  if (!haveMesh)
137  {
138  const bool oldParRun = Pstream::parRun(false);
139 
140 
141  // Create dummy mesh. Only used on procs that don't have mesh.
142  IOobject noReadIO(io);
143  noReadIO.readOpt(IOobject::NO_READ);
144  noReadIO.writeOpt(IOobject::AUTO_WRITE);
145  fvMesh dummyMesh(noReadIO, Zero, false);
146 
147  // Add patches
148  List<polyPatch*> patches(patchEntries.size());
149  label nPatches = 0;
150 
151  forAll(patchEntries, patchi)
152  {
153  const entry& e = patchEntries[patchi];
154  const word type(e.dict().get<word>("type"));
155  const word& name = e.keyword();
156 
157  if
158  (
159  type != processorPolyPatch::typeName
160  && type != processorCyclicPolyPatch::typeName
161  )
162  {
163  dictionary patchDict(e.dict());
164  patchDict.set("nFaces", 0);
165  patchDict.set("startFace", 0);
166 
167  patches[patchi] = polyPatch::New
168  (
169  name,
170  patchDict,
171  nPatches++,
172  dummyMesh.boundaryMesh()
173  ).ptr();
174  }
175  }
176  patches.setSize(nPatches);
177  dummyMesh.addFvPatches(patches, false); // no parallel comms
178 
179 
180  // Add some dummy zones so upon reading it does not read them
181  // from the undecomposed case. Should be done as extra argument to
182  // regIOobject::readStream?
183  List<pointZone*> pz
184  (
185  1,
186  new pointZone
187  (
188  "dummyPointZone",
189  0,
190  dummyMesh.pointZones()
191  )
192  );
193  List<faceZone*> fz
194  (
195  1,
196  new faceZone
197  (
198  "dummyFaceZone",
199  0,
200  dummyMesh.faceZones()
201  )
202  );
203  List<cellZone*> cz
204  (
205  1,
206  new cellZone
207  (
208  "dummyCellZone",
209  0,
210  dummyMesh.cellZones()
211  )
212  );
213  dummyMesh.addZones(pz, fz, cz);
214  dummyMesh.pointZones().clear();
215  dummyMesh.faceZones().clear();
216  dummyMesh.cellZones().clear();
217  //Pout<< "Writing dummy mesh to " << dummyMesh.polyMesh::objectPath()
218  // << endl;
219  dummyMesh.write();
220 
221  Pstream::parRun(oldParRun); // Restore parallel state
222  }
223 
224 
225 
226  // Read mesh
227  // ~~~~~~~~~
228  // Now all processors have a (possibly zero size) mesh so read in
229  // parallel
230 
231  //Pout<< "Reading mesh from " << io.objectPath() << endl;
232  auto meshPtr = autoPtr<fvMesh>::New(io);
233  fvMesh& mesh = *meshPtr;
234 
235  // Make sure to use a non-parallel geometry calculation method
236  {
237  tmp<fvGeometryScheme> basicGeometry
238  (
240  (
241  mesh,
242  dictionary(),
243  basicFvGeometryScheme::typeName
244  )
245  );
246  mesh.geometry(basicGeometry);
247  }
248 
249 
250  // Sync patches
251  // ~~~~~~~~~~~~
252 
253  if (!Pstream::master() && haveMesh)
254  {
255  // Check master names against mine
256 
257  const polyBoundaryMesh& patches = mesh.boundaryMesh();
258 
259  forAll(patchEntries, patchi)
260  {
261  const entry& e = patchEntries[patchi];
262  const word type(e.dict().get<word>("type"));
263  const word& name = e.keyword();
264 
265  if
266  (
267  type == processorPolyPatch::typeName
268  || type == processorCyclicPolyPatch::typeName
269  )
270  {
271  break;
272  }
273 
274  if (patchi >= patches.size())
275  {
277  << "Non-processor patches not synchronised."
278  << endl
279  << "Processor " << Pstream::myProcNo()
280  << " has only " << patches.size()
281  << " patches, master has "
282  << patchi
283  << exit(FatalError);
284  }
285 
286  if
287  (
288  type != patches[patchi].type()
289  || name != patches[patchi].name()
290  )
291  {
293  << "Non-processor patches not synchronised."
294  << endl
295  << "Master patch " << patchi
296  << " name:" << type
297  << " type:" << type << endl
298  << "Processor " << Pstream::myProcNo()
299  << " patch " << patchi
300  << " has name:" << patches[patchi].name()
301  << " type:" << patches[patchi].type()
302  << exit(FatalError);
303  }
304  }
305  }
306 
307 
308  // Determine zones
309  // ~~~~~~~~~~~~~~~
310 
311  wordList pointZoneNames(mesh.pointZones().names());
312  Pstream::scatter(pointZoneNames);
313  wordList faceZoneNames(mesh.faceZones().names());
314  Pstream::scatter(faceZoneNames);
315  wordList cellZoneNames(mesh.cellZones().names());
316  Pstream::scatter(cellZoneNames);
317 
318  if (!haveMesh)
319  {
320  // Add the zones. Make sure to remove the old dummy ones first
321  mesh.pointZones().clear();
322  mesh.faceZones().clear();
323  mesh.cellZones().clear();
324 
325  List<pointZone*> pz(pointZoneNames.size());
326  forAll(pointZoneNames, i)
327  {
328  pz[i] = new pointZone
329  (
330  pointZoneNames[i],
331  i,
332  mesh.pointZones()
333  );
334  }
335  List<faceZone*> fz(faceZoneNames.size());
336  forAll(faceZoneNames, i)
337  {
338  fz[i] = new faceZone
339  (
340  faceZoneNames[i],
341  i,
342  mesh.faceZones()
343  );
344  }
345  List<cellZone*> cz(cellZoneNames.size());
346  forAll(cellZoneNames, i)
347  {
348  cz[i] = new cellZone
349  (
350  cellZoneNames[i],
351  i,
352  mesh.cellZones()
353  );
354  }
355  mesh.addZones(pz, fz, cz);
356  }
357 
358 
359  // Determine sets
360  // ~~~~~~~~~~~~~~
361 
362  wordList pointSetNames;
363  wordList faceSetNames;
364  wordList cellSetNames;
365  if (Pstream::master())
366  {
367  // Read sets
368  const bool oldParRun = Pstream::parRun(false);
369  IOobjectList objects(mesh, mesh.facesInstance(), "polyMesh/sets");
370  Pstream::parRun(oldParRun);
371 
372  pointSetNames = objects.sortedNames(pointSet::typeName);
373  faceSetNames = objects.sortedNames(faceSet::typeName);
374  cellSetNames = objects.sortedNames(cellSet::typeName);
375  }
376  Pstream::scatter(pointSetNames);
377  Pstream::scatter(faceSetNames);
378  Pstream::scatter(cellSetNames);
379 
380  if (!haveMesh)
381  {
382  for (const word& setName : pointSetNames)
383  {
384  pointSet(mesh, setName, 0).write();
385  }
386  for (const word& setName : faceSetNames)
387  {
388  faceSet(mesh, setName, 0).write();
389  }
390  for (const word& setName : cellSetNames)
391  {
392  cellSet(mesh, setName, 0).write();
393  }
394  }
395 
396 
397  // Force recreation of globalMeshData.
398  mesh.globalData();
399 
400 
401  // Do some checks.
402 
403  // Check if the boundary definition is unique
404  mesh.boundaryMesh().checkDefinition(true);
405  // Check if the boundary processor patches are correct
406  mesh.boundaryMesh().checkParallelSync(true);
407  // Check names of zones are equal
408  mesh.cellZones().checkDefinition(true);
409  mesh.cellZones().checkParallelSync(true);
410  mesh.faceZones().checkDefinition(true);
411  mesh.faceZones().checkParallelSync(true);
412  mesh.pointZones().checkDefinition(true);
413  mesh.pointZones().checkParallelSync(true);
414 
415  return meshPtr;
416 }
417 
418 
419 bool Foam::removeEmptyDir(const fileName& path)
420 {
421  // Return true if empty directory. Note bypass of fileHandler to be
422  // consistent with polyMesh.removeFiles for now.
423 
424  {
425  fileNameList files
426  (
428  (
429  path,
430  fileName::FILE,
431  false, // filterGz
432  false // followLink
433  )
434  );
435  if (files.size())
436  {
437  return false;
438  }
439  }
440  {
441  fileNameList dirs
442  (
444  (
445  path,
446  fileName::DIRECTORY,
447  false, // filterGz
448  false // followLink
449  )
450  );
451  if (dirs.size())
452  {
453  return false;
454  }
455  }
456  {
457  fileNameList links
458  (
460  (
461  path,
462  fileName::LINK,
463  false, // filterGz
464  false // followLink
465  )
466  );
467  if (links.size())
468  {
469  return false;
470  }
471  }
472  {
473  fileNameList other
474  (
476  (
477  path,
478  fileName::UNDEFINED,
479  false, // filterGz
480  false // followLink
481  )
482  );
483  if (other.size())
484  {
485  return false;
486  }
487  }
488 
489  // Avoid checking success of deletion since initial path might not
490  // exist (e.g. contain 'region0'). Will stop when trying to delete
491  // parent directory anyway since now not empty.
492  Foam::rm(path);
493  return true;
494 }
495 
496 
497 void Foam::removeEmptyDirs(const fileName& meshPath)
498 {
499  // Delete resulting directory if empty
500  fileName path(meshPath);
501  path.clean();
502 
503  // Do subdirectories
504  {
505  const fileNameList dirs
506  (
508  (
509  path,
510  fileName::DIRECTORY,
511  false, // filterGz
512  false // followLink
513  )
514  );
515  for (const auto& dir : dirs)
516  {
517  removeEmptyDirs(path/dir);
518  }
519  }
520 
522 }
523 
524 
525 // ************************************************************************* //
meshPtr
Foam::autoPtr< Foam::fvMesh > meshPtr(nullptr)
Foam::loadOrCreateMesh
autoPtr< fvMesh > loadOrCreateMesh(const bool decompose, const IOobject &io)
Load (if it exists) or create zero cell mesh given an IOobject:
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::rm
bool rm(const fileName &file)
Remove a file (or its gz equivalent), returning true if successful.
Definition: MSwindows.C:1004
Foam::fileHandler
const fileOperation & fileHandler()
Get current file handler.
Definition: fileOperation.C:1485
IOobjectList.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
polyBoundaryMeshEntries.H
processorCyclicPolyPatch.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:62
nPatches
const label nPatches
Definition: printMeshSummary.H:30
faceSet.H
Foam::FatalError
error FatalError
processorPolyPatch.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::fileNameList
List< fileName > fileNameList
A List of fileNames.
Definition: fileNameList.H:58
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Time.H
Foam::autoPtr< Foam::fvMesh >
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::removeEmptyDirs
void removeEmptyDirs(const fileName &path)
Remove empty directories from bottom up.
Foam::fileOperation::isFile
virtual bool isFile(const fileName &, const bool checkGzip=true, const bool followLink=true) const =0
Does the name exist as a FILE in the file system?
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
path
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
Foam::removeEmptyDir
bool removeEmptyDir(const fileName &path)
Remove empty directory. Return true if removed.
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
cellSet.H
basicFvGeometryScheme.H
Foam::readDir
fileNameList readDir(const fileName &directory, const fileName::Type type=fileName::FILE, const bool filtergz=true, const bool followLink=true)
Read a directory and return the entries as a fileName List.
Definition: MSwindows.C:707
loadOrCreateMesh.H
Load or create (0 size) a mesh. Used in distributing meshes to a larger number of processors.
pointSet.H