loadOrCreateMesh.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) 2012-2021 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 "loadOrCreateMesh.H"
27 #include "processorPolyPatch.H"
29 #include "IOPtrList.H"
30 #include "OSspecific.H"
31 
32 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTemplateTypeNameAndDebug(IOPtrList<entry>, 0);
37 }
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 // Read mesh if available. Otherwise create empty mesh with same non-proc
42 // patches as proc0 mesh. Requires all processors to have all patches
43 // (and in same order).
45 (
46  const IOobject& io
47 )
48 {
49  fileName meshSubDir;
50 
51  if (io.name() == polyMesh::defaultRegion)
52  {
53  meshSubDir = polyMesh::meshSubDir;
54  }
55  else
56  {
57  meshSubDir = io.name()/polyMesh::meshSubDir;
58  }
59 
60 
61  // Scatter master patches
62  PtrList<entry> patchEntries;
63  if (Pstream::master())
64  {
65  // Read PtrList of dictionary as dictionary.
66  const word oldTypeName = IOPtrList<entry>::typeName;
67  const_cast<word&>(IOPtrList<entry>::typeName) = word::null;
68  IOPtrList<entry> dictList
69  (
70  IOobject
71  (
72  "boundary",
73  io.time().findInstance
74  (
75  meshSubDir,
76  "boundary",
77  IOobject::MUST_READ
78  ),
79  meshSubDir,
80  io.db(),
81  IOobject::MUST_READ,
82  IOobject::NO_WRITE,
83  false
84  )
85  );
86  const_cast<word&>(IOPtrList<entry>::typeName) = oldTypeName;
87  // Fake type back to what was in field
88  const_cast<word&>(dictList.type()) = dictList.headerClassName();
89 
90  patchEntries.transfer(dictList);
91 
92  // Send patches
93  for
94  (
95  int slave=Pstream::firstSlave();
96  slave<=Pstream::lastSlave();
97  slave++
98  )
99  {
100  OPstream toSlave(Pstream::commsTypes::scheduled, slave);
101  toSlave << patchEntries;
102  }
103  }
104  else
105  {
106  // Receive patches
107  IPstream fromMaster
108  (
109  Pstream::commsTypes::scheduled,
110  Pstream::masterNo()
111  );
112  fromMaster >> patchEntries;
113  }
114 
115 
116 
117  // Check who has a mesh
118  const bool haveMesh = fileHandler().isDir
119  (
120  fileHandler().filePath(io.time().path()/io.instance()/meshSubDir)
121  );
122 
123  if (!haveMesh)
124  {
125  bool oldParRun = Pstream::parRun();
126  Pstream::parRun() = false;
127 
128 
129  // Create dummy mesh. Only used on procs that don't have mesh.
130  IOobject noReadIO(io);
131  noReadIO.readOpt() = IOobject::NO_READ;
132  fvMesh dummyMesh
133  (
134  noReadIO,
135  pointField(),
136  faceList(),
137  labelList(),
138  labelList(),
139  false
140  );
141 
142  // Add patches
143  List<polyPatch*> patches(patchEntries.size());
144  label nPatches = 0;
145 
146  forAll(patchEntries, patchi)
147  {
148  const entry& e = patchEntries[patchi];
149  const word type(e.dict().lookup("type"));
150  const word& name = e.keyword();
151 
152  if
153  (
154  type != processorPolyPatch::typeName
155  && type != processorCyclicPolyPatch::typeName
156  )
157  {
158  dictionary patchDict(e.dict());
159  patchDict.set("nFaces", 0);
160  patchDict.set("startFace", 0);
161 
163  (
164  name,
165  patchDict,
166  nPatches++,
167  dummyMesh.boundaryMesh()
168  ).ptr();
169  }
170  }
171  patches.setSize(nPatches);
172  dummyMesh.addFvPatches(patches, false); // no parallel comms
173 
174  // Add some dummy zones so upon reading it does not read them
175  // from the undecomposed case. Should be done as extra argument to
176  // regIOobject::readStream?
177  List<pointZone*> pz
178  (
179  1,
180  new pointZone
181  (
182  "dummyPointZone",
183  labelList(0),
184  0,
185  dummyMesh.pointZones()
186  )
187  );
188  List<faceZone*> fz
189  (
190  1,
191  new faceZone
192  (
193  "dummyFaceZone",
194  labelList(0),
195  boolList(0),
196  0,
197  dummyMesh.faceZones()
198  )
199  );
200  List<cellZone*> cz
201  (
202  1,
203  new cellZone
204  (
205  "dummyCellZone",
206  labelList(0),
207  0,
208  dummyMesh.cellZones()
209  )
210  );
211  dummyMesh.addZones(pz, fz, cz);
212  dummyMesh.write();
213 
214  Pstream::parRun() = oldParRun;
215  }
216 
217  autoPtr<fvMesh> meshPtr(new fvMesh(io, false));
218  fvMesh& mesh = meshPtr();
219 
220 
221  // Sync patches
222  // ~~~~~~~~~~~~
223 
224  if (!Pstream::master() && haveMesh)
225  {
226  // Check master names against mine
227 
228  const polyBoundaryMesh& patches = mesh.boundaryMesh();
229 
230  forAll(patchEntries, patchi)
231  {
232  const entry& e = patchEntries[patchi];
233  const word type(e.dict().lookup("type"));
234  const word& name = e.keyword();
235 
236  if (type == processorPolyPatch::typeName)
237  {
238  break;
239  }
240 
241  if (patchi >= patches.size())
242  {
244  << "Non-processor patches not synchronised."
245  << endl
246  << "Processor " << Pstream::myProcNo()
247  << " has only " << patches.size()
248  << " patches, master has "
249  << patchi
250  << exit(FatalError);
251  }
252 
253  if
254  (
255  type != patches[patchi].type()
256  || name != patches[patchi].name()
257  )
258  {
260  << "Non-processor patches not synchronised."
261  << endl
262  << "Master patch " << patchi
263  << " name:" << type
264  << " type:" << type << endl
265  << "Processor " << Pstream::myProcNo()
266  << " patch " << patchi
267  << " has name:" << patches[patchi].name()
268  << " type:" << patches[patchi].type()
269  << exit(FatalError);
270  }
271  }
272  }
273 
274 
275  // Determine zones
276  // ~~~~~~~~~~~~~~~
277 
278  wordList pointZoneNames(mesh.pointZones().names());
279  Pstream::scatter(pointZoneNames);
280  wordList faceZoneNames(mesh.faceZones().names());
281  Pstream::scatter(faceZoneNames);
282  wordList cellZoneNames(mesh.cellZones().names());
283  Pstream::scatter(cellZoneNames);
284 
285  if (!haveMesh)
286  {
287  // Add the zones. Make sure to remove the old dummy ones first
288  mesh.pointZones().clear();
289  mesh.faceZones().clear();
290  mesh.cellZones().clear();
291 
292  List<pointZone*> pz(pointZoneNames.size());
293  forAll(pointZoneNames, i)
294  {
295  pz[i] = new pointZone
296  (
297  pointZoneNames[i],
298  labelList(0),
299  i,
300  mesh.pointZones()
301  );
302  }
303  List<faceZone*> fz(faceZoneNames.size());
304  forAll(faceZoneNames, i)
305  {
306  fz[i] = new faceZone
307  (
308  faceZoneNames[i],
309  labelList(0),
310  boolList(0),
311  i,
312  mesh.faceZones()
313  );
314  }
315  List<cellZone*> cz(cellZoneNames.size());
316  forAll(cellZoneNames, i)
317  {
318  cz[i] = new cellZone
319  (
320  cellZoneNames[i],
321  labelList(0),
322  i,
323  mesh.cellZones()
324  );
325  }
326  mesh.addZones(pz, fz, cz);
327  }
328 
329 
330  if (!haveMesh)
331  {
332  // We created a dummy mesh file above. Delete it.
333  const fileName meshFiles = io.time().path()/io.instance()/meshSubDir;
334  // Pout<< "Removing dummy mesh " << meshFiles << endl;
335  mesh.removeFiles();
336  rmDir(meshFiles);
337  }
338 
339  // Force recreation of globalMeshData.
340  mesh.clearOut();
341  mesh.globalData();
342 
343 
344  // Do some checks.
345 
346  // Check if the boundary definition is unique
347  mesh.boundaryMesh().checkDefinition(true);
348  // Check if the boundary processor patches are correct
349  mesh.boundaryMesh().checkParallelSync(true);
350  // Check names of zones are equal
351  mesh.cellZones().checkDefinition(true);
352  mesh.cellZones().checkParallelSync(true);
353  mesh.faceZones().checkDefinition(true);
354  mesh.faceZones().checkParallelSync(true);
355  mesh.pointZones().checkDefinition(true);
356  mesh.pointZones().checkParallelSync(true);
357 
358  return meshPtr;
359 }
360 
361 
362 // ************************************************************************* //
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual bool isDir(const fileName &, const bool followLink=true) const =0
Does the name exist as a directory in the file system?
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
label patchi
static fvMesh * meshPtr
Definition: globalFoam.H:52
const fvPatchList & patches
Load or create (0 size) a mesh. Used in distributing meshes to a larger number of processors.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const dimensionedScalar e
Elementary charge.
Namespace for OpenFOAM.
const fileOperation & fileHandler()
Get current file handler.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
List< word > wordList
A List of words.
Definition: fileName.H:54
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
bool rmDir(const fileName &)
Remove a directory and its contents.
Definition: POSIX.C:1047
error FatalError
List< face > faceList
Definition: faceListFwd.H:43
defineTemplateTypeNameAndDebug(IOPtrList< ensightPart >, 0)
autoPtr< fvMesh > loadOrCreateMesh(const IOobject &io)
Load (if it exists) or create zero cell mesh given an IOobject:
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
label nPatches
Definition: readKivaGrid.H:396