blockMeshTopology.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2019-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "blockMesh.H"
30#include "blockMeshTools.H"
31#include "Time.H"
32#include "preservePatchTypes.H"
33#include "emptyPolyPatch.H"
34#include "cyclicPolyPatch.H"
35
36// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
37
38namespace Foam
39{
40
41// Replace (<block> <face>) face description
42// with the corresponding block face
43// - otherwise check point labels are in range
44
45template<class Source>
47(
48 const Source& source,
49 const word& patchName,
50 const PtrList<block>& blocks,
51 const label nPoints,
52 faceList& patchFaces
53)
54{
55 const label nBlocks = blocks.size();
56
57 forAll(patchFaces, facei)
58 {
59 face& f = patchFaces[facei];
60
61 // Replace (<block> <face>) face description
62 // with the corresponding block face
63 if (f.size() == 2)
64 {
65 const label bi = f[0];
66 const label fi = f[1];
67
68 if (bi < 0 || bi >= nBlocks)
69 {
71 << "Block index out of range."
72 << " Patch face (" << bi << ' ' << fi << ")\n"
73 << " Number of blocks = " << nBlocks
74 << ", block index = " << bi << nl
75 << " on patch " << patchName << ", face " << facei
77 }
78 else if (fi >= blocks[bi].blockShape().nFaces())
79 {
81 << "Block face index out of range."
82 << " Patch face (" << bi << ' ' << fi << ")\n"
83 << " Number of block faces = "
84 << blocks[bi].blockShape().nFaces()
85 << ", face index = " << fi << nl
86 << " on patch " << patchName << ", face " << facei
88 }
89 else
90 {
91 f = blocks[bi].blockShape().face(fi);
92 }
93 }
94 else
95 {
96 for (const label pointi : f)
97 {
98 if (pointi < 0 || pointi >= nPoints)
99 {
101 << "Point label " << pointi
102 << " out of range 0.." << (nPoints - 1) << nl
103 << " on patch " << patchName << ", face " << facei
104 << exit(FatalIOError);
105 }
106 }
107 }
108 }
109}
110
111} // End namespace Foam
112
113
114// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
115
116void Foam::blockMesh::readPatches
117(
118 const dictionary& meshDescription,
119 faceListList& tmpBlocksPatches,
120 wordList& patchNames,
121 wordList& patchTypes,
122 wordList& nbrPatchNames
123)
124{
125 // Collect all variables
126 dictionary varDict(meshDescription.subOrEmptyDict("namedVertices"));
127 varDict.merge(meshDescription.subOrEmptyDict("namedBlocks"));
128
129
130 ITstream& patchStream = meshDescription.lookup("patches");
131
132 // Read number of patches in mesh
133 label nPatches = 0;
134
135 if (patchStream.peek().isLabel())
136 {
137 patchStream >> nPatches;
138
139 tmpBlocksPatches.setSize(nPatches);
142 nbrPatchNames.setSize(nPatches);
143 }
144
145 // Read beginning of blocks
146 patchStream.readBegin("patches");
147
148 nPatches = 0;
149
150 token lastToken(patchStream);
151 while (!lastToken.isPunctuation(token::END_LIST))
152 {
153 if (tmpBlocksPatches.size() <= nPatches)
154 {
155 tmpBlocksPatches.setSize(nPatches + 1);
158 nbrPatchNames.setSize(nPatches + 1);
159 }
160
161 patchStream.putBack(lastToken);
162
163 patchStream
166
167 // Read patch faces
168 tmpBlocksPatches[nPatches] = blockMeshTools::read<face>
169 (
170 patchStream,
171 varDict
172 );
173
174
175 // Check for multiple patches
176 for (label i = 0; i < nPatches; i++)
177 {
178 if (patchNames[nPatches] == patchNames[i])
179 {
180 FatalIOErrorInFunction(patchStream)
181 << "Duplicate patch " << patchNames[nPatches]
182 << " at line " << patchStream.lineNumber()
183 << exit(FatalIOError);
184 }
185 }
186
188 (
189 patchStream,
191 *this,
192 vertices_.size(),
193 tmpBlocksPatches[nPatches]
194 );
195
196 ++nPatches;
197
198
199 // Split old style cyclics
200
202 {
203 word halfA = patchNames[nPatches-1] + "_half0";
204 word halfB = patchNames[nPatches-1] + "_half1";
205
206 FatalIOErrorInFunction(patchStream)
207 << "Old-style cyclic definition."
208 << " Splitting patch "
209 << patchNames[nPatches-1] << " into two halves "
210 << halfA << " and " << halfB << endl
211 << " Alternatively use new 'boundary' dictionary syntax."
212 << exit(FatalIOError);
213
214 // Add extra patch
215 if (tmpBlocksPatches.size() <= nPatches)
216 {
217 tmpBlocksPatches.setSize(nPatches + 1);
220 nbrPatchNames.setSize(nPatches + 1);
221 }
222
223 // Update halfA info
224 patchNames[nPatches-1] = halfA;
225 nbrPatchNames[nPatches-1] = halfB;
226
227 // Update halfB info
229 patchNames[nPatches] = halfB;
230 nbrPatchNames[nPatches] = halfA;
231
232 // Split faces
233 if ((tmpBlocksPatches[nPatches-1].size() % 2) != 0)
234 {
235 FatalIOErrorInFunction(patchStream)
236 << "Size of cyclic faces is not a multiple of 2 :"
237 << tmpBlocksPatches[nPatches-1]
238 << exit(FatalIOError);
239 }
240 label sz = tmpBlocksPatches[nPatches-1].size()/2;
241 faceList unsplitFaces(tmpBlocksPatches[nPatches-1], true);
242 tmpBlocksPatches[nPatches-1] = faceList
243 (
244 SubList<face>(unsplitFaces, sz)
245 );
246 tmpBlocksPatches[nPatches] = faceList
247 (
248 SubList<face>(unsplitFaces, sz, sz)
249 );
250
251 nPatches++;
252 }
253
254 patchStream >> lastToken;
255 }
256 patchStream.putBack(lastToken);
257
258 // Read end of blocks
259 patchStream.readEnd("patches");
260}
261
262
263void Foam::blockMesh::readBoundary
264(
265 const dictionary& meshDescription,
267 faceListList& tmpBlocksPatches,
268 PtrList<dictionary>& patchDicts
269)
270{
271 // Collect all variables
272 dictionary varDict(meshDescription.subOrEmptyDict("namedVertices"));
273 varDict.merge(meshDescription.subOrEmptyDict("namedBlocks"));
274
275
276 // Read like boundary file
277 const PtrList<entry> patchesInfo
278 (
279 meshDescription.lookup("boundary")
280 );
281
282 patchNames.setSize(patchesInfo.size());
283 tmpBlocksPatches.setSize(patchesInfo.size());
284 patchDicts.setSize(patchesInfo.size());
285
286 forAll(tmpBlocksPatches, patchi)
287 {
288 const entry& patchInfo = patchesInfo[patchi];
289
290 if (!patchInfo.isDict())
291 {
292 FatalIOErrorInFunction(meshDescription)
293 << "Entry " << patchInfo << " in boundary section is not a"
294 << " valid dictionary."
295 << exit(FatalIOError);
296 }
297
298 patchNames[patchi] = patchInfo.keyword();
299
300 // Construct patch dictionary
301 patchDicts.set(patchi, new dictionary(patchInfo.dict()));
302
303 // Read block faces
304 tmpBlocksPatches[patchi] = blockMeshTools::read<face>
305 (
306 patchDicts[patchi].lookup("faces"),
307 varDict
308 );
309
310
312 (
313 patchInfo.dict(),
314 patchNames[patchi],
315 *this,
316 vertices_.size(),
317 tmpBlocksPatches[patchi]
318 );
319 }
320}
321
322
323Foam::cellShapeList Foam::blockMesh::getBlockShapes() const
324{
325 const blockMesh& blocks = *this;
326
327 cellShapeList shapes(blocks.size());
328
329 forAll(blocks, blocki)
330 {
331 shapes[blocki] = blocks[blocki].blockShape();
332 }
333
334 return shapes;
335}
336
337
338// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
339
341Foam::blockMesh::createTopology
342(
343 const IOdictionary& meshDescription,
344 const word& regionName
345)
346{
347 word defaultPatchName = "defaultFaces";
348 word defaultPatchType = emptyPolyPatch::typeName;
349
350 // Read the names/types for the unassigned patch faces
351 // this is a bit heavy handed (and ugly), but there is currently
352 // no easy way to rename polyMesh patches subsequently
353 if (const dictionary* dictptr = meshDescription.findDict("defaultPatch"))
354 {
355 dictptr->readIfPresent("name", defaultPatchName);
356 dictptr->readIfPresent("type", defaultPatchType);
357 }
358
359
360 // Scaling, transformations
361 readPointTransforms(meshDescription);
362
363 // Read the block edges
364 if (meshDescription.found("edges"))
365 {
366 if (verbose_)
367 {
368 Info<< "Creating block edges" << endl;
369 }
370
371 blockEdgeList edges
372 (
373 meshDescription.lookup("edges"),
374 blockEdge::iNew(meshDescription, geometry_, vertices_)
375 );
376
377 edges_.transfer(edges);
378 }
379 else
380 {
381 edges_.clear();
382 if (verbose_)
383 {
384 Info<< "No non-linear block edges defined" << endl;
385 }
386 }
387
388
389 // Read the block faces
390 if (meshDescription.found("faces"))
391 {
392 if (verbose_)
393 {
394 Info<< "Creating block faces" << endl;
395 }
396
397 blockFaceList faces
398 (
399 meshDescription.lookup("faces"),
400 blockFace::iNew(meshDescription, geometry_)
401 );
402
403 faces_.transfer(faces);
404 }
405 else
406 {
407 faces_.clear();
408 if (verbose_)
409 {
410 Info<< "No non-planar block faces defined" << endl;
411 }
412 }
413
414
415 // Create the blocks
416 if (verbose_)
417 {
418 Info<< "Creating topology blocks" << endl;
419 }
420 {
421 blockList blocks
422 (
423 meshDescription.lookup("blocks"),
424 block::iNew(meshDescription, vertices_, edges_, faces_)
425 );
426
427 transfer(blocks);
428 }
429
430
431 // Create patch information
432
433 faceListList tmpBlocksPatches;
435 PtrList<dictionary> patchDicts;
436
437
438 if (verbose_)
439 {
440 Info<< nl << "Creating topology patches - ";
441 }
442
443 if (meshDescription.found("patches"))
444 {
445 if (verbose_)
446 {
447 Info<< "from patches section" << endl;
448 }
449
451 wordList nbrPatchNames;
452
453 readPatches
454 (
455 meshDescription,
456 tmpBlocksPatches,
459 nbrPatchNames
460 );
461
462 if (verbose_)
463 {
464 Info<< nl
465 << "Reading physicalType from existing boundary file" << endl;
466 }
467
468 patchDicts.resize(patchNames.size());
469
471 (
472 meshDescription.time(),
473 meshDescription.time().constant(),
477 defaultPatchName,
478 defaultPatchType
479 );
480
481
482 // Add cyclic info (might not be present from older file)
483 forAll(patchDicts, patchi)
484 {
485 if (!patchDicts.set(patchi))
486 {
487 patchDicts.set(patchi, new dictionary());
488 }
489
490 dictionary& dict = patchDicts[patchi];
491
492 // Add but not override type
493 if (!dict.found("type"))
494 {
495 dict.add("type", patchTypes[patchi], false);
496 }
497 else if (dict.get<word>("type") != patchTypes[patchi])
498 {
499 FatalIOErrorInFunction(meshDescription)
500 << "For patch " << patchNames[patchi]
501 << " overriding type '" << patchTypes[patchi]
502 << "' with '" << dict.get<word>("type")
503 << "' (read from boundary file)"
504 << exit(FatalIOError);
505 }
506
507 // Override neighbour patch name
508 if (!nbrPatchNames[patchi].empty())
509 {
510 dict.set("neighbourPatch", nbrPatchNames[patchi]);
511 }
512 }
513 }
514 else if (meshDescription.found("boundary"))
515 {
516 if (verbose_)
517 {
518 Info<< "from boundary section" << endl;
519 }
520
521 readBoundary
522 (
523 meshDescription,
525 tmpBlocksPatches,
527 );
528 }
529 else
530 {
531 if (verbose_)
532 {
533 Info<< "with default boundary only!!" << endl;
534 }
535 }
536
537
538 if (verbose_)
539 {
540 Info<< nl << "Creating block mesh topology";
541 if (hasPointTransforms())
542 {
543 Info<< " - scaling/transform applied later";
544 }
545 Info<< endl;
546 }
547
548 auto blockMeshPtr = autoPtr<polyMesh>::New
549 (
550 IOobject
551 (
553 meshDescription.time().constant(),
554 meshDescription.time(),
557 false
558 ),
559 pointField(vertices_), // Use a copy of vertices
560 getBlockShapes(),
561 tmpBlocksPatches,
564 defaultPatchName,
565 defaultPatchType
566 );
567
568 check(*blockMeshPtr, meshDescription);
569
570 return blockMeshPtr;
571}
572
573
574// ************************************************************************* //
void setSize(const label n)
Alias for resize()
Definition: List.H:218
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
wordList patchNames() const
Return the patch names.
Definition: blockMesh.C:399
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:640
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition: dictionary.C:780
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:324
@ END_LIST
End list [isseparator].
Definition: token.H:156
A class for handling words, derived from Foam::string.
Definition: word.H:68
Foam::word regionName(Foam::polyMesh::defaultRegion)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
const label nPatches
label nPoints
Namespace for OpenFOAM.
List< cellShape > cellShapeList
List of cellShapes and PtrList of List of cellShape.
Definition: cellShapeList.H:45
List< word > wordList
A List of words.
Definition: fileName.H:63
PtrList< blockFace > blockFaceList
A PtrList of blockFaces.
Definition: blockFaceList.H:47
messageStream Info
Information stream (stdout output on master, null elsewhere)
static void check(const int retVal, const char *what)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
PtrList< block > blockList
A PtrList of blocks.
Definition: blockList.H:47
PtrList< blockEdge > blockEdgeList
A PtrList of blockEdges.
Definition: blockEdgeList.H:47
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
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.
IOerror FatalIOError
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
List< faceList > faceListList
A List of faceList.
Definition: faceListFwd.H:49
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
static void rewritePatchLabels(const Source &source, const word &patchName, const PtrList< block > &blocks, const label nPoints, faceList &patchFaces)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
preservePatchTypes
wordList patchTypes(nPatches)
wordList patchNames(nPatches)
PtrList< dictionary > patchDicts
Definition: readKivaGrid.H:532
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
static const char *const typeName
The type name used in ensight case files.