blockMeshTopology.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
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 "blockMesh.H"
27 #include "blockMeshTools.H"
28 #include "Time.H"
29 #include "preservePatchTypes.H"
30 #include "emptyPolyPatch.H"
31 #include "cyclicPolyPatch.H"
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 template<class Source>
36 void Foam::blockMesh::checkPatchLabels
37 (
38  const Source& source,
39  const word& patchName,
40  const pointField& points,
41  faceList& patchFaces
42 ) const
43 {
44  forAll(patchFaces, facei)
45  {
46  face& f = patchFaces[facei];
47 
48  // Replace (<block> <face>) face description
49  // with the corresponding block face
50  if (f.size() == 2)
51  {
52  const label bi = f[0];
53  const label fi = f[1];
54 
55  if (bi >= size())
56  {
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
62  << exit(FatalIOError);
63  }
64  else if (fi >= operator[](bi).blockShape().faces().size())
65  {
67  << "Block face index out of range for patch face " << f
68  << nl
69  << " Number of block faces = "
70  << operator[](bi).blockShape().faces().size()
71  << ", index = " << f[1] << nl
72  << " on patch " << patchName << ", face " << facei
73  << exit(FatalIOError);
74  }
75  else
76  {
77  f = operator[](bi).blockShape().faces()[fi];
78  }
79  }
80  else
81  {
82  forAll(f, fp)
83  {
84  if (f[fp] < 0)
85  {
87  << "Negative point label " << f[fp] << nl
88  << " on patch " << patchName << ", face " << facei
89  << exit(FatalIOError);
90  }
91  else if (f[fp] >= points.size())
92  {
94  << "Point label " << f[fp]
95  << " out of range 0.." << points.size() - 1 << nl
96  << " on patch " << patchName << ", face " << facei
97  << exit(FatalIOError);
98  }
99  }
100  }
101  }
102 }
103 
104 
105 void Foam::blockMesh::readPatches
106 (
107  const dictionary& meshDescription,
108  faceListList& tmpBlocksPatches,
111  wordList& nbrPatchNames
112 )
113 {
114  // Collect all variables
115  dictionary varDict(meshDescription.subOrEmptyDict("namedVertices"));
116  varDict.merge(meshDescription.subOrEmptyDict("namedBlocks"));
117 
118 
119  ITstream& patchStream(meshDescription.lookup("patches"));
120 
121  // Read number of patches in mesh
122  label nPatches = 0;
123 
124  token firstToken(patchStream);
125 
126  if (firstToken.isLabel())
127  {
128  nPatches = firstToken.labelToken();
129 
130  tmpBlocksPatches.setSize(nPatches);
133  nbrPatchNames.setSize(nPatches);
134  }
135  else
136  {
137  patchStream.putBack(firstToken);
138  }
139 
140  // Read beginning of blocks
141  patchStream.readBegin("patches");
142 
143  nPatches = 0;
144 
145  token lastToken(patchStream);
146  while
147  (
148  !(
149  lastToken.isPunctuation()
150  && lastToken.pToken() == token::END_LIST
151  )
152  )
153  {
154  if (tmpBlocksPatches.size() <= nPatches)
155  {
156  tmpBlocksPatches.setSize(nPatches + 1);
159  nbrPatchNames.setSize(nPatches + 1);
160  }
161 
162  patchStream.putBack(lastToken);
163 
164  patchStream
165  >> patchTypes[nPatches]
166  >> patchNames[nPatches];
167 
168  // Read patch faces
169  tmpBlocksPatches[nPatches] = blockMeshTools::read<face>
170  (
171  patchStream,
172  varDict
173  );
174 
175 
176  // Check for multiple patches
177  for (label i = 0; i < nPatches; i++)
178  {
179  if (patchNames[nPatches] == patchNames[i])
180  {
181  FatalIOErrorInFunction(patchStream)
182  << "Duplicate patch " << patchNames[nPatches]
183  << " at line " << patchStream.lineNumber()
184  << exit(FatalIOError);
185  }
186  }
187 
188  checkPatchLabels
189  (
190  patchStream,
192  vertices_,
193  tmpBlocksPatches[nPatches]
194  );
195 
196  nPatches++;
197 
198 
199  // Split old style cyclics
200 
201  if (patchTypes[nPatches-1] == cyclicPolyPatch::typeName)
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 
263 void 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.lookupOrDefault("boundary", PtrList<entry>())
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 
311  checkPatchLabels
312  (
313  patchInfo.dict(),
315  vertices_,
316  tmpBlocksPatches[patchi]
317  );
318  }
319 }
320 
321 
322 void Foam::blockMesh::createCellShapes
323 (
324  cellShapeList& tmpBlockCells
325 )
326 {
327  const blockMesh& blocks = *this;
328 
329  tmpBlockCells.setSize(blocks.size());
330  forAll(blocks, blocki)
331  {
332  tmpBlockCells[blocki] = blocks[blocki].blockShape();
333  }
334 }
335 
336 
337 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
338 
339 void Foam::blockMesh::defaultPatchError
340 (
341  const word& defaultPatchName,
342  const dictionary& meshDescription
343 ) const
344 {
345  // This will be made a fatal error in a future release
346  // to ensure that the defaultPatch is handled explicitly
347  // avoiding common user errors
348  // FatalIOErrorInFunction(meshDescription)
349  IOWarningInFunction(meshDescription)
350  << "The 'defaultPatch' type must be specified"
351  " for the '" << defaultPatchName << "' patch, e.g. for snappyHexMesh"
352  << nl << nl
353  << " defaultPatch\n"
354  " {\n"
355  " name default; // optional\n"
356  " type patch;\n"
357  " }\n\n"
358  << "or for 2D meshes" << nl << nl
359  << " defaultPatch\n"
360  " {\n"
361  " name frontAndBack; // optional\n"
362  " type empty;\n"
363  " }\n"
364  << endl;
365  // << exit(FatalIOError);
366 }
367 
368 
369 Foam::polyMesh* Foam::blockMesh::createTopology
370 (
371  const IOdictionary& meshDescription,
372  const word& regionName
373 )
374 {
375  checkBlockFaceOrientation = meshDescription.lookupOrDefault<Switch>
376  (
377  "checkBlockFaceOrientation",
378  true
379  );
380 
381  blockList& blocks = *this;
382 
383  word defaultPatchName = "defaultFaces";
384  word defaultPatchType = emptyPolyPatch::typeName;
385  bool defaultPatchTypeSet = false;
386 
387  // Read the names/types for the unassigned patch faces
388  // this is a bit heavy handed (and ugly), but there is currently
389  // no easy way to rename polyMesh patches subsequently
390  if (const dictionary* dictPtr = meshDescription.subDictPtr("defaultPatch"))
391  {
392  dictPtr->readIfPresent("name", defaultPatchName);
393  defaultPatchType = dictPtr->lookup<word>("type");
394  defaultPatchTypeSet = true;
395  }
396 
397  // Optional 'convertToMeters' or 'scale' scaling factor
398  if (!meshDescription.readIfPresent("convertToMeters", scaleFactor_))
399  {
400  meshDescription.readIfPresent("scale", scaleFactor_);
401  }
402 
403  // Read the block edges
404  if (meshDescription.found("edges"))
405  {
406  if (verboseOutput)
407  {
408  Info<< "Creating block edges" << endl;
409  }
410 
411  blockEdgeList edges
412  (
413  meshDescription.lookup("edges"),
414  blockEdge::iNew(meshDescription, geometry_, vertices_)
415  );
416 
417  edges_.transfer(edges);
418  }
419  else if (verboseOutput)
420  {
421  Info<< "No non-linear block edges defined" << endl;
422  }
423 
424 
425  // Read the block faces
426  if (meshDescription.found("faces"))
427  {
428  if (verboseOutput)
429  {
430  Info<< "Creating block faces" << endl;
431  }
432 
433  blockFaceList faces
434  (
435  meshDescription.lookup("faces"),
436  blockFace::iNew(meshDescription, geometry_)
437  );
438 
439  faces_.transfer(faces);
440  }
441  else if (verboseOutput)
442  {
443  Info<< "No non-planar block faces defined" << endl;
444  }
445 
446 
447  // Create the blocks
448  if (verboseOutput)
449  {
450  Info<< "Creating topology blocks" << endl;
451  }
452  {
453  blockList blocks
454  (
455  meshDescription.lookup("blocks"),
456  block::iNew(meshDescription, vertices_, edges_, faces_)
457  );
458 
459  transfer(blocks);
460  }
461 
462 
463 
464  polyMesh* blockMeshPtr = nullptr;
465 
466  // Create the patches
467 
468  if (verboseOutput)
469  {
470  Info<< "Creating topology patches" << endl;
471  }
472 
473  if (meshDescription.found("patches"))
474  {
475  IOWarningInFunction(meshDescription)
476  << "Reading patches section" << nl
477  << " The 'patches' entry is deprecated and has been superseded"
478  << " by the more consistent and general 'boundary' entry."
479  << endl;
480 
481  faceListList tmpBlocksPatches;
484  wordList nbrPatchNames;
485 
486  readPatches
487  (
488  meshDescription,
489  tmpBlocksPatches,
490  patchNames,
491  patchTypes,
492  nbrPatchNames
493  );
494 
495  Info<< nl << "Creating block mesh topology" << endl;
496 
497  cellShapeList tmpBlockCells(blocks.size());
498  createCellShapes(tmpBlockCells);
499 
500  Info<< nl << "Reading physicalType from existing boundary file" << endl;
501 
502  PtrList<dictionary> patchDicts(patchNames.size());
503  word defaultFacesType;
504 
506  (
507  meshDescription.time(),
508  meshDescription.time().constant(),
510  patchNames,
511  patchDicts,
512  defaultPatchName,
513  defaultPatchType
514  );
515 
516  // Add cyclic info (might not be present from older file)
518  {
519  if (!patchDicts.set(patchi))
520  {
521  patchDicts.set(patchi, new dictionary());
522  }
523 
524  dictionary& dict = patchDicts[patchi];
525 
526  // Add but not override type
527  if (!dict.found("type"))
528  {
529  dict.add("type", patchTypes[patchi], false);
530  }
531  else if (word(dict.lookup("type")) != patchTypes[patchi])
532  {
533  FatalIOErrorInFunction(meshDescription)
534  << "For patch " << patchNames[patchi]
535  << " overriding type '" << patchTypes[patchi]
536  << "' with '" << word(dict.lookup("type"))
537  << "' (read from boundary file)"
538  << exit(FatalIOError);
539  }
540 
541  // Override neighbourpatch name
542  if (nbrPatchNames[patchi] != word::null)
543  {
544  dict.set("neighbourPatch", nbrPatchNames[patchi]);
545  }
546  }
547 
548  blockMeshPtr = new polyMesh
549  (
550  IOobject
551  (
552  regionName,
553  meshDescription.time().constant(),
554  meshDescription.time(),
557  false
558  ),
559  clone(vertices_), // Copy these points, do NOT move
560  tmpBlockCells,
561  tmpBlocksPatches,
562  patchNames,
563  patchDicts,
564  defaultPatchName,
565  defaultPatchType
566  );
567 
568  if
569  (
570  !defaultPatchTypeSet
571  && blockMeshPtr->boundaryMesh().findPatchID(defaultPatchName) != -1
572  )
573  {
574  defaultPatchError(defaultPatchName, meshDescription);
575  }
576  }
577  else
578  {
580  faceListList tmpBlocksPatches;
581  PtrList<dictionary> patchDicts;
582 
583  readBoundary
584  (
585  meshDescription,
586  patchNames,
587  tmpBlocksPatches,
588  patchDicts
589  );
590 
591  Info<< nl << "Creating block mesh topology" << endl;
592 
593  cellShapeList tmpBlockCells(blocks.size());
594  createCellShapes(tmpBlockCells);
595 
596  blockMeshPtr = new polyMesh
597  (
598  IOobject
599  (
600  regionName,
601  meshDescription.time().constant(),
602  meshDescription.time(),
605  false
606  ),
607  clone(vertices_), // Copy these points, do NOT move
608  tmpBlockCells,
609  tmpBlocksPatches,
610  patchNames,
611  patchDicts,
612  defaultPatchName,
613  defaultPatchType
614  );
615 
616  if
617  (
618  !defaultPatchTypeSet
619  && blockMeshPtr->boundaryMesh().findPatchID(defaultPatchName) != -1
620  )
621  {
622  defaultPatchError(defaultPatchName, meshDescription);
623  }
624  }
625 
626  check(*blockMeshPtr, meshDescription);
627 
628  return blockMeshPtr;
629 }
630 
631 
632 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
const T & operator[](const label) const
Return element const reference.
Definition: UPtrListI.H:96
const blockFaceList & faces() const
Return the curved faces.
Definition: blockMesh.H:257
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.
Definition: dictionary.C:860
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:1307
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:1169
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:659
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:271
@ END_LIST
Definition: token.H:107
static const word null
An empty word.
Definition: word.H:77
Foam::word regionName
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
label patchi
const pointField & points
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
List< cellShape > cellShapeList
List of cellShapes and PtrList of List of cellShape.
Definition: cellShapeList.H:43
List< word > wordList
A List of words.
Definition: fileName.H:54
PtrList< blockFace > blockFaceList
A PtrList of blockFaces.
Definition: blockFaceList.H:45
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
PtrList< block > blockList
A PtrList of blocks.
Definition: blockList.H:45
PtrList< blockEdge > blockEdgeList
A PtrList of blockEdges.
Definition: blockEdgeList.H:45
List< faceList > faceListList
Definition: faceListFwd.H:45
T clone(const T &t)
Definition: List.H:55
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
Definition: faceListFwd.H:43
static const char nl
Definition: Ostream.H:260
preservePatchTypes
wordList patchTypes(nPatches)
wordList patchNames(nPatches)
PtrList< dictionary > patchDicts
Definition: readKivaGrid.H:531
labelList f(nPoints)
word defaultFacesType
Definition: readKivaGrid.H:455
label nPatches
Definition: readKivaGrid.H:396
dictionary dict