35 template<
class Source>
36 void Foam::blockMesh::checkPatchLabels
39 const word& patchName,
46 face&
f = patchFaces[facei];
58 <<
"Block index out of range for patch face " <<
f <<
nl
59 <<
" Number of blocks = " <<
size()
60 <<
", index = " <<
f[0] <<
nl
61 <<
" on patch " << patchName <<
", face " << facei
64 else if (fi >=
operator[](bi).blockShape().
faces().
size())
67 <<
"Block face index out of range for patch face " <<
f
69 <<
" Number of block faces = "
71 <<
", index = " <<
f[1] <<
nl
72 <<
" on patch " << patchName <<
", face " << facei
87 <<
"Negative point label " <<
f[fp] <<
nl
88 <<
" on patch " << patchName <<
", face " << facei
94 <<
"Point label " <<
f[fp]
96 <<
" on patch " << patchName <<
", face " << facei
105 void Foam::blockMesh::readPatches
107 const dictionary& meshDescription,
115 dictionary varDict(meshDescription.subOrEmptyDict(
"namedVertices"));
116 varDict.merge(meshDescription.subOrEmptyDict(
"namedBlocks"));
119 ITstream& patchStream(meshDescription.lookup(
"patches"));
124 token firstToken(patchStream);
126 if (firstToken.isLabel())
137 patchStream.putBack(firstToken);
141 patchStream.readBegin(
"patches");
145 token lastToken(patchStream);
149 lastToken.isPunctuation()
154 if (tmpBlocksPatches.size() <=
nPatches)
156 tmpBlocksPatches.setSize(
nPatches + 1);
159 nbrPatchNames.setSize(
nPatches + 1);
162 patchStream.putBack(lastToken);
169 tmpBlocksPatches[
nPatches] = blockMeshTools::read<face>
183 <<
" at line " << patchStream.lineNumber()
207 <<
"Old-style cyclic definition."
208 <<
" Splitting patch "
210 << halfA <<
" and " << halfB <<
endl
211 <<
" Alternatively use new 'boundary' dictionary syntax."
215 if (tmpBlocksPatches.size() <=
nPatches)
217 tmpBlocksPatches.setSize(
nPatches + 1);
220 nbrPatchNames.setSize(
nPatches + 1);
233 if ((tmpBlocksPatches[
nPatches-1].size() % 2) != 0)
236 <<
"Size of cyclic faces is not a multiple of 2 :"
244 SubList<face>(unsplitFaces, sz)
248 SubList<face>(unsplitFaces, sz, sz)
254 patchStream >> lastToken;
256 patchStream.putBack(lastToken);
259 patchStream.readEnd(
"patches");
263 void Foam::blockMesh::readBoundary
265 const dictionary& meshDescription,
272 dictionary varDict(meshDescription.subOrEmptyDict(
"namedVertices"));
273 varDict.merge(meshDescription.subOrEmptyDict(
"namedBlocks"));
277 const PtrList<entry> patchesInfo
279 meshDescription.lookupOrDefault(
"boundary", PtrList<entry>())
283 tmpBlocksPatches.setSize(patchesInfo.size());
288 const entry& patchInfo = patchesInfo[
patchi];
290 if (!patchInfo.isDict())
293 <<
"Entry " << patchInfo <<
" in boundary section is not a"
294 <<
" valid dictionary."
304 tmpBlocksPatches[
patchi] = blockMeshTools::read<face>
322 void Foam::blockMesh::createCellShapes
327 const blockMesh& blocks = *
this;
329 tmpBlockCells.setSize(blocks.size());
332 tmpBlockCells[blocki] = blocks[blocki].blockShape();
339 void Foam::blockMesh::defaultPatchError
341 const word& defaultPatchName,
342 const dictionary& meshDescription
350 <<
"The 'defaultPatch' type must be specified"
351 " for the '" << defaultPatchName <<
"' patch, e.g. for snappyHexMesh"
355 " name default; // optional\n"
358 <<
"or for 2D meshes" <<
nl <<
nl
361 " name frontAndBack; // optional\n"
371 const IOdictionary& meshDescription,
376 checkBlockFaceOrientation = meshDescription.lookupOrDefault<Switch>
378 "checkBlockFaceOrientation",
384 word defaultPatchName =
"defaultFaces";
385 word defaultPatchType = emptyPolyPatch::typeName;
386 bool defaultPatchTypeSet =
false;
391 if (
const dictionary* dictPtr = meshDescription.subDictPtr(
"defaultPatch"))
393 dictPtr->readIfPresent(
"name", defaultPatchName);
394 defaultPatchType = dictPtr->lookup<word>(
"type");
395 defaultPatchTypeSet =
true;
399 if (!meshDescription.readIfPresent(
"convertToMeters", scaleFactor_))
401 meshDescription.readIfPresent(
"scale", scaleFactor_);
405 if (meshDescription.found(
"edges"))
409 Info<<
"Creating block edges" <<
endl;
414 meshDescription.lookup(
"edges"),
415 blockEdge::iNew(meshDescription, geometry_, vertices_)
418 edges_.transfer(edges);
420 else if (verboseOutput)
422 Info<<
"No non-linear block edges defined" <<
endl;
427 if (meshDescription.found(
"faces"))
431 Info<<
"Creating block faces" <<
endl;
436 meshDescription.lookup(
"faces"),
437 blockFace::iNew(meshDescription, geometry_)
440 faces_.transfer(faces);
442 else if (verboseOutput)
444 Info<<
"No non-planar block faces defined" <<
endl;
451 Info<<
"Creating topology blocks" <<
endl;
456 meshDescription.lookup(
"blocks"),
457 block::iNew(meshDescription, vertices_, edges_, faces_)
465 polyMesh* blockMeshPtr =
nullptr;
471 Info<<
"Creating topology patches" <<
endl;
474 if (meshDescription.found(
"patches"))
477 <<
"Reading patches section" <<
nl
478 <<
" The 'patches' entry is deprecated and has been superseded"
479 <<
" by the more consistent and general 'boundary' entry."
496 Info<<
nl <<
"Creating block mesh topology" <<
endl;
499 createCellShapes(tmpBlockCells);
501 Info<<
nl <<
"Reading physicalType from existing boundary file" <<
endl;
508 meshDescription.time(),
538 <<
"' (read from boundary file)"
545 dict.set(
"neighbourPatch", nbrPatchNames[
patchi]);
549 blockMeshPtr =
new polyMesh
555 meshDescription.time(),
572 && blockMeshPtr->boundaryMesh().findIndex(defaultPatchName) != -1
575 defaultPatchError(defaultPatchName, meshDescription);
592 Info<<
nl <<
"Creating block mesh topology" <<
endl;
595 createCellShapes(tmpBlockCells);
597 blockMeshPtr =
new polyMesh
603 meshDescription.time(),
620 && blockMeshPtr->boundaryMesh().findIndex(defaultPatchName) != -1
623 defaultPatchError(defaultPatchName, meshDescription);
627 check(*blockMeshPtr, meshDescription);
#define forAll(list, i)
Loop across all elements in list.
void size(const label)
Override size to be inconsistent with allocated storage.
void setSize(const label)
Reset size of List.
label size() const
Return the number of elements in the UPtrList.
const T & operator[](const label) const
Return element const reference.
const blockFaceList & faces() const
Return the curved faces.
const pointField & points() const
The points for the entire mesh.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
bool add(entry *, bool mergeEntry=false)
Add a new entry.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Mesh consisting of general polyhedral cells.
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
static const word null
An empty word.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
errorManipArg< error, int > exit(error &err, const int errNo=1)
List< cellShape > cellShapeList
List of cellShapes and PtrList of List of cellShape.
List< word > wordList
A List of words.
PtrList< blockFace > blockFaceList
A PtrList of blockFaces.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const word & regionName(const solver ®ion)
vectorField pointField
pointField is a vectorField.
PtrList< block > blockList
A PtrList of blocks.
PtrList< blockEdge > blockEdgeList
A PtrList of blockEdges.
List< faceList > faceListList
void preservePatchTypes(const objectRegistry &obr, const word &meshInstance, const fileName &meshDir, const wordList &patchNames, PtrList< dictionary > &patchDicts, const word &defaultFacesName, word &defaultFacesType)
Preserve patch types.
wordList patchTypes(nPatches)
wordList patchNames(nPatches)
PtrList< dictionary > patchDicts