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-2024 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 "argList.H"
52 #include "Time.H"
53 #include "systemDict.H"
54 #include "blockMesh.H"
55 #include "mergePatchPairs.H"
56 #include "polyTopoChange.H"
57 #include "emptyPolyPatch.H"
58 #include "cyclicPolyPatch.H"
59 #include "OFstream.H"
60 
61 using namespace Foam;
62 
63 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
65 int main(int argc, char *argv[])
66 {
68  Foam::argList::removeOption("noFunctionObjects");
69  #include "addDictOption.H"
71  (
72  "blockTopology",
73  "write block edges and centres as .obj files"
74  );
76  (
77  "noClean",
78  "keep the existing files in the polyMesh"
79  );
80 
82  (
83  "Block description\n"
84  "\n"
85  " For a given block, the correspondence between the ordering of\n"
86  " vertex labels and face labels is shown below.\n"
87  " For vertex numbering in the sequence 0 to 7 (block, centre):\n"
88  " faces 0 (f0) and 1 are left and right, respectively;\n"
89  " faces 2 and 3 are front and back; \n"
90  " and faces 4 and 5 are bottom and top::\n"
91  "\n"
92  " 7 ---- 6\n"
93  " f5 |\\ |\\ f3\n"
94  " | | 4 ---- 5 \\\n"
95  " | 3 |--- 2 | \\\n"
96  " | \\| \\| f2\n"
97  " f4 0 ---- 1\n"
98  "\n"
99  " Z f0 ----- f1\n"
100  " | Y\n"
101  " | /\n"
102  " O --- X\n"
103  );
104 
105  #include "addMeshOption.H"
106  #include "addRegionOption.H"
107  #include "setRootCase.H"
108  #include "setMeshPath.H"
109  #include "createTime.H"
110 
111  const word dictName("blockMeshDict");
112 
113 
114  // Check if the mesh is specified otherwise mesh the default mesh
115  if (!meshPath.empty())
116  {
117  Info<< nl << "Generating mesh " << meshPath << endl;
118  }
119 
121  word regionPath;
122 
123  // Check if the region is specified otherwise mesh the default region
125  {
126  Info<< nl << "Generating mesh for region " << regionName << endl;
127  regionPath = regionName;
128  }
129 
130  if (!args.optionFound("noClean"))
131  {
132  const fileName polyMeshPath
133  (
134  runTime.path()/runTime.constant()
135  /meshPath/regionPath/polyMesh::meshSubDir
136  );
137 
138  if (exists(polyMeshPath))
139  {
140  if (exists(polyMeshPath/dictName))
141  {
142  Info<< "Not deleting polyMesh directory " << nl
143  << " " << polyMeshPath << nl
144  << " because it contains " << dictName << endl;
145  }
146  else
147  {
148  Info<< "Deleting polyMesh directory" << nl
149  << " " << polyMeshPath << endl;
150  rmDir(polyMeshPath);
151  }
152  }
153  }
154 
155  typeIOobject<IOdictionary> meshDictIO
156  (
157  systemDictIO(dictName, args, runTime, regionName, meshPath)
158  );
159 
160  if (!meshDictIO.headerOk())
161  {
163  << "Cannot find file " << meshDictIO.relativeObjectPath()
164  << nl
165  << exit(FatalError);
166  }
167 
168  Info<< "Creating block mesh from\n "
169  << meshDictIO.relativeObjectPath() << endl;
170 
171  IOdictionary meshDict(meshDictIO);
172  blockMesh blocks(meshDict, runTime.constant()/meshPath, regionName);
173 
174 
175  if (args.optionFound("blockTopology"))
176  {
177  // Write mesh as edges.
178  {
179  fileName objMeshFile("blockTopology.obj");
180 
181  OFstream str(runTime.path()/objMeshFile);
182 
183  Info<< nl << "Dumping block structure as Lightwave obj format"
184  << " to " << objMeshFile << endl;
185 
186  blocks.writeTopology(str);
187  }
188 
189  // Write centres of blocks
190  {
191  fileName objCcFile("blockCentres.obj");
192 
193  OFstream str(runTime.path()/objCcFile);
194 
195  Info<< nl << "Dumping block centres as Lightwave obj format"
196  << " to " << objCcFile << endl;
197 
198  const polyMesh& topo = blocks.topology();
199 
200  const pointField& cellCentres = topo.cellCentres();
201 
202  forAll(cellCentres, celli)
203  {
204  // point cc = b.blockShape().centre(b.points());
205  const point& cc = cellCentres[celli];
206 
207  str << "v " << cc.x() << ' ' << cc.y() << ' ' << cc.z() << nl;
208  }
209  }
210 
211  Info<< nl << "end" << endl;
212 
213  return 0;
214  }
215 
216 
217  Info<< nl << "Creating polyMesh from blockMesh" << endl;
218 
219  word defaultFacesName = "defaultFaces";
220  word defaultFacesType = emptyPolyPatch::typeName;
221  polyMesh mesh
222  (
223  IOobject
224  (
225  regionName,
226  runTime.constant(),
227  meshPath,
228  runTime
229  ),
230  clone(blocks.points()), // could we reuse space?
231  blocks.cells(),
232  blocks.patches(),
233  blocks.patchNames(),
234  blocks.patchDicts(),
237  );
238 
239 
240  // Read in a list of dictionaries for the merge patch pairs
241  if (meshDict.found("mergePatchPairs"))
242  {
243  List<Pair<word>> patchPairNames
244  (
245  meshDict.lookup("mergePatchPairs")
246  );
247 
248  if (patchPairNames.size())
249  {
250  const word oldInstance = mesh.pointsInstance();
251  mergePatchPairs(mesh, patchPairNames, 1e-4);
252  mesh.setInstance(oldInstance);
253  }
254  }
255  else
256  {
257  Info<< nl << "No patch pairs to merge" << endl;
258  }
259 
260  label nZones = blocks.numZonedBlocks();
261 
262  if (nZones > 0)
263  {
264  Info<< nl << "Adding cell zones" << endl;
265 
266  // Map from zoneName to cellZone index
267  HashTable<label> zoneMap(nZones);
268 
269  // Cells per zone.
270  List<DynamicList<label>> zoneCells(nZones);
271 
272  // Running cell counter
273  label celli = 0;
274 
275  // Largest zone so far
276  label freeZoneI = 0;
277 
278  forAll(blocks, blockI)
279  {
280  const block& b = blocks[blockI];
281  const List<FixedList<label, 8>> blockCells = b.cells();
282  const word& zoneName = b.zoneName();
283 
284  if (zoneName.size())
285  {
286  HashTable<label>::const_iterator iter = zoneMap.find(zoneName);
287 
288  label zoneI;
289 
290  if (iter == zoneMap.end())
291  {
292  zoneI = freeZoneI++;
293 
294  Info<< " " << zoneI << '\t' << zoneName << endl;
295 
296  zoneMap.insert(zoneName, zoneI);
297  }
298  else
299  {
300  zoneI = iter();
301  }
302 
303  forAll(blockCells, i)
304  {
305  zoneCells[zoneI].append(celli++);
306  }
307  }
308  else
309  {
310  celli += blockCells.size();
311  }
312  }
313 
314 
315  List<cellZone*> cz(zoneMap.size());
316 
317  forAllConstIter(HashTable<label>, zoneMap, iter)
318  {
319  label zoneI = iter();
320 
321  cz[zoneI] = new cellZone
322  (
323  iter.key(),
324  zoneCells[zoneI].shrink(),
325  mesh.cellZones()
326  );
327  }
328 
329  mesh.pointZones().setSize(0);
330  mesh.faceZones().setSize(0);
331  mesh.cellZones().setSize(0);
333  }
334 
335 
336  // Detect any cyclic patches and force re-ordering of the faces
337  {
339  bool hasCyclic = false;
341  {
342  if (isA<cyclicPolyPatch>(patches[patchi]))
343  {
344  hasCyclic = true;
345  break;
346  }
347  }
348 
349  if (hasCyclic)
350  {
351  Info<< nl << "Detected cyclic patches; ordering boundary faces"
352  << endl;
353  const word oldInstance = mesh.instance();
354  polyTopoChange meshMod(mesh);
355  meshMod.changeMesh(mesh);
356  mesh.setInstance(oldInstance);
357  }
358  }
359 
360 
361  // Set the precision of the points data to 10
363 
364  Info<< nl << "Writing polyMesh" << endl;
365  mesh.removeFiles();
366  if (!mesh.write())
367  {
369  << "Failed writing polyMesh."
370  << exit(FatalError);
371  }
372 
373 
374  // Write summary
375  {
377 
378  Info<< "----------------" << nl
379  << "Mesh Information" << nl
380  << "----------------" << nl
381  << " " << "boundingBox: " << boundBox(mesh.points()) << nl
382  << " " << "nPoints: " << mesh.nPoints() << nl
383  << " " << "nCells: " << mesh.nCells() << nl
384  << " " << "nFaces: " << mesh.nFaces() << nl
385  << " " << "nInternalFaces: " << mesh.nInternalFaces() << nl;
386 
387  Info<< "----------------" << nl
388  << "Patches" << nl
389  << "----------------" << nl;
390 
392  {
393  const polyPatch& p = patches[patchi];
394 
395  Info<< " " << "patch " << patchi
396  << " (start: " << p.start()
397  << " size: " << p.size()
398  << ") name: " << p.name()
399  << nl;
400  }
401  }
402 
403  Info<< "\nEnd\n" << endl;
404 
405  return 0;
406 }
407 
408 
409 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
void shrink()
Shrink the allocated table to approx. twice number of elements.
Definition: HashTable.C:574
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:167
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
fileName & instance() const
Return the instance directory, constant, system, <time> etc.
Definition: IOobject.C:352
const word & name() const
Return name.
Definition: IOobject.H:307
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:473
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:87
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
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
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:255
static void noParallel()
Remove the parallel options.
Definition: argList.C:175
static void removeOption(const word &opt)
Remove option from validOptions and from optionUsage.
Definition: argList.C:168
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
Named list of cell indices representing a sub-set of the mesh.
Definition: cellZone.H:61
A class for handling file names.
Definition: fileName.H:82
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1785
Class to stitch mesh by merging patch-pairs.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const pointZoneList & pointZones() const
Return point zones.
Definition: polyMesh.H:437
const cellZoneList & cellZones() const
Return cell zones.
Definition: polyMesh.H:449
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:270
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:988
const faceZoneList & faceZones() const
Return face zones.
Definition: polyMesh.H:443
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
void removeFiles(const fileName &instanceDir) const
Remove all files from mesh instance.
Definition: polyMesh.C:1551
void addZones(const List< pointZone * > &pz, const List< faceZone * > &fz, const List< cellZone * > &cz)
Add mesh zones.
Definition: polyMesh.C:1145
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1331
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:273
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:91
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.
label nInternalFaces() const
const vectorField & cellCentres() const
label nCells() const
label nPoints() const
label nFaces() const
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:530
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
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
const doubleScalar e
Definition: doubleScalar.H:106
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:258
const word & regionName(const solver &region)
Definition: solver.H:218
IOobject systemDictIO(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion, const fileName &path=fileName::null)
Definition: systemDict.C:33
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
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:267
word defaultFacesName
Definition: readKivaGrid.H:454
word defaultFacesType
Definition: readKivaGrid.H:455
word meshPath
Definition: setMeshPath.H:1
Foam::argList args(argc, argv)
volScalarField & p