blockMesh.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-2021 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 #include "systemDict.H"
55 
56 #include "blockMesh.H"
57 #include "attachPolyTopoChanger.H"
58 #include "polyTopoChange.H"
59 #include "emptyPolyPatch.H"
60 #include "cyclicPolyPatch.H"
61 
62 #include "argList.H"
63 #include "OSspecific.H"
64 #include "OFstream.H"
65 
66 #include "Pair.H"
67 #include "slidingInterface.H"
68 
69 using namespace Foam;
70 
71 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
72 
73 int main(int argc, char *argv[])
74 {
76  #include "addDictOption.H"
78  (
79  "blockTopology",
80  "write block edges and centres as .obj files"
81  );
83  (
84  "noClean",
85  "keep the existing files in the polyMesh"
86  );
87 
89  (
90  "Block description\n"
91  "\n"
92  " For a given block, the correspondence between the ordering of\n"
93  " vertex labels and face labels is shown below.\n"
94  " For vertex numbering in the sequence 0 to 7 (block, centre):\n"
95  " faces 0 (f0) and 1 are left and right, respectively;\n"
96  " faces 2 and 3 are front and back; \n"
97  " and faces 4 and 5 are bottom and top::\n"
98  "\n"
99  " 7 ---- 6\n"
100  " f5 |\\ |\\ f3\n"
101  " | | 4 ---- 5 \\\n"
102  " | 3 |--- 2 | \\\n"
103  " | \\| \\| f2\n"
104  " f4 0 ---- 1\n"
105  "\n"
106  " Z f0 ----- f1\n"
107  " | Y\n"
108  " | /\n"
109  " O --- X\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
123  {
124  Info<< nl << "Generating mesh for region " << regionName << endl;
125  regionPath = regionName;
126  }
127 
128  if (!args.optionFound("noClean"))
129  {
130  fileName polyMeshPath
131  (
132  runTime.path()/runTime.constant()/regionPath/polyMesh::meshSubDir
133  );
134 
135  if (exists(polyMeshPath))
136  {
137  if (exists(polyMeshPath/dictName))
138  {
139  Info<< "Not deleting polyMesh directory " << nl
140  << " " << polyMeshPath << nl
141  << " because it contains " << dictName << endl;
142  }
143  else
144  {
145  Info<< "Deleting polyMesh directory" << nl
146  << " " << polyMeshPath << endl;
147  rmDir(polyMeshPath);
148  }
149  }
150  }
151 
152  typeIOobject<IOdictionary> meshDictIO
153  (
154  systemDictIO(dictName, args, runTime, regionName)
155  );
156 
157  if (!meshDictIO.headerOk())
158  {
160  << "Cannot find file " << meshDictIO.relativeObjectPath()
161  << nl
162  << exit(FatalError);
163  }
164 
165  Info<< "Creating block mesh from\n "
166  << meshDictIO.relativeObjectPath() << endl;
167 
168  IOdictionary meshDict(meshDictIO);
169  blockMesh blocks(meshDict, regionName);
170 
171 
172  if (args.optionFound("blockTopology"))
173  {
174  // Write mesh as edges.
175  {
176  fileName objMeshFile("blockTopology.obj");
177 
178  OFstream str(runTime.path()/objMeshFile);
179 
180  Info<< nl << "Dumping block structure as Lightwave obj format"
181  << " to " << objMeshFile << endl;
182 
183  blocks.writeTopology(str);
184  }
185 
186  // Write centres of blocks
187  {
188  fileName objCcFile("blockCentres.obj");
189 
190  OFstream str(runTime.path()/objCcFile);
191 
192  Info<< nl << "Dumping block centres as Lightwave obj format"
193  << " to " << objCcFile << endl;
194 
195  const polyMesh& topo = blocks.topology();
196 
197  const pointField& cellCentres = topo.cellCentres();
198 
199  forAll(cellCentres, celli)
200  {
201  // point cc = b.blockShape().centre(b.points());
202  const point& cc = cellCentres[celli];
203 
204  str << "v " << cc.x() << ' ' << cc.y() << ' ' << cc.z() << nl;
205  }
206  }
207 
208  Info<< nl << "end" << endl;
209 
210  return 0;
211  }
212 
213 
214  Info<< nl << "Creating polyMesh from blockMesh" << endl;
215 
216  word defaultFacesName = "defaultFaces";
217  word defaultFacesType = emptyPolyPatch::typeName;
218  polyMesh mesh
219  (
220  IOobject
221  (
222  regionName,
223  runTime.constant(),
224  runTime
225  ),
226  clone(blocks.points()), // could we re-use space?
227  blocks.cells(),
228  blocks.patches(),
229  blocks.patchNames(),
230  blocks.patchDicts(),
233  );
234 
235 
236  // Read in a list of dictionaries for the merge patch pairs
237  if (meshDict.found("mergePatchPairs"))
238  {
239  List<Pair<word>> mergePatchPairs
240  (
241  meshDict.lookup("mergePatchPairs")
242  );
243 
244  #include "mergePatchPairs.H"
245  }
246  else
247  {
248  Info<< nl << "There are no merge patch pairs edges" << endl;
249  }
250 
251 
252  // Set any cellZones (note: cell labelling unaffected by above
253  // mergePatchPairs)
254 
255  label nZones = blocks.numZonedBlocks();
256 
257  if (nZones > 0)
258  {
259  Info<< nl << "Adding cell zones" << endl;
260 
261  // Map from zoneName to cellZone index
262  HashTable<label> zoneMap(nZones);
263 
264  // Cells per zone.
265  List<DynamicList<label>> zoneCells(nZones);
266 
267  // Running cell counter
268  label celli = 0;
269 
270  // Largest zone so far
271  label freeZoneI = 0;
272 
273  forAll(blocks, blockI)
274  {
275  const block& b = blocks[blockI];
276  const List<FixedList<label, 8>> blockCells = b.cells();
277  const word& zoneName = b.zoneName();
278 
279  if (zoneName.size())
280  {
281  HashTable<label>::const_iterator iter = zoneMap.find(zoneName);
282 
283  label zoneI;
284 
285  if (iter == zoneMap.end())
286  {
287  zoneI = freeZoneI++;
288 
289  Info<< " " << zoneI << '\t' << zoneName << endl;
290 
291  zoneMap.insert(zoneName, zoneI);
292  }
293  else
294  {
295  zoneI = iter();
296  }
297 
298  forAll(blockCells, i)
299  {
300  zoneCells[zoneI].append(celli++);
301  }
302  }
303  else
304  {
305  celli += blockCells.size();
306  }
307  }
308 
309 
310  List<cellZone*> cz(zoneMap.size());
311 
312  forAllConstIter(HashTable<label>, zoneMap, iter)
313  {
314  label zoneI = iter();
315 
316  cz[zoneI] = new cellZone
317  (
318  iter.key(),
319  zoneCells[zoneI].shrink(),
320  zoneI,
321  mesh.cellZones()
322  );
323  }
324 
325  mesh.pointZones().setSize(0);
326  mesh.faceZones().setSize(0);
327  mesh.cellZones().setSize(0);
328  mesh.addZones(List<pointZone*>(0), List<faceZone*>(0), cz);
329  }
330 
331 
332  // Detect any cyclic patches and force re-ordering of the faces
333  {
334  const polyPatchList& patches = mesh.boundaryMesh();
335  bool hasCyclic = false;
337  {
338  if (isA<cyclicPolyPatch>(patches[patchi]))
339  {
340  hasCyclic = true;
341  break;
342  }
343  }
344 
345  if (hasCyclic)
346  {
347  Info<< nl << "Detected cyclic patches; ordering boundary faces"
348  << endl;
349  const word oldInstance = mesh.instance();
350  polyTopoChange meshMod(mesh);
351  meshMod.changeMesh(mesh, false);
352  mesh.setInstance(oldInstance);
353  }
354  }
355 
356 
357  // Set the precision of the points data to 10
359 
360  Info<< nl << "Writing polyMesh" << endl;
361  mesh.removeFiles();
362  if (!mesh.write())
363  {
365  << "Failed writing polyMesh."
366  << exit(FatalError);
367  }
368 
369 
370  // Write summary
371  {
372  const polyPatchList& patches = mesh.boundaryMesh();
373 
374  Info<< "----------------" << nl
375  << "Mesh Information" << nl
376  << "----------------" << nl
377  << " " << "boundingBox: " << boundBox(mesh.points()) << nl
378  << " " << "nPoints: " << mesh.nPoints() << nl
379  << " " << "nCells: " << mesh.nCells() << nl
380  << " " << "nFaces: " << mesh.nFaces() << nl
381  << " " << "nInternalFaces: " << mesh.nInternalFaces() << nl;
382 
383  Info<< "----------------" << nl
384  << "Patches" << nl
385  << "----------------" << nl;
386 
388  {
389  const polyPatch& p = patches[patchi];
390 
391  Info<< " " << "patch " << patchi
392  << " (start: " << p.start()
393  << " size: " << p.size()
394  << ") name: " << p.name()
395  << nl;
396  }
397  }
398 
399  Info<< "\nEnd\n" << endl;
400 
401  return 0;
402 }
403 
404 
405 // ************************************************************************* //
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
void shrink()
Shrink the allocated table to approx. twice number of elements.
Definition: HashTable.C:500
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:310
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:458
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Output to file stream.
Definition: OFstream.H:86
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
const Cmpt & z() const
Definition: VectorI.H:87
const Cmpt & y() const
Definition: VectorI.H:81
const Cmpt & x() const
Definition: VectorI.H:75
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:204
static void noParallel()
Remove the parallel options.
Definition: argList.C:175
A multi-block mesh generator.
Definition: blockMesh.H:65
Creates a single block of cells from point coordinates, numbers of cells in each direction and an exp...
Definition: block.H:66
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:59
A subset of mesh cells.
Definition: cellZone.H:64
A class for handling file names.
Definition: fileName.H:82
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:269
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:266
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
Direct mesh changes based on v1.3 polyTopoChange syntax.
const vectorField & cellCentres() const
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:531
A class for handling words, derived from string.
Definition: word.H:62
Foam::word regionName
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
int main(int argc, char *argv[])
Definition: financialFoam.C:44
label patchi
const fvPatchList & patches
volScalarField & b
Definition: createFields.H:27
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
bool exists(const fileName &, const bool checkVariants=true, const bool followLink=true)
Does the name exist (as directory or file) in the file system?
Definition: POSIX.C:520
IOobject systemDictIO(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion)
Definition: systemDict.C:33
messageStream Info
bool rmDir(const fileName &)
Remove a directory and its contents.
Definition: POSIX.C:1047
T clone(const T &t)
Definition: List.H:55
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
error FatalError
static const char nl
Definition: Ostream.H:260
word defaultFacesName
Definition: readKivaGrid.H:454
word defaultFacesType
Definition: readKivaGrid.H:455
Foam::argList args(argc, argv)
volScalarField & p