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 | Copyright (C) 2012-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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 "Time.H"
30 //#include "IOPtrList.H"
32 
33 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
34 
35 //namespace Foam
36 //{
37 // defineTemplateTypeNameAndDebug(IOPtrList<entry>, 0);
38 //}
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 // Read mesh if available. Otherwise create empty mesh with same non-proc
43 // patches as proc0 mesh. Requires all processors to have all patches
44 // (and in same order).
46 (
47  const IOobject& io
48 )
49 {
50  // Region name
51  // ~~~~~~~~~~~
52 
53  fileName meshSubDir;
54 
55  if (io.name() == polyMesh::defaultRegion)
56  {
57  meshSubDir = polyMesh::meshSubDir;
58  }
59  else
60  {
61  meshSubDir = io.name()/polyMesh::meshSubDir;
62  }
63 
64 
65  // Patch types
66  // ~~~~~~~~~~~
67  // Read and scatter master patches (without reading master mesh!)
68 
69  PtrList<entry> patchEntries;
70  if (Pstream::master())
71  {
73  //const word oldTypeName = IOPtrList<entry>::typeName;
74  //const_cast<word&>(IOPtrList<entry>::typeName) = word::null;
75  //IOPtrList<entry> dictList
76  //(
77  // IOobject
78  // (
79  // "boundary",
80  // io.time().findInstance
81  // (
82  // meshSubDir,
83  // "boundary",
84  // IOobject::MUST_READ
85  // ),
86  // meshSubDir,
87  // io.db(),
88  // IOobject::MUST_READ,
89  // IOobject::NO_WRITE,
90  // false
91  // )
92  //);
93  //const_cast<word&>(IOPtrList<entry>::typeName) = oldTypeName;
95  //const_cast<word&>(dictList.type()) = dictList.headerClassName();
96  //
97  //patchEntries.transfer(dictList);
98  const fileName facesInstance = io.time().findInstance
99  (
100  meshSubDir,
101  "faces",
102  IOobject::MUST_READ
103  );
104 
105  patchEntries = polyBoundaryMeshEntries
106  (
107  IOobject
108  (
109  "boundary",
110  facesInstance,
111  meshSubDir,
112  io.db(),
113  IOobject::MUST_READ,
114  IOobject::NO_WRITE,
115  false
116  )
117  );
118 
119  // Send patches
120  for
121  (
122  int slave=Pstream::firstSlave();
123  slave<=Pstream::lastSlave();
124  slave++
125  )
126  {
127  OPstream toSlave(Pstream::scheduled, slave);
128  toSlave << patchEntries;
129  }
130  }
131  else
132  {
133  // Receive patches
134  IPstream fromMaster(Pstream::scheduled, Pstream::masterNo());
135  fromMaster >> patchEntries;
136  }
137 
138 
139 
140  // Dummy meshes
141  // ~~~~~~~~~~~~
142 
143  // Check who has a mesh
144  //const bool haveMesh = isDir(io.time().path()/io.instance()/meshSubDir);
145  const bool haveMesh = isFile
146  (
147  io.time().path()/io.instance()/meshSubDir/"faces"
148  );
149 
150 
151  if (!haveMesh)
152  {
153  bool oldParRun = Pstream::parRun();
154  Pstream::parRun() = false;
155 
156 
157  // Create dummy mesh. Only used on procs that don't have mesh.
158  IOobject noReadIO(io);
159  noReadIO.readOpt() = IOobject::NO_READ;
160  fvMesh dummyMesh
161  (
162  noReadIO,
163  xferCopy(pointField()),
164  xferCopy(faceList()),
165  xferCopy(labelList()),
166  xferCopy(labelList()),
167  false
168  );
169 
170  // Add patches
171  List<polyPatch*> patches(patchEntries.size());
172  label nPatches = 0;
173 
174  forAll(patchEntries, patchI)
175  {
176  const entry& e = patchEntries[patchI];
177  const word type(e.dict().lookup("type"));
178  const word& name = e.keyword();
179 
180  if
181  (
182  type != processorPolyPatch::typeName
183  && type != processorCyclicPolyPatch::typeName
184  )
185  {
186  dictionary patchDict(e.dict());
187  patchDict.set("nFaces", 0);
188  patchDict.set("startFace", 0);
189 
190  patches[patchI] = polyPatch::New
191  (
192  name,
193  patchDict,
194  nPatches++,
195  dummyMesh.boundaryMesh()
196  ).ptr();
197  }
198  }
199  patches.setSize(nPatches);
200  dummyMesh.addFvPatches(patches, false); // no parallel comms
201 
202 
203  // Add some dummy zones so upon reading it does not read them
204  // from the undecomposed case. Should be done as extra argument to
205  // regIOobject::readStream?
207  (
208  1,
209  new pointZone
210  (
211  "dummyPointZone",
212  labelList(0),
213  0,
214  dummyMesh.pointZones()
215  )
216  );
217  List<faceZone*> fz
218  (
219  1,
220  new faceZone
221  (
222  "dummyFaceZone",
223  labelList(0),
224  boolList(0),
225  0,
226  dummyMesh.faceZones()
227  )
228  );
229  List<cellZone*> cz
230  (
231  1,
232  new cellZone
233  (
234  "dummyCellZone",
235  labelList(0),
236  0,
237  dummyMesh.cellZones()
238  )
239  );
240  dummyMesh.addZones(pz, fz, cz);
241  dummyMesh.pointZones().clear();
242  dummyMesh.faceZones().clear();
243  dummyMesh.cellZones().clear();
244  //Pout<< "Writing dummy mesh to " << dummyMesh.polyMesh::objectPath()
245  // << endl;
246  dummyMesh.write();
247 
248  Pstream::parRun() = oldParRun;
249  }
250 
251 
252 
253  // Read mesh
254  // ~~~~~~~~~
255  // Now all processors have a (possibly zero size) mesh so read in
256  // parallel
257 
258  //Pout<< "Reading mesh from " << io.objectPath() << endl;
259  autoPtr<fvMesh> meshPtr(new fvMesh(io));
260  fvMesh& mesh = meshPtr();
261 
262 
263 
264  // Sync patches
265  // ~~~~~~~~~~~~
266 
267  if (!Pstream::master() && haveMesh)
268  {
269  // Check master names against mine
270 
271  const polyBoundaryMesh& patches = mesh.boundaryMesh();
272 
273  forAll(patchEntries, patchI)
274  {
275  const entry& e = patchEntries[patchI];
276  const word type(e.dict().lookup("type"));
277  const word& name = e.keyword();
278 
279  if (type == processorPolyPatch::typeName)
280  {
281  break;
282  }
283 
284  if (patchI >= patches.size())
285  {
287  << "Non-processor patches not synchronised."
288  << endl
289  << "Processor " << Pstream::myProcNo()
290  << " has only " << patches.size()
291  << " patches, master has "
292  << patchI
293  << exit(FatalError);
294  }
295 
296  if
297  (
298  type != patches[patchI].type()
299  || name != patches[patchI].name()
300  )
301  {
303  << "Non-processor patches not synchronised."
304  << endl
305  << "Master patch " << patchI
306  << " name:" << type
307  << " type:" << type << endl
308  << "Processor " << Pstream::myProcNo()
309  << " patch " << patchI
310  << " has name:" << patches[patchI].name()
311  << " type:" << patches[patchI].type()
312  << exit(FatalError);
313  }
314  }
315  }
316 
317 
318  // Determine zones
319  // ~~~~~~~~~~~~~~~
320 
321  wordList pointZoneNames(mesh.pointZones().names());
322  Pstream::scatter(pointZoneNames);
323  wordList faceZoneNames(mesh.faceZones().names());
324  Pstream::scatter(faceZoneNames);
325  wordList cellZoneNames(mesh.cellZones().names());
326  Pstream::scatter(cellZoneNames);
327 
328  if (!haveMesh)
329  {
330  // Add the zones. Make sure to remove the old dummy ones first
331  mesh.pointZones().clear();
332  mesh.faceZones().clear();
333  mesh.cellZones().clear();
334 
335  List<pointZone*> pz(pointZoneNames.size());
336  forAll(pointZoneNames, i)
337  {
338  pz[i] = new pointZone
339  (
340  pointZoneNames[i],
341  labelList(0),
342  i,
343  mesh.pointZones()
344  );
345  }
346  List<faceZone*> fz(faceZoneNames.size());
347  forAll(faceZoneNames, i)
348  {
349  fz[i] = new faceZone
350  (
351  faceZoneNames[i],
352  labelList(0),
353  boolList(0),
354  i,
355  mesh.faceZones()
356  );
357  }
358  List<cellZone*> cz(cellZoneNames.size());
359  forAll(cellZoneNames, i)
360  {
361  cz[i] = new cellZone
362  (
363  cellZoneNames[i],
364  labelList(0),
365  i,
366  mesh.cellZones()
367  );
368  }
369  mesh.addZones(pz, fz, cz);
370  }
371 
372 
373 // if (!haveMesh)
374 // {
375 // // We created a dummy mesh file above. Delete it.
376 // const fileName meshFiles = io.time().path()/io.instance()/meshSubDir;
377 // //Pout<< "Removing dummy mesh " << meshFiles << endl;
378 // mesh.removeFiles();
379 // rmDir(meshFiles);
380 // }
381 //
382  // Force recreation of globalMeshData.
383 // mesh.clearOut();
384  mesh.globalData();
385 
386 
387  // Do some checks.
388 
389  // Check if the boundary definition is unique
390  mesh.boundaryMesh().checkDefinition(true);
391  // Check if the boundary processor patches are correct
392  mesh.boundaryMesh().checkParallelSync(true);
393  // Check names of zones are equal
394  mesh.cellZones().checkDefinition(true);
395  mesh.cellZones().checkParallelSync(true);
396  mesh.faceZones().checkDefinition(true);
397  mesh.faceZones().checkParallelSync(true);
398  mesh.pointZones().checkDefinition(true);
399  mesh.pointZones().checkParallelSync(true);
400 
401  return meshPtr;
402 }
403 
404 
405 // ************************************************************************* //
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:65
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Foam::ZoneMesh::clear
void clear()
Clear the zones.
Definition: ZoneMesh.C:408
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
Foam::polyBoundaryMeshEntries
Foam::polyBoundaryMeshEntries.
Definition: polyBoundaryMeshEntries.H:48
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
loadOrCreateMesh.H
Load or create (0 size) a mesh. Used in distributing meshes to a larger number of processors.
nPatches
label nPatches
Definition: readKivaGrid.H:402
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::OPstream
Output inter-processor communications stream.
Definition: OPstream.H:50
Foam::fvMesh::addFvPatches
void addFvPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:462
Foam::pointZone
A subset of mesh points. The labels of points in the zone can be obtained from the addressing() list.
Definition: pointZone.H:62
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:469
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
polyBoundaryMeshEntries.H
Foam::IOobject::instance
const fileName & instance() const
Definition: IOobject.H:350
Foam::IOobject::time
const Time & time() const
Return time.
Definition: IOobject.C:245
processorCyclicPolyPatch.H
Foam::cellZone
A subset of mesh cells.
Definition: cellZone.H:61
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
Foam::IOobject::db
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:239
Foam::fvMesh::write
virtual bool write() const
Write mesh using IO settings from time.
Definition: fvMesh.C:873
Foam::loadOrCreateMesh
autoPtr< fvMesh > loadOrCreateMesh(const IOobject &io)
Load (if it exists) or create zero cell mesh given an IOobject:
Definition: loadOrCreateMesh.C:46
Foam::label
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
Foam::IOobject::readOpt
readOption readOpt() const
Definition: IOobject.H:317
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
Foam::polyMesh::pointZones
const pointZoneMesh & pointZones() const
Return point zone mesh.
Definition: polyMesh.H:457
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobject.H:273
Foam::PtrList< entry >
meshPtr
static fvMesh * meshPtr
Definition: globalFoam.H:52
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
processorPolyPatch.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::isFile
bool isFile(const fileName &, const bool checkGzip=true)
Does the name exist as a FILE in the file system?
Definition: POSIX.C:622
Foam::boolList
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
Foam::Time::findInstance
word findInstance(const fileName &dir, const word &name=word::null, const IOobject::readOption rOpt=IOobject::MUST_READ, const word &stopInstance=word::null) const
Return the location of "dir" containing the file "name".
Definition: findInstance.C:38
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::autoPtr< Foam::fvMesh >
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::Time::path
fileName path() const
Return path.
Definition: Time.H:281
Foam::faceList
List< face > faceList
Definition: faceListFwd.H:43
Foam::xferCopy
Xfer< T > xferCopy(const T &)
Construct by copying the contents of the arg.
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:94
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
patches
patches[0]
Definition: createSingleCellMesh.H:36
Foam::IPstream
Input inter-processor communications stream.
Definition: IPstream.H:50
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::polyMesh::addZones
void addZones(const List< pointZone * > &pz, const List< faceZone * > &fz, const List< cellZone * > &cz)
Add mesh zones.
Definition: polyMesh.C:922
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::dictionary::set
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:856