Go to the documentation of this file.
55 if (io.
name() == polyMesh::defaultRegion)
57 meshSubDir = polyMesh::meshSubDir;
61 meshSubDir = io.
name()/polyMesh::meshSubDir;
70 if (Pstream::master())
122 int slave=Pstream::firstSlave();
123 slave<=Pstream::lastSlave();
127 OPstream toSlave(Pstream::scheduled, slave);
128 toSlave << patchEntries;
134 IPstream fromMaster(Pstream::scheduled, Pstream::masterNo());
135 fromMaster >> patchEntries;
145 const bool haveMesh =
isFile
153 bool oldParRun = Pstream::parRun();
154 Pstream::parRun() =
false;
159 noReadIO.
readOpt() = IOobject::NO_READ;
174 forAll(patchEntries, patchI)
176 const entry&
e = patchEntries[patchI];
177 const word type(
e.dict().lookup(
"type"));
182 type != processorPolyPatch::typeName
183 &&
type != processorCyclicPolyPatch::typeName
187 patchDict.
set(
"nFaces", 0);
188 patchDict.
set(
"startFace", 0);
248 Pstream::parRun() = oldParRun;
267 if (!Pstream::master() && haveMesh)
273 forAll(patchEntries, patchI)
275 const entry&
e = patchEntries[patchI];
276 const word type(
e.dict().lookup(
"type"));
279 if (
type == processorPolyPatch::typeName)
287 <<
"Non-processor patches not synchronised."
289 <<
"Processor " << Pstream::myProcNo()
290 <<
" has only " <<
patches.size()
291 <<
" patches, master has "
303 <<
"Non-processor patches not synchronised."
305 <<
"Master patch " << patchI
308 <<
"Processor " << Pstream::myProcNo()
309 <<
" patch " << patchI
310 <<
" has name:" <<
patches[patchI].name()
311 <<
" type:" <<
patches[patchI].type()
322 Pstream::scatter(pointZoneNames);
324 Pstream::scatter(faceZoneNames);
326 Pstream::scatter(cellZoneNames);
331 mesh.pointZones().clear();
332 mesh.faceZones().clear();
333 mesh.cellZones().clear();
369 mesh.addZones(pz, fz, cz);
390 mesh.boundaryMesh().checkDefinition(
true);
392 mesh.boundaryMesh().checkParallelSync(
true);
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);
A keyword and a list of tokens is an 'entry'.
vectorField pointField
pointField is a vectorField.
void clear()
Clear the zones.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
A class for handling words, derived from string.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
A class for handling file names.
List< label > labelList
A List of labels.
Foam::polyBoundaryMeshEntries.
Load or create (0 size) a mesh. Used in distributing meshes to a larger number of processors.
#define forAll(list, i)
Loop across all elements in list.
Output inter-processor communications stream.
void addFvPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches. Constructor helper.
A subset of mesh points. The labels of points in the zone can be obtained from the addressing() list.
const cellZoneMesh & cellZones() const
Return cell zone mesh.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const fileName & instance() const
const Time & time() const
Return time.
const faceZoneMesh & faceZones() const
Return face zone mesh.
const objectRegistry & db() const
Return the local objectRegistry.
virtual bool write() const
Write mesh using IO settings from time.
autoPtr< fvMesh > loadOrCreateMesh(const IOobject &io)
Load (if it exists) or create zero cell mesh given an IOobject:
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
readOption readOpt() const
A subset of mesh faces organised as a primitive patch.
const pointZoneMesh & pointZones() const
Return point zone mesh.
const word & name() const
Return name.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Mesh data needed to do the Finite Volume discretisation.
bool isFile(const fileName &, const bool checkGzip=true)
Does the name exist as a FILE in the file system?
List< bool > boolList
Bool container classes.
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".
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
fileName path() const
Return path.
Xfer< T > xferCopy(const T &)
Construct by copying the contents of the arg.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
const dimensionedScalar e
Elementary charge.
label size() const
Return the number of elements in the PtrList.
Input inter-processor communications stream.
void size(const label)
Override size to be inconsistent with allocated storage.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
void addZones(const List< pointZone * > &pz, const List< faceZone * > &fz, const List< cellZone * > &cz)
Add mesh zones.
word name(const complex &)
Return a string representation of a complex.
void set(entry *)
Assign a new entry, overwrite any existing entry.