blockMesh.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-2017 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 Application
25  blockMesh
26 
27 Description
28  A multi-block mesh generator.
29 
30  Uses the block mesh description found in
31  - \c system/blockMeshDict
32  - \c system/<region>/blockMeshDict
33  - \c constant/polyMesh/blockMeshDict
34  - \c constant/<region>/polyMesh/blockMeshDict
35 
36 Usage
37  \b blockMesh [OPTION]
38 
39  Options:
40  - \par -blockTopology
41  Write the topology as a set of edges in OBJ format.
42 
43  - \par -region <name>
44  Specify an alternative mesh region.
45 
46  - \par -dict <filename>
47  Specify alternative dictionary for the block mesh description.
48 
49 \*---------------------------------------------------------------------------*/
50 
51 #include "Time.H"
52 #include "IOdictionary.H"
53 #include "IOPtrList.H"
54 
55 #include "blockMesh.H"
56 #include "attachPolyTopoChanger.H"
57 #include "emptyPolyPatch.H"
58 #include "cellSet.H"
59 
60 #include "argList.H"
61 #include "OSspecific.H"
62 #include "OFstream.H"
63 
64 #include "Pair.H"
65 #include "slidingInterface.H"
66 
67 using namespace Foam;
68 
69 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
70 
71 int main(int argc, char *argv[])
72 {
75  (
76  "blockTopology",
77  "write block edges and centres as .obj files"
78  );
80  (
81  "noClean",
82  "keep the existing files in the polyMesh"
83  );
85  (
86  "dict",
87  "file",
88  "specify alternative dictionary for the blockMesh description"
89  );
90 
92  (
93  "Block description\n"
94  "\n"
95  " For a given block, the correspondence between the ordering of\n"
96  " vertex labels and face labels is shown below.\n"
97  " For vertex numbering in the sequence 0 to 7 (block, centre):\n"
98  " faces 0 (f0) and 1 are left and right, respectively;\n"
99  " faces 2 and 3 are bottom and top;\n"
100  " and faces 4 and 5 are front the back:\n"
101  "\n"
102  " 4 ---- 5\n"
103  " f3 |\\ |\\ f5\n"
104  " | | 7 ---- 6 \\\n"
105  " | 0 |--- 1 | \\\n"
106  " | \\| \\| f4\n"
107  " f2 3 ---- 2\n"
108  "\n"
109  " f0 ----- f1\n"
110  );
111 
112  #include "addRegionOption.H"
113  #include "setRootCase.H"
114  #include "createTime.H"
115 
116  const word dictName("blockMeshDict");
117 
119  word regionPath;
120 
121  // Check if the region is specified otherwise mesh the default region
122  if (args.optionReadIfPresent("region", regionName, polyMesh::defaultRegion))
123  {
124  Info<< nl << "Generating mesh for region " << regionName << endl;
125  regionPath = regionName;
126  }
127 
128  // Search for the appropriate blockMesh dictionary....
129 
131 
132  // Check if the dictionary is specified on the command-line
133  if (args.optionFound("dict"))
134  {
135  dictPath = args["dict"];
136 
137  dictPath =
138  (
139  isDir(dictPath)
140  ? dictPath/dictName
141  : dictPath
142  );
143  }
144  // Check if dictionary is present in the constant directory
145  else if
146  (
147  exists
148  (
149  runTime.path()/runTime.constant()
150  /regionPath/polyMesh::meshSubDir/dictName
151  )
152  )
153  {
154  dictPath =
155  runTime.constant()
156  /regionPath/polyMesh::meshSubDir/dictName;
157  }
158  // Otherwise assume the dictionary is present in the system directory
159  else
160  {
161  dictPath = runTime.system()/regionPath/dictName;
162  }
163 
164  if (!args.optionFound("noClean"))
165  {
166  fileName polyMeshPath
167  (
168  runTime.path()/runTime.constant()/regionPath/polyMesh::meshSubDir
169  );
170 
171  if (exists(polyMeshPath))
172  {
173  if (exists(polyMeshPath/dictName))
174  {
175  Info<< "Not deleting polyMesh directory " << nl
176  << " " << polyMeshPath << nl
177  << " because it contains " << dictName << endl;
178  }
179  else
180  {
181  Info<< "Deleting polyMesh directory" << nl
182  << " " << polyMeshPath << endl;
183  rmDir(polyMeshPath);
184  }
185  }
186  }
187 
188  IOobject meshDictIO
189  (
190  dictPath,
191  runTime,
194  false
195  );
196 
197  if (!meshDictIO.typeHeaderOk<IOdictionary>(true))
198  {
200  << meshDictIO.objectPath()
201  << nl
202  << exit(FatalError);
203  }
204 
205  Info<< "Creating block mesh from\n "
206  << meshDictIO.objectPath() << endl;
207 
208  IOdictionary meshDict(meshDictIO);
209  blockMesh blocks(meshDict, regionName);
210 
211 
212  if (args.optionFound("blockTopology"))
213  {
214  // Write mesh as edges.
215  {
216  fileName objMeshFile("blockTopology.obj");
217 
218  OFstream str(runTime.path()/objMeshFile);
219 
220  Info<< nl << "Dumping block structure as Lightwave obj format"
221  << " to " << objMeshFile << endl;
222 
223  blocks.writeTopology(str);
224  }
225 
226  // Write centres of blocks
227  {
228  fileName objCcFile("blockCentres.obj");
229 
230  OFstream str(runTime.path()/objCcFile);
231 
232  Info<< nl << "Dumping block centres as Lightwave obj format"
233  << " to " << objCcFile << endl;
234 
235  const polyMesh& topo = blocks.topology();
236 
237  const pointField& cellCentres = topo.cellCentres();
238 
239  forAll(cellCentres, celli)
240  {
241  //point cc = b.blockShape().centre(b.points());
242  const point& cc = cellCentres[celli];
243 
244  str << "v " << cc.x() << ' ' << cc.y() << ' ' << cc.z() << nl;
245  }
246  }
247 
248  Info<< nl << "end" << endl;
249 
250  return 0;
251  }
252 
253 
254  Info<< nl << "Creating polyMesh from blockMesh" << endl;
255 
256  word defaultFacesName = "defaultFaces";
257  word defaultFacesType = emptyPolyPatch::typeName;
258  polyMesh mesh
259  (
260  IOobject
261  (
262  regionName,
263  runTime.constant(),
264  runTime
265  ),
266  xferCopy(blocks.points()), // could we re-use space?
267  blocks.cells(),
268  blocks.patches(),
269  blocks.patchNames(),
270  blocks.patchDicts(),
273  );
274 
275 
276  // Read in a list of dictionaries for the merge patch pairs
277  if (meshDict.found("mergePatchPairs"))
278  {
279  List<Pair<word>> mergePatchPairs
280  (
281  meshDict.lookup("mergePatchPairs")
282  );
283 
284  #include "mergePatchPairs.H"
285  }
286  else
287  {
288  Info<< nl << "There are no merge patch pairs edges" << endl;
289  }
290 
291 
292  // Set any cellZones (note: cell labelling unaffected by above
293  // mergePatchPairs)
294 
295  label nZones = blocks.numZonedBlocks();
296 
297  if (nZones > 0)
298  {
299  Info<< nl << "Adding cell zones" << endl;
300 
301  // Map from zoneName to cellZone index
302  HashTable<label> zoneMap(nZones);
303 
304  // Cells per zone.
305  List<DynamicList<label>> zoneCells(nZones);
306 
307  // Running cell counter
308  label celli = 0;
309 
310  // Largest zone so far
311  label freeZoneI = 0;
312 
313  forAll(blocks, blockI)
314  {
315  const block& b = blocks[blockI];
316  const List<FixedList<label, 8>> blockCells = b.cells();
317  const word& zoneName = b.zoneName();
318 
319  if (zoneName.size())
320  {
321  HashTable<label>::const_iterator iter = zoneMap.find(zoneName);
322 
323  label zoneI;
324 
325  if (iter == zoneMap.end())
326  {
327  zoneI = freeZoneI++;
328 
329  Info<< " " << zoneI << '\t' << zoneName << endl;
330 
331  zoneMap.insert(zoneName, zoneI);
332  }
333  else
334  {
335  zoneI = iter();
336  }
337 
338  forAll(blockCells, i)
339  {
340  zoneCells[zoneI].append(celli++);
341  }
342  }
343  else
344  {
345  celli += b.cells().size();
346  }
347  }
348 
349 
350  List<cellZone*> cz(zoneMap.size());
351 
352  Info<< nl << "Writing cell zones as cellSets" << endl;
353 
354  forAllConstIter(HashTable<label>, zoneMap, iter)
355  {
356  label zoneI = iter();
357 
358  cz[zoneI] = new cellZone
359  (
360  iter.key(),
361  zoneCells[zoneI].shrink(),
362  zoneI,
363  mesh.cellZones()
364  );
365 
366  // Write as cellSet for ease of processing
367  cellSet cset(mesh, iter.key(), zoneCells[zoneI].shrink());
368  cset.write();
369  }
370 
371  mesh.pointZones().setSize(0);
372  mesh.faceZones().setSize(0);
373  mesh.cellZones().setSize(0);
375  }
376 
377  // Set the precision of the points data to 10
379 
380  Info<< nl << "Writing polyMesh" << endl;
381  mesh.removeFiles();
382  if (!mesh.write())
383  {
385  << "Failed writing polyMesh."
386  << exit(FatalError);
387  }
388 
389 
390  //
391  // write some information
392  //
393  {
394  const polyPatchList& patches = mesh.boundaryMesh();
395 
396  Info<< "----------------" << nl
397  << "Mesh Information" << nl
398  << "----------------" << nl
399  << " " << "boundingBox: " << boundBox(mesh.points()) << nl
400  << " " << "nPoints: " << mesh.nPoints() << nl
401  << " " << "nCells: " << mesh.nCells() << nl
402  << " " << "nFaces: " << mesh.nFaces() << nl
403  << " " << "nInternalFaces: " << mesh.nInternalFaces() << nl;
404 
405  Info<< "----------------" << nl
406  << "Patches" << nl
407  << "----------------" << nl;
408 
409  forAll(patches, patchi)
410  {
411  const polyPatch& p = patches[patchi];
412 
413  Info<< " " << "patch " << patchi
414  << " (start: " << p.start()
415  << " size: " << p.size()
416  << ") name: " << p.name()
417  << nl;
418  }
419  }
420 
421  Info<< "\nEnd\n" << endl;
422 
423  return 0;
424 }
425 
426 
427 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
word defaultFacesType
Definition: readKivaGrid.H:461
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
const word & name() const
Return name.
A class for handling file names.
Definition: fileName.H:69
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:466
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:110
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
label nInternalFaces() const
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
virtual bool write(const bool valid=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:877
label nFaces() const
Foam::word regionName
Output to file stream.
Definition: OFstream.H:82
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:312
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:461
List< FixedList< label, 8 > > cells() const
Return the cells for filling the block.
Definition: blockCreate.C:349
label nCells() const
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:309
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
static void noParallel()
Remove the parallel options.
Definition: argList.C:148
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
const Cmpt & z() const
Definition: VectorI.H:87
Xfer< T > xferCopy(const T &)
Construct by copying the contents of the arg.
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:198
const Cmpt & y() const
Definition: VectorI.H:81
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1011
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
dynamicFvMesh & mesh
A multi-block mesh generator.
Definition: blockMesh.H:61
void removeFiles(const fileName &instanceDir) const
Remove all files from mesh instance.
Definition: polyMesh.C:1212
const pointZoneMesh & pointZones() const
Return point zone mesh.
Definition: polyMesh.H:460
bool isDir(const fileName &, const bool followLink=true)
Does the name exist as a DIRECTORY in the file system?
Definition: POSIX.C:536
A class for handling words, derived from string.
Definition: word.H:59
fileName dictPath
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:472
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:95
const word dictName("particleTrackDict")
void addZones(const List< pointZone *> &pz, const List< faceZone *> &fz, const List< cellZone *> &cz)
Add mesh zones.
Definition: polyMesh.C:954
const vectorField & cellCentres() const
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
const Cmpt & x() const
Definition: VectorI.H:75
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
Creates a single block of cells from point coordinates, numbers of cells in each direction and an exp...
Definition: block.H:63
bool rmDir(const fileName &)
Remove a dirctory and its contents.
Definition: POSIX.C:1008
static const char nl
Definition: Ostream.H:262
A subset of mesh cells.
Definition: cellZone.H:61
word defaultFacesName
Definition: readKivaGrid.H:460
label patchi
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
A collection of cell labels.
Definition: cellSet.H:48
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:300
messageStream Info
label nPoints() const
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:85
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:126
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
virtual bool write(const bool valid=true) const
Write using setting from DB.
const word & zoneName() const
Return the (optional) zone name.
Foam::argList args(argc, argv)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
bool exists(const fileName &, const bool checkGzip=true, const bool followLink=true)
Does the name exist (as DIRECTORY or FILE) in the file system?
Definition: POSIX.C:517
void shrink()
Shrink the allocated table to approx. twice number of elements.
Definition: HashTable.C:500
Namespace for OpenFOAM.