All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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
122  if (args.optionReadIfPresent("region", regionName, polyMesh::defaultRegion))
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  (
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  IOobject meshDictIO(systemDictIO(dictName, args, runTime, regionName));
153 
154  if (!meshDictIO.typeHeaderOk<IOdictionary>(true))
155  {
157  << "Cannot find file " << meshDictIO.localObjectPath()
158  << nl
159  << exit(FatalError);
160  }
161 
162  Info<< "Creating block mesh from\n "
163  << meshDictIO.localObjectPath() << endl;
164 
165  IOdictionary meshDict(meshDictIO);
166  blockMesh blocks(meshDict, regionName);
167 
168 
169  if (args.optionFound("blockTopology"))
170  {
171  // Write mesh as edges.
172  {
173  fileName objMeshFile("blockTopology.obj");
174 
175  OFstream str(runTime.path()/objMeshFile);
176 
177  Info<< nl << "Dumping block structure as Lightwave obj format"
178  << " to " << objMeshFile << endl;
179 
180  blocks.writeTopology(str);
181  }
182 
183  // Write centres of blocks
184  {
185  fileName objCcFile("blockCentres.obj");
186 
187  OFstream str(runTime.path()/objCcFile);
188 
189  Info<< nl << "Dumping block centres as Lightwave obj format"
190  << " to " << objCcFile << endl;
191 
192  const polyMesh& topo = blocks.topology();
193 
194  const pointField& cellCentres = topo.cellCentres();
195 
196  forAll(cellCentres, celli)
197  {
198  // point cc = b.blockShape().centre(b.points());
199  const point& cc = cellCentres[celli];
200 
201  str << "v " << cc.x() << ' ' << cc.y() << ' ' << cc.z() << nl;
202  }
203  }
204 
205  Info<< nl << "end" << endl;
206 
207  return 0;
208  }
209 
210 
211  Info<< nl << "Creating polyMesh from blockMesh" << endl;
212 
213  word defaultFacesName = "defaultFaces";
214  word defaultFacesType = emptyPolyPatch::typeName;
215  polyMesh mesh
216  (
217  IOobject
218  (
219  regionName,
220  runTime.constant(),
221  runTime
222  ),
223  clone(blocks.points()), // could we re-use space?
224  blocks.cells(),
225  blocks.patches(),
226  blocks.patchNames(),
227  blocks.patchDicts(),
230  );
231 
232 
233  // Read in a list of dictionaries for the merge patch pairs
234  if (meshDict.found("mergePatchPairs"))
235  {
236  List<Pair<word>> mergePatchPairs
237  (
238  meshDict.lookup("mergePatchPairs")
239  );
240 
241  #include "mergePatchPairs.H"
242  }
243  else
244  {
245  Info<< nl << "There are no merge patch pairs edges" << endl;
246  }
247 
248 
249  // Set any cellZones (note: cell labelling unaffected by above
250  // mergePatchPairs)
251 
252  label nZones = blocks.numZonedBlocks();
253 
254  if (nZones > 0)
255  {
256  Info<< nl << "Adding cell zones" << endl;
257 
258  // Map from zoneName to cellZone index
259  HashTable<label> zoneMap(nZones);
260 
261  // Cells per zone.
262  List<DynamicList<label>> zoneCells(nZones);
263 
264  // Running cell counter
265  label celli = 0;
266 
267  // Largest zone so far
268  label freeZoneI = 0;
269 
270  forAll(blocks, blockI)
271  {
272  const block& b = blocks[blockI];
273  const List<FixedList<label, 8>> blockCells = b.cells();
274  const word& zoneName = b.zoneName();
275 
276  if (zoneName.size())
277  {
278  HashTable<label>::const_iterator iter = zoneMap.find(zoneName);
279 
280  label zoneI;
281 
282  if (iter == zoneMap.end())
283  {
284  zoneI = freeZoneI++;
285 
286  Info<< " " << zoneI << '\t' << zoneName << endl;
287 
288  zoneMap.insert(zoneName, zoneI);
289  }
290  else
291  {
292  zoneI = iter();
293  }
294 
295  forAll(blockCells, i)
296  {
297  zoneCells[zoneI].append(celli++);
298  }
299  }
300  else
301  {
302  celli += b.cells().size();
303  }
304  }
305 
306 
307  List<cellZone*> cz(zoneMap.size());
308 
309  forAllConstIter(HashTable<label>, zoneMap, iter)
310  {
311  label zoneI = iter();
312 
313  cz[zoneI] = new cellZone
314  (
315  iter.key(),
316  zoneCells[zoneI].shrink(),
317  zoneI,
318  mesh.cellZones()
319  );
320  }
321 
322  mesh.pointZones().setSize(0);
323  mesh.faceZones().setSize(0);
324  mesh.cellZones().setSize(0);
326  }
327 
328 
329  // Detect any cyclic patches and force re-ordering of the faces
330  {
331  const polyPatchList& patches = mesh.boundaryMesh();
332  bool hasCyclic = false;
333  forAll(patches, patchi)
334  {
335  if (isA<cyclicPolyPatch>(patches[patchi]))
336  {
337  hasCyclic = true;
338  break;
339  }
340  }
341 
342  if (hasCyclic)
343  {
344  Info<< nl << "Detected cyclic patches; ordering boundary faces"
345  << endl;
346  const word oldInstance = mesh.instance();
347  polyTopoChange meshMod(mesh);
348  meshMod.changeMesh(mesh, false);
349  mesh.setInstance(oldInstance);
350  }
351  }
352 
353 
354  // Set the precision of the points data to 10
356 
357  Info<< nl << "Writing polyMesh" << endl;
358  mesh.removeFiles();
359  if (!mesh.write())
360  {
362  << "Failed writing polyMesh."
363  << exit(FatalError);
364  }
365 
366 
367  // Write summary
368  {
369  const polyPatchList& patches = mesh.boundaryMesh();
370 
371  Info<< "----------------" << nl
372  << "Mesh Information" << nl
373  << "----------------" << nl
374  << " " << "boundingBox: " << boundBox(mesh.points()) << nl
375  << " " << "nPoints: " << mesh.nPoints() << nl
376  << " " << "nCells: " << mesh.nCells() << nl
377  << " " << "nFaces: " << mesh.nFaces() << nl
378  << " " << "nInternalFaces: " << mesh.nInternalFaces() << nl;
379 
380  Info<< "----------------" << nl
381  << "Patches" << nl
382  << "----------------" << nl;
383 
384  forAll(patches, patchi)
385  {
386  const polyPatch& p = patches[patchi];
387 
388  Info<< " " << "patch " << patchi
389  << " (start: " << p.start()
390  << " size: " << p.size()
391  << ") name: " << p.name()
392  << nl;
393  }
394  }
395 
396  Info<< "\nEnd\n" << endl;
397 
398  return 0;
399 }
400 
401 
402 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
word defaultFacesType
Definition: readKivaGrid.H:461
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
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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:79
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:112
const meshCellZones & cellZones() const
Return cell zones.
Definition: polyMesh.H:482
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
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:59
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:458
List< FixedList< label, 8 > > cells() const
Return the cells for filling the block.
Definition: blockCreate.C:349
label nCells() const
engineTime & runTime
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:309
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
static void noParallel()
Remove the parallel options.
Definition: argList.C:175
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
const Cmpt & z() const
Definition: VectorI.H:87
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1093
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:204
const Cmpt & y() const
Definition: VectorI.H:81
const meshPointZones & pointZones() const
Return point zones.
Definition: polyMesh.H:470
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1131
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...
T clone(const T &t)
Definition: List.H:54
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:1424
A class for handling words, derived from string.
Definition: word.H:59
const word & constant() const
Return constant name.
Definition: TimePaths.H:123
void addZones(const List< pointZone *> &pz, const List< faceZone *> &fz, const List< cellZone *> &cz)
Add mesh zones.
Definition: polyMesh.C:953
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 directory and its contents.
Definition: POSIX.C:1047
static const char nl
Definition: Ostream.H:260
A subset of mesh cells.
Definition: cellZone.H:61
word defaultFacesName
Definition: readKivaGrid.H:460
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:32
const fileName & instance() const
Definition: IOobject.H:390
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:70
const meshFaceZones & faceZones() const
Return face zones.
Definition: polyMesh.H:476
Direct mesh changes based on v1.3 polyTopoChange syntax.
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:309
messageStream Info
label nPoints() const
IOobject systemDictIO(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion)
Definition: systemDict.C:33
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:118
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
const word & zoneName() const
Return the (optional) zone name.
Foam::argList args(argc, argv)
fileName path() const
Return path.
Definition: TimePaths.H:139
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
const word dictName("noiseDict")
void shrink()
Shrink the allocated table to approx. twice number of elements.
Definition: HashTable.C:500
Namespace for OpenFOAM.