Go to the documentation of this file.
42 bool topologyOK =
true;
49 token firstToken(patchStream);
62 patchStream.
putBack(firstToken);
70 token lastToken(patchStream);
75 && lastToken.
pToken() == token::END_LIST
94 patchStream >> tmpBlocksPatches[
nPatches];
105 <<
". Exiting !" <<
nl
110 topologyOK = topologyOK && patchLabelsOK
128 <<
"Old-style cyclic definition."
129 <<
" Splitting patch "
131 << halfA <<
" and " << halfB <<
endl
132 <<
" Alternatively use new 'boundary' dictionary syntax."
136 if (tmpBlocksPatches.size() <=
nPatches)
138 tmpBlocksPatches.setSize(
nPatches + 1);
153 if ((tmpBlocksPatches[
nPatches-1].size() % 2) != 0)
156 <<
"Size of cyclic faces is not a multiple of 2 :"
174 patchStream >> lastToken;
176 patchStream.
putBack(lastToken);
179 patchStream.
readEnd(
"patches");
193 bool topologyOK =
true;
198 meshDescription.
lookup(
"boundary")
205 forAll(tmpBlocksPatches, patchI)
207 const entry& patchInfo = patchesInfo[patchI];
212 <<
"Entry " << patchInfo <<
" in boundary section is not a"
220 patchDicts[patchI].lookup(
"faces") >> tmpBlocksPatches[patchI];
222 topologyOK = topologyOK && patchLabelsOK
226 tmpBlocksPatches[patchI]
244 tmpBlockCells[blockI] =
cellShape(blocks[blockI].blockShape());
246 if (tmpBlockCells[blockI].
mag(blockPointField_) < 0.0)
249 <<
"negative volume block : " << blockI
250 <<
", probably defined inside-out" <<
endl;
264 bool topologyOK =
true;
268 word defaultPatchName =
"defaultFaces";
269 word defaultPatchType = emptyPolyPatch::typeName;
276 dictPtr->readIfPresent(
"name", defaultPatchName);
277 dictPtr->readIfPresent(
"type", defaultPatchType);
281 if (!meshDescription.
readIfPresent(
"convertToMeters", scaleFactor_))
290 if (meshDescription.
found(
"edges"))
294 Info<<
"Creating curved edges" <<
endl;
302 token firstToken(is);
307 edges_.setSize(nEdges);
324 && lastToken.
pToken() == token::END_LIST
328 if (edges_.size() <= nEdges)
330 edges_.setSize(nEdges + 1);
350 else if (verboseOutput)
352 Info<<
"No non-linear edges defined" <<
endl;
361 Info<<
"Creating topology blocks" <<
endl;
370 token firstToken(is);
392 && lastToken.
pToken() == token::END_LIST
396 if (blocks.
size() <= nBlocks)
414 topologyOK = topologyOK && blockLabelsOK
418 blocks[nBlocks].blockShape()
439 Info<<
"Creating topology patches" <<
endl;
442 if (meshDescription.
found(
"patches"))
444 Info<<
nl <<
"Reading patches section" <<
endl;
451 topologyOK = topologyOK && readPatches
463 <<
"Cannot create mesh due to errors in topology, exiting !"
467 Info<<
nl <<
"Creating block mesh topology" <<
endl;
470 createCellShapes(tmpBlockCells);
473 Info<<
nl <<
"Reading physicalType from existing boundary file" <<
endl;
480 meshDescription.
time(),
482 polyMesh::meshSubDir,
501 if (!
dict.found(
"type"))
512 <<
"' with '" <<
word(
dict.lookup(
"type"))
513 <<
"' (read from boundary file)"
518 if (nbrPatchNames[patchI] != word::null)
520 dict.set(
"neighbourPatch", nbrPatchNames[patchI]);
531 meshDescription.
time(),
545 else if (meshDescription.
found(
"boundary"))
551 topologyOK = topologyOK && readBoundary
562 <<
"Cannot create mesh due to errors in topology, exiting !"
567 Info<<
nl <<
"Creating block mesh topology" <<
endl;
570 createCellShapes(tmpBlockCells);
580 meshDescription.
time(),
595 checkBlockMesh(*blockMeshPtr);
A keyword and a list of tokens is an 'entry'.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Istream & readEnd(const char *funcName)
Creates a single block of cells from point coordinates, numbers of cells in each direction and an exp...
bool readPatches(const dictionary &meshDescription, faceListList &tmpBlocksPatches, wordList &patchNames, wordList &patchTypes, wordList &nbrPatchNames)
preservePatchTypes(runTime, runTime.constant(), polyMesh::meshSubDir, patchNames, patchDicts, defaultFacesName, defaultFacesType)
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)
const keyType & keyword() const
Return keyword.
#define forAll(list, i)
Loop across all elements in list.
bool readBoundary(const dictionary &meshDescription, wordList &patchNames, faceListList &tmpBlocksPatches, PtrList< dictionary > &patchDicts)
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
A List obtained as a section of another List.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Ostream & endl(Ostream &os)
Add newline and flush stream.
A token holds items read from Istream.
dimensioned< scalar > mag(const dimensioned< Type > &)
const Time & time() const
Return time.
Mesh consisting of general polyhedral cells.
PtrList< dictionary > patchDicts
wordList patchTypes(nPatches)
void createCellShapes(cellShapeList &tmpBlockCells)
bool set(const label) const
Is element set.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
virtual bool isDict() const
Return true if this entry is a dictionary.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
wordList patchNames(nPatches)
punctuationToken pToken() const
A list of keyword definitions, which are a keyword followed by any number of values (e....
An analytical geometric cellShape.
virtual const dictionary & dict() const =0
Return dictionary if this entry is a dictionary.
errorManipArg< error, int > exit(error &err, const int errNo=1)
void setSize(const label)
Reset size of List.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Xfer< T > xferCopy(const T &)
Construct by copying the contents of the arg.
void putBack(const token &)
Put back token.
const word & constant() const
Return constant name.
label size() const
Return the number of elements in the PtrList.
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
label lineNumber() const
Return current stream line number.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
void size(const label)
Override size to be inconsistent with allocated storage.
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
Istream & readBegin(const char *funcName)
bool isPunctuation() const
#define WarningInFunction
Report a warning using Foam::Warning.
polyMesh * createTopology(const IOdictionary &, const word ®ionName)
const dictionary * subDictPtr(const word &) const
Find and return a sub-dictionary pointer if present.
A multi-block mesh generator.