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-2016 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  "dict",
82  "file",
83  "specify alternative dictionary for the blockMesh description"
84  );
85 
86  #include "addRegionOption.H"
87  #include "setRootCase.H"
88  #include "createTime.H"
89 
90  const word dictName("blockMeshDict");
91 
93  word regionPath;
94 
95  // Check if the region is specified otherwise mesh the default region
96  if (args.optionReadIfPresent("region", regionName, polyMesh::defaultRegion))
97  {
98  Info<< nl << "Generating mesh for region " << regionName << endl;
99  regionPath = regionName;
100  }
101 
102  // Search for the appropriate blockMesh dictionary....
103 
105 
106  // Check if the dictionary is specified on the command-line
107  if (args.optionFound("dict"))
108  {
109  dictPath = args["dict"];
110 
111  dictPath =
112  (
113  isDir(dictPath)
114  ? dictPath/dictName
115  : dictPath
116  );
117  }
118  // Check if dictionary is present in the constant directory
119  else if
120  (
121  exists
122  (
123  runTime.path()/runTime.constant()
124  /regionPath/polyMesh::meshSubDir/dictName
125  )
126  )
127  {
128  dictPath =
129  runTime.constant()
130  /regionPath/polyMesh::meshSubDir/dictName;
131  }
132  // Otherwise assume the dictionary is present in the system directory
133  else
134  {
135  dictPath = runTime.system()/regionPath/dictName;
136  }
137 
138  IOobject meshDictIO
139  (
140  dictPath,
141  runTime,
144  false
145  );
146 
147  if (!meshDictIO.headerOk())
148  {
150  << meshDictIO.objectPath()
151  << nl
152  << exit(FatalError);
153  }
154 
155  Info<< "Creating block mesh from\n "
156  << meshDictIO.objectPath() << endl;
157 
158  blockMesh::verbose(true);
159 
160  IOdictionary meshDict(meshDictIO);
161  blockMesh blocks(meshDict, regionName);
162 
163 
164  if (args.optionFound("blockTopology"))
165  {
166  // Write mesh as edges.
167  {
168  fileName objMeshFile("blockTopology.obj");
169 
170  OFstream str(runTime.path()/objMeshFile);
171 
172  Info<< nl << "Dumping block structure as Lightwave obj format"
173  << " to " << objMeshFile << endl;
174 
175  blocks.writeTopology(str);
176  }
177 
178  // Write centres of blocks
179  {
180  fileName objCcFile("blockCentres.obj");
181 
182  OFstream str(runTime.path()/objCcFile);
183 
184  Info<< nl << "Dumping block centres as Lightwave obj format"
185  << " to " << objCcFile << endl;
186 
187  const polyMesh& topo = blocks.topology();
188 
189  const pointField& cellCentres = topo.cellCentres();
190 
191  forAll(cellCentres, celli)
192  {
193  //point cc = b.blockShape().centre(b.points());
194  const point& cc = cellCentres[celli];
195 
196  str << "v " << cc.x() << ' ' << cc.y() << ' ' << cc.z() << nl;
197  }
198  }
199 
200  Info<< nl << "end" << endl;
201 
202  return 0;
203  }
204 
205 
206  Info<< nl << "Creating polyMesh from blockMesh" << endl;
207 
208  word defaultFacesName = "defaultFaces";
209  word defaultFacesType = emptyPolyPatch::typeName;
210  polyMesh mesh
211  (
212  IOobject
213  (
214  regionName,
215  runTime.constant(),
216  runTime
217  ),
218  xferCopy(blocks.points()), // could we re-use space?
219  blocks.cells(),
220  blocks.patches(),
221  blocks.patchNames(),
222  blocks.patchDicts(),
225  );
226 
227 
228  // Read in a list of dictionaries for the merge patch pairs
229  if (meshDict.found("mergePatchPairs"))
230  {
231  List<Pair<word>> mergePatchPairs
232  (
233  meshDict.lookup("mergePatchPairs")
234  );
235 
236  #include "mergePatchPairs.H"
237  }
238  else
239  {
240  Info<< nl << "There are no merge patch pairs edges" << endl;
241  }
242 
243 
244  // Set any cellZones (note: cell labelling unaffected by above
245  // mergePatchPairs)
246 
247  label nZones = blocks.numZonedBlocks();
248 
249  if (nZones > 0)
250  {
251  Info<< nl << "Adding cell zones" << endl;
252 
253  // Map from zoneName to cellZone index
254  HashTable<label> zoneMap(nZones);
255 
256  // Cells per zone.
257  List<DynamicList<label>> zoneCells(nZones);
258 
259  // Running cell counter
260  label celli = 0;
261 
262  // Largest zone so far
263  label freeZoneI = 0;
264 
265  forAll(blocks, blockI)
266  {
267  const block& b = blocks[blockI];
268  const labelListList& blockCells = b.cells();
269  const word& zoneName = b.blockDef().zoneName();
270 
271  if (zoneName.size())
272  {
273  HashTable<label>::const_iterator iter = zoneMap.find(zoneName);
274 
275  label zoneI;
276 
277  if (iter == zoneMap.end())
278  {
279  zoneI = freeZoneI++;
280 
281  Info<< " " << zoneI << '\t' << zoneName << endl;
282 
283  zoneMap.insert(zoneName, zoneI);
284  }
285  else
286  {
287  zoneI = iter();
288  }
289 
290  forAll(blockCells, i)
291  {
292  zoneCells[zoneI].append(celli++);
293  }
294  }
295  else
296  {
297  celli += b.cells().size();
298  }
299  }
300 
301 
302  List<cellZone*> cz(zoneMap.size());
303 
304  Info<< nl << "Writing cell zones as cellSets" << endl;
305 
306  forAllConstIter(HashTable<label>, zoneMap, iter)
307  {
308  label zoneI = iter();
309 
310  cz[zoneI] = new cellZone
311  (
312  iter.key(),
313  zoneCells[zoneI].shrink(),
314  zoneI,
315  mesh.cellZones()
316  );
317 
318  // Write as cellSet for ease of processing
319  cellSet cset(mesh, iter.key(), zoneCells[zoneI].shrink());
320  cset.write();
321  }
322 
323  mesh.pointZones().setSize(0);
324  mesh.faceZones().setSize(0);
325  mesh.cellZones().setSize(0);
327  }
328 
329  // Set the precision of the points data to 10
331 
332  Info<< nl << "Writing polyMesh" << endl;
333  mesh.removeFiles();
334  if (!mesh.write())
335  {
337  << "Failed writing polyMesh."
338  << exit(FatalError);
339  }
340 
341 
342  //
343  // write some information
344  //
345  {
346  const polyPatchList& patches = mesh.boundaryMesh();
347 
348  Info<< "----------------" << nl
349  << "Mesh Information" << nl
350  << "----------------" << nl
351  << " " << "boundingBox: " << boundBox(mesh.points()) << nl
352  << " " << "nPoints: " << mesh.nPoints() << nl
353  << " " << "nCells: " << mesh.nCells() << nl
354  << " " << "nFaces: " << mesh.nFaces() << nl
355  << " " << "nInternalFaces: " << mesh.nInternalFaces() << nl;
356 
357  Info<< "----------------" << nl
358  << "Patches" << nl
359  << "----------------" << nl;
360 
361  forAll(patches, patchi)
362  {
363  const polyPatch& p = patches[patchi];
364 
365  Info<< " " << "patch " << patchi
366  << " (start: " << p.start()
367  << " size: " << p.size()
368  << ") name: " << p.name()
369  << nl;
370  }
371  }
372 
373  Info<< "\nEnd\n" << endl;
374 
375  return 0;
376 }
377 
378 
379 // ************************************************************************* //
word defaultFacesType
Definition: readKivaGrid.H:461
const Cmpt & z() const
Definition: VectorI.H:87
#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
A class for handling file names.
Definition: fileName.H:69
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:106
const Cmpt & x() const
Definition: VectorI.H:75
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::word regionName
Output to file stream.
Definition: OFstream.H:81
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:309
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:461
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:198
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:306
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
static void noParallel()
Remove the parallel options.
Definition: argList.C:146
const blockDescriptor & blockDef() const
Return the block definition.
Definition: blockI.H:41
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
Xfer< T > xferCopy(const T &)
Construct by copying the contents of the arg.
bool isDir(const fileName &)
Does the name exist as a DIRECTORY in the file system?
Definition: POSIX.C:486
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
label nCells() const
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:138
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:58
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:300
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
A class for handling words, derived from string.
Definition: word.H:59
fileName dictPath
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:93
const Cmpt & y() const
Definition: VectorI.H:81
const vectorField & cellCentres() const
const word & name() const
Return name.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
static void verbose(const bool on=true)
Enable/disable verbose information about the progress.
Definition: blockMesh.C:70
Creates a single block of cells from point coordinates, numbers of cells in each direction and an exp...
Definition: block.H:63
static const char nl
Definition: Ostream.H:262
const word & zoneName() const
Return the (optional) zone name.
const pointZoneMesh & pointZones() const
Return point zone mesh.
Definition: polyMesh.H:457
A subset of mesh cells.
Definition: cellZone.H:61
word defaultFacesName
Definition: readKivaGrid.H:460
label nFaces() const
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:469
bool exists(const fileName &, const bool checkGzip=true)
Does the name exist (as DIRECTORY or FILE) in the file system?
Definition: POSIX.C:480
label patchi
word dictName("noiseDict")
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:62
A collection of cell labels.
Definition: cellSet.H:48
virtual bool write() const
Write using setting from DB.
label nPoints() const
messageStream Info
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
void addZones(const List< pointZone * > &pz, const List< faceZone * > &fz, const List< cellZone * > &cz)
Add mesh zones.
Definition: polyMesh.C:922
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:83
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::argList args(argc, argv)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
label nInternalFaces() const
virtual bool write() const
Write mesh using IO settings from time.
Definition: fvMesh.C:870
void removeFiles(const fileName &instanceDir) const
Remove all files from mesh instance.
Definition: polyMesh.C:1170
void shrink()
Shrink the allocated table to approx. twice number of elements.
Definition: HashTable.C:496
Namespace for OpenFOAM.
const labelListList & cells() const
Return the cells for filling the block.
Definition: block.C:73