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 | Copyright (C) 2011-2015 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 "Time.H"
28 #include "preservePatchTypes.H"
29 #include "emptyPolyPatch.H"
30 #include "cyclicPolyPatch.H"
31 
32 
34 (
35  const dictionary& meshDescription,
36  faceListList& tmpBlocksPatches,
39  wordList& nbrPatchNames
40 )
41 {
42  bool topologyOK = true;
43 
44  ITstream& patchStream(meshDescription.lookup("patches"));
45 
46  // read number of patches in mesh
47  label nPatches = 0;
48 
49  token firstToken(patchStream);
50 
51  if (firstToken.isLabel())
52  {
53  nPatches = firstToken.labelToken();
54 
55  tmpBlocksPatches.setSize(nPatches);
56  patchNames.setSize(nPatches);
57  patchTypes.setSize(nPatches);
58  nbrPatchNames.setSize(nPatches);
59  }
60  else
61  {
62  patchStream.putBack(firstToken);
63  }
64 
65  // Read beginning of blocks
66  patchStream.readBegin("patches");
67 
68  nPatches = 0;
69 
70  token lastToken(patchStream);
71  while
72  (
73  !(
74  lastToken.isPunctuation()
75  && lastToken.pToken() == token::END_LIST
76  )
77  )
78  {
79  if (tmpBlocksPatches.size() <= nPatches)
80  {
81  tmpBlocksPatches.setSize(nPatches + 1);
82  patchNames.setSize(nPatches + 1);
83  patchTypes.setSize(nPatches + 1);
84  nbrPatchNames.setSize(nPatches + 1);
85  }
86 
87  patchStream.putBack(lastToken);
88 
89  patchStream
91  >> patchNames[nPatches];
92 
93  // Read patch faces
94  patchStream >> tmpBlocksPatches[nPatches];
95 
96 
97  // Catch multiple patches asap.
98  for (label i = 0; i < nPatches; i++)
99  {
100  if (patchNames[nPatches] == patchNames[i])
101  {
103  << "Duplicate patch " << patchNames[nPatches]
104  << " at line " << patchStream.lineNumber()
105  << ". Exiting !" << nl
106  << exit(FatalError);
107  }
108  }
109 
110  topologyOK = topologyOK && patchLabelsOK
111  (
112  nPatches,
113  blockPointField_,
114  tmpBlocksPatches[nPatches]
115  );
116 
117  nPatches++;
118 
119 
120  // Split old style cyclics
121 
122  if (patchTypes[nPatches-1] == cyclicPolyPatch::typeName)
123  {
124  word halfA = patchNames[nPatches-1] + "_half0";
125  word halfB = patchNames[nPatches-1] + "_half1";
126 
128  << "Old-style cyclic definition."
129  << " Splitting patch "
130  << patchNames[nPatches-1] << " into two halves "
131  << halfA << " and " << halfB << endl
132  << " Alternatively use new 'boundary' dictionary syntax."
133  << endl;
134 
135  // Add extra patch
136  if (tmpBlocksPatches.size() <= nPatches)
137  {
138  tmpBlocksPatches.setSize(nPatches + 1);
139  patchNames.setSize(nPatches + 1);
140  patchTypes.setSize(nPatches + 1);
141  nbrPatchNames.setSize(nPatches + 1);
142  }
143 
144  // Update halfA info
145  patchNames[nPatches-1] = halfA;
146  nbrPatchNames[nPatches-1] = halfB;
147  // Update halfB info
149  patchNames[nPatches] = halfB;
150  nbrPatchNames[nPatches] = halfA;
151 
152  // Split faces
153  if ((tmpBlocksPatches[nPatches-1].size() % 2) != 0)
154  {
156  << "Size of cyclic faces is not a multiple of 2 :"
157  << tmpBlocksPatches[nPatches-1]
158  << exit(FatalError);
159  }
160  label sz = tmpBlocksPatches[nPatches-1].size()/2;
161  faceList unsplitFaces(tmpBlocksPatches[nPatches-1], true);
162  tmpBlocksPatches[nPatches-1] = faceList
163  (
164  SubList<face>(unsplitFaces, sz)
165  );
166  tmpBlocksPatches[nPatches] = faceList
167  (
168  SubList<face>(unsplitFaces, sz, sz)
169  );
170 
171  nPatches++;
172  }
173 
174  patchStream >> lastToken;
175  }
176  patchStream.putBack(lastToken);
177 
178  // Read end of blocks
179  patchStream.readEnd("patches");
180 
181  return topologyOK;
182 }
183 
184 
186 (
187  const dictionary& meshDescription,
189  faceListList& tmpBlocksPatches,
191 )
192 {
193  bool topologyOK = true;
194 
195  // Read like boundary file
196  const PtrList<entry> patchesInfo
197  (
198  meshDescription.lookup("boundary")
199  );
200 
201  patchNames.setSize(patchesInfo.size());
202  tmpBlocksPatches.setSize(patchesInfo.size());
203  patchDicts.setSize(patchesInfo.size());
204 
205  forAll(tmpBlocksPatches, patchI)
206  {
207  const entry& patchInfo = patchesInfo[patchI];
208 
209  if (!patchInfo.isDict())
210  {
211  FatalIOErrorInFunction(meshDescription)
212  << "Entry " << patchInfo << " in boundary section is not a"
213  << " valid dictionary." << exit(FatalIOError);
214  }
215 
216  patchNames[patchI] = patchInfo.keyword();
217  // Construct dictionary
218  patchDicts.set(patchI, new dictionary(patchInfo.dict()));
219  // Read block faces
220  patchDicts[patchI].lookup("faces") >> tmpBlocksPatches[patchI];
221 
222  topologyOK = topologyOK && patchLabelsOK
223  (
224  patchI,
225  blockPointField_,
226  tmpBlocksPatches[patchI]
227  );
228  }
229 
230  return topologyOK;
231 }
232 
233 
235 (
236  cellShapeList& tmpBlockCells
237 )
238 {
239  const blockMesh& blocks = *this;
240 
241  tmpBlockCells.setSize(blocks.size());
242  forAll(blocks, blockI)
243  {
244  tmpBlockCells[blockI] = cellShape(blocks[blockI].blockShape());
245 
246  if (tmpBlockCells[blockI].mag(blockPointField_) < 0.0)
247  {
249  << "negative volume block : " << blockI
250  << ", probably defined inside-out" << endl;
251  }
252  }
253 }
254 
255 
256 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
257 
259 (
260  const IOdictionary& meshDescription,
261  const word& regionName
262 )
263 {
264  bool topologyOK = true;
265 
266  blockList& blocks = *this;
267 
268  word defaultPatchName = "defaultFaces";
269  word defaultPatchType = emptyPolyPatch::typeName;
270 
271  // get names/types for the unassigned patch faces
272  // this is a bit heavy handed (and ugly), but there is currently
273  // no easy way to rename polyMesh patches subsequently
274  if (const dictionary* dictPtr = meshDescription.subDictPtr("defaultPatch"))
275  {
276  dictPtr->readIfPresent("name", defaultPatchName);
277  dictPtr->readIfPresent("type", defaultPatchType);
278  }
279 
280  // optional 'convertToMeters' or 'scale' scaling factor
281  if (!meshDescription.readIfPresent("convertToMeters", scaleFactor_))
282  {
283  meshDescription.readIfPresent("scale", scaleFactor_);
284  }
285 
286 
287  //
288  // get the non-linear edges in mesh
289  //
290  if (meshDescription.found("edges"))
291  {
292  if (verboseOutput)
293  {
294  Info<< "Creating curved edges" << endl;
295  }
296 
297  ITstream& is(meshDescription.lookup("edges"));
298 
299  // read number of edges in mesh
300  label nEdges = 0;
301 
302  token firstToken(is);
303 
304  if (firstToken.isLabel())
305  {
306  nEdges = firstToken.labelToken();
307  edges_.setSize(nEdges);
308  }
309  else
310  {
311  is.putBack(firstToken);
312  }
313 
314  // Read beginning of edges
315  is.readBegin("edges");
316 
317  nEdges = 0;
318 
319  token lastToken(is);
320  while
321  (
322  !(
323  lastToken.isPunctuation()
324  && lastToken.pToken() == token::END_LIST
325  )
326  )
327  {
328  if (edges_.size() <= nEdges)
329  {
330  edges_.setSize(nEdges + 1);
331  }
332 
333  is.putBack(lastToken);
334 
335  edges_.set
336  (
337  nEdges,
338  curvedEdge::New(blockPointField_, is)
339  );
340 
341  nEdges++;
342 
343  is >> lastToken;
344  }
345  is.putBack(lastToken);
346 
347  // Read end of edges
348  is.readEnd("edges");
349  }
350  else if (verboseOutput)
351  {
352  Info<< "No non-linear edges defined" << endl;
353  }
354 
355 
356  //
357  // Create the blocks
358  //
359  if (verboseOutput)
360  {
361  Info<< "Creating topology blocks" << endl;
362  }
363 
364  {
365  ITstream& is(meshDescription.lookup("blocks"));
366 
367  // read number of blocks in mesh
368  label nBlocks = 0;
369 
370  token firstToken(is);
371 
372  if (firstToken.isLabel())
373  {
374  nBlocks = firstToken.labelToken();
375  blocks.setSize(nBlocks);
376  }
377  else
378  {
379  is.putBack(firstToken);
380  }
381 
382  // Read beginning of blocks
383  is.readBegin("blocks");
384 
385  nBlocks = 0;
386 
387  token lastToken(is);
388  while
389  (
390  !(
391  lastToken.isPunctuation()
392  && lastToken.pToken() == token::END_LIST
393  )
394  )
395  {
396  if (blocks.size() <= nBlocks)
397  {
398  blocks.setSize(nBlocks + 1);
399  }
400 
401  is.putBack(lastToken);
402 
403  blocks.set
404  (
405  nBlocks,
406  new block
407  (
408  blockPointField_,
409  edges_,
410  is
411  )
412  );
413 
414  topologyOK = topologyOK && blockLabelsOK
415  (
416  nBlocks,
417  blockPointField_,
418  blocks[nBlocks].blockShape()
419  );
420 
421  nBlocks++;
422 
423  is >> lastToken;
424  }
425  is.putBack(lastToken);
426 
427  // Read end of blocks
428  is.readEnd("blocks");
429  }
430 
431 
432  polyMesh* blockMeshPtr = NULL;
433 
434  //
435  // Create the patches
436  //
437  if (verboseOutput)
438  {
439  Info<< "Creating topology patches" << endl;
440  }
441 
442  if (meshDescription.found("patches"))
443  {
444  Info<< nl << "Reading patches section" << endl;
445 
446  faceListList tmpBlocksPatches;
449  wordList nbrPatchNames;
450 
451  topologyOK = topologyOK && readPatches
452  (
453  meshDescription,
454  tmpBlocksPatches,
455  patchNames,
456  patchTypes,
457  nbrPatchNames
458  );
459 
460  if (!topologyOK)
461  {
463  << "Cannot create mesh due to errors in topology, exiting !"
464  << nl << exit(FatalError);
465  }
466 
467  Info<< nl << "Creating block mesh topology" << endl;
468 
469  cellShapeList tmpBlockCells(blocks.size());
470  createCellShapes(tmpBlockCells);
471 
472 
473  Info<< nl << "Reading physicalType from existing boundary file" << endl;
474 
477 
479  (
480  meshDescription.time(),
481  meshDescription.time().constant(),
482  polyMesh::meshSubDir,
483  patchNames,
484  patchDicts,
485  defaultPatchName,
486  defaultPatchType
487  );
488 
489 
490  // Add cyclic info (might not be present from older file)
491  forAll(patchDicts, patchI)
492  {
493  if (!patchDicts.set(patchI))
494  {
495  patchDicts.set(patchI, new dictionary());
496  }
497 
498  dictionary& dict = patchDicts[patchI];
499 
500  // Add but not override type
501  if (!dict.found("type"))
502  {
503  dict.add("type", patchTypes[patchI], false);
504  }
505  else if (word(dict.lookup("type")) != patchTypes[patchI])
506  {
508  (
509  meshDescription
510  ) << "For patch " << patchNames[patchI]
511  << " overriding type '" << patchTypes[patchI]
512  << "' with '" << word(dict.lookup("type"))
513  << "' (read from boundary file)"
514  << endl;
515  }
516 
517  // Override neighbourpatch name
518  if (nbrPatchNames[patchI] != word::null)
519  {
520  dict.set("neighbourPatch", nbrPatchNames[patchI]);
521  }
522  }
523 
524 
525  blockMeshPtr = new polyMesh
526  (
527  IOobject
528  (
529  regionName,
530  meshDescription.time().constant(),
531  meshDescription.time(),
532  IOobject::NO_READ,
533  IOobject::NO_WRITE,
534  false
535  ),
536  xferCopy(blockPointField_), // copy these points, do NOT move
537  tmpBlockCells,
538  tmpBlocksPatches,
539  patchNames,
540  patchDicts,
541  defaultPatchName,
542  defaultPatchType
543  );
544  }
545  else if (meshDescription.found("boundary"))
546  {
548  faceListList tmpBlocksPatches;
550 
551  topologyOK = topologyOK && readBoundary
552  (
553  meshDescription,
554  patchNames,
555  tmpBlocksPatches,
556  patchDicts
557  );
558 
559  if (!topologyOK)
560  {
562  << "Cannot create mesh due to errors in topology, exiting !"
563  << nl << exit(FatalError);
564  }
565 
566 
567  Info<< nl << "Creating block mesh topology" << endl;
568 
569  cellShapeList tmpBlockCells(blocks.size());
570  createCellShapes(tmpBlockCells);
571 
572  // Extract
573 
574  blockMeshPtr = new polyMesh
575  (
576  IOobject
577  (
578  regionName,
579  meshDescription.time().constant(),
580  meshDescription.time(),
581  IOobject::NO_READ,
582  IOobject::NO_WRITE,
583  false
584  ),
585  xferCopy(blockPointField_), // copy these points, do NOT move
586  tmpBlockCells,
587  tmpBlocksPatches,
588  patchNames,
589  patchDicts,
590  defaultPatchName,
591  defaultPatchType
592  );
593  }
594 
595  checkBlockMesh(*blockMeshPtr);
596 
597  return blockMeshPtr;
598 }
599 
600 
601 // ************************************************************************* //
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:65
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::Istream::readEnd
Istream & readEnd(const char *funcName)
Definition: Istream.C:105
Foam::block
Creates a single block of cells from point coordinates, numbers of cells in each direction and an exp...
Definition: block.H:63
Foam::blockMesh::readPatches
bool readPatches(const dictionary &meshDescription, faceListList &tmpBlocksPatches, wordList &patchNames, wordList &patchTypes, wordList &nbrPatchNames)
Definition: blockMeshTopology.C:34
preservePatchTypes
preservePatchTypes(runTime, runTime.constant(), polyMesh::meshSubDir, patchNames, patchDicts, defaultFacesName, defaultFacesType)
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::entry::keyword
const keyType & keyword() const
Return keyword.
Definition: entry.H:120
nPatches
label nPatches
Definition: readKivaGrid.H:402
cyclicPolyPatch.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::blockMesh::readBoundary
bool readBoundary(const dictionary &meshDescription, wordList &patchNames, faceListList &tmpBlocksPatches, PtrList< dictionary > &patchDicts)
Definition: blockMeshTopology.C:186
Foam::dictionary::readIfPresent
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
Definition: dictionaryTemplates.C:94
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
preservePatchTypes.H
preservePatchTypes
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::token
A token holds items read from Istream.
Definition: token.H:67
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
blockMesh.H
Foam::IOobject::time
const Time & time() const
Return time.
Definition: IOobject.C:245
defaultFacesType
word defaultFacesType
Definition: readKivaGrid.H:461
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
patchDicts
PtrList< dictionary > patchDicts
Definition: readKivaGrid.H:537
Foam::token::isLabel
bool isLabel() const
Definition: tokenI.H:262
patchTypes
wordList patchTypes(nPatches)
Foam::blockMesh::createCellShapes
void createCellShapes(cellShapeList &tmpBlockCells)
Definition: blockMeshTopology.C:235
regionName
Foam::word regionName
Definition: createNamedDynamicFvMesh.H:1
Foam::PtrList::set
bool set(const label) const
Is element set.
Foam::label
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
Foam::dictionary::found
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:304
Foam::entry::isDict
virtual bool isDict() const
Return true if this entry is a dictionary.
Definition: entry.H:153
Foam::ITstream
Input token stream.
Definition: ITstream.H:49
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
patchNames
wordList patchNames(nPatches)
Foam::token::pToken
punctuationToken pToken() const
Definition: tokenI.H:208
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::cellShape
An analytical geometric cellShape.
Definition: cellShape.H:69
Foam::entry::dict
virtual const dictionary & dict() const =0
Return dictionary if this entry is a dictionary.
emptyPolyPatch.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::List::setSize
void setSize(const label)
Reset size of List.
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::faceList
List< face > faceList
Definition: faceListFwd.H:43
Foam::xferCopy
Xfer< T > xferCopy(const T &)
Construct by copying the contents of the arg.
Foam::List< faceList >
Foam::Istream::putBack
void putBack(const token &)
Put back token.
Definition: Istream.C:30
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePaths.H:130
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
Foam::PtrList::setSize
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:142
Foam::IOstream::lineNumber
label lineNumber() const
Return current stream line number.
Definition: IOstream.H:438
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:330
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
IOWarningInFunction
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
Definition: messageStream.H:271
Foam::Istream::readBegin
Istream & readBegin(const char *funcName)
Definition: Istream.C:88
Foam::token::isPunctuation
bool isPunctuation() const
Definition: tokenI.H:203
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::blockMesh::createTopology
polyMesh * createTopology(const IOdictionary &, const word &regionName)
Definition: blockMeshTopology.C:259
Foam::token::labelToken
label labelToken() const
Definition: tokenI.H:267
Foam::dictionary::subDictPtr
const dictionary * subDictPtr(const word &) const
Find and return a sub-dictionary pointer if present.
Definition: dictionary.C:616
Foam::blockMesh
A multi-block mesh generator.
Definition: blockMesh.H:58