initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / applications / utilities / mesh / generation / blockMesh / blockMeshApp.C
blobe73a5262b9b33702b83ecccb9a730254b56a1776
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Application
26     blockMesh
28 Description
29     A multi-block mesh generator.
31     Uses the block mesh description found in
32     @a constant/polyMesh/blockMeshDict
33     (or @a constant/\<region\>/polyMesh/blockMeshDict).
35 Usage
37     - blockMesh [OPTION]
39     @param -blockTopology \n
40     Write the topology as a set of edges in OBJ format.
42     @param -region \<name\> \n
43     Specify an alternative mesh region.
45     @param -dict \<dictionary\> \n
46     Specify an alternative dictionary for the block mesh description.
48 \*---------------------------------------------------------------------------*/
50 #include "Time.H"
51 #include "IOdictionary.H"
52 #include "IOPtrList.H"
54 #include "blockMesh.H"
55 #include "attachPolyTopoChanger.H"
56 #include "preservePatchTypes.H"
57 #include "emptyPolyPatch.H"
58 #include "cellSet.H"
60 #include "argList.H"
61 #include "OSspecific.H"
62 #include "OFstream.H"
64 #include "Pair.H"
65 #include "slidingInterface.H"
67 using namespace Foam;
69 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
70 // Main program:
72 int main(int argc, char *argv[])
74     argList::noParallel();
76 #   include "addOptions.H"
77 #   include "setRootCase.H"
78 #   include "createTime.H"
79 #   include "checkOptions.H"
81     word regionName;
82     fileName polyMeshDir;
83     word dictName("blockMeshDict");
84     fileName dictPath(runTime.constant());
86     if (args.options().found("region"))
87     {
88         // constant/<region>/polyMesh/blockMeshDict
89         regionName  = args.options()["region"];
90         polyMeshDir = regionName/polyMesh::meshSubDir;
92         Info<< nl << "Generating mesh for region " << regionName << endl;
93     }
94     else
95     {
96         // constant/polyMesh/blockMeshDict
97         regionName  = polyMesh::defaultRegion;
98         polyMeshDir = polyMesh::meshSubDir;
99     }
101     fileName dictLocal = polyMeshDir;
103     if (args.options().found("dict"))
104     {
105         wordList elems(fileName(args.options()["dict"]).components());
106         dictName = elems[elems.size()-1];
107         dictPath = elems[0];
108         dictLocal = "";
110         if (elems.size() == 1)
111         {
112             dictPath = ".";
113         }
114         else if (elems.size() > 2)
115         {
116             dictLocal = fileName(SubList<word>(elems, elems.size()-2, 1));
117         }
118     }
121     IOobject meshDictIo
122     (
123         dictName,
124         dictPath,
125         dictLocal,
126         runTime,
127         IOobject::MUST_READ,
128         IOobject::NO_WRITE,
129         false
130     );
132     if (!meshDictIo.headerOk())
133     {
134         FatalErrorIn(args.executable())
135             << "Cannot open mesh description file\n    "
136             << meshDictIo.objectPath()
137             << nl
138             << exit(FatalError);
139     }
141     Info<< nl << "Creating block mesh from\n    "
142         << meshDictIo.objectPath() << endl;
144     IOdictionary meshDict(meshDictIo);
146     blockMesh blocks(meshDict);
148     if (writeTopo)
149     {
150         // Write mesh as edges.
151         {
152             fileName objMeshFile("blockTopology.obj");
154             OFstream str(runTime.path()/objMeshFile);
156             Info<< nl << "Dumping block structure as Lightwave obj format"
157                 << " to " << objMeshFile << endl;
159             blocks.writeTopology(str);
160         }
162         // Write centres of blocks
163         {
164             fileName objCcFile("blockCentres.obj");
166             OFstream str(runTime.path()/objCcFile);
168             Info<< nl << "Dumping block centres as Lightwave obj format"
169                 << " to " << objCcFile << endl;
171             const polyMesh& topo = blocks.topology();
173             const pointField& cellCentres = topo.cellCentres();
175             forAll(cellCentres, cellI)
176             {
177                 //point cc = b.blockShape().centre(b.points());
178                 const point& cc = cellCentres[cellI];
180                 str << "v " << cc.x() << ' ' << cc.y() << ' ' << cc.z() << nl;
181             }
182         }
184         Info<< nl << "end" << endl;
186         return 0;
187     }
191     Info<< nl << "Creating mesh from block mesh" << endl;
193     wordList patchNames = blocks.patchNames();
194     wordList patchTypes = blocks.patchTypes();
195     word defaultFacesName = "defaultFaces";
196     word defaultFacesType = emptyPolyPatch::typeName;
197     wordList patchPhysicalTypes = blocks.patchPhysicalTypes();
199     preservePatchTypes
200     (
201         runTime,
202         runTime.constant(),
203         polyMeshDir,
204         patchNames,
205         patchTypes,
206         defaultFacesName,
207         defaultFacesType,
208         patchPhysicalTypes
209     );
211     polyMesh mesh
212     (
213         IOobject
214         (
215             regionName,
216             runTime.constant(),
217             runTime
218         ),
219         blocks.points(),
220         blocks.cells(),
221         blocks.patches(),
222         patchNames,
223         patchTypes,
224         defaultFacesName,
225         defaultFacesType,
226         patchPhysicalTypes
227     );
230     // Read in a list of dictionaries for the merge patch pairs
231     if (meshDict.found("mergePatchPairs"))
232     {
233         List<Pair<word> > mergePatchPairs
234         (
235             meshDict.lookup("mergePatchPairs")
236         );
238         if (mergePatchPairs.size())
239         {
240             FatalErrorIn(args.executable())
241                 << "mergePatchPairs not currently supported."
242                 << exit(FatalError);
243         }
245         //#include "mergePatchPairs.H"
246     }
247     else
248     {
249         Info<< nl << "There are no merge patch pairs edges" << endl;
250     }
253     // Set any cellZones (note: cell labelling unaffected by above
254     // mergePatchPairs)
256     label nZones = blocks.numZonedBlocks();
258     if (nZones > 0)
259     {
260         Info<< nl << "Adding cell zones" << endl;
262         // Map from zoneName to cellZone index
263         HashTable<label> zoneMap(nZones);
265         // Cells per zone.
266         List<DynamicList<label> > zoneCells(nZones);
268         // Running cell counter
269         label cellI = 0;
271         // Largest zone so far
272         label freeZoneI = 0;
274         forAll(blocks, blockI)
275         {
276             const block& b = blocks[blockI];
277             const labelListList& blockCells = b.cells();
278             const word& zoneName = b.blockDef().zoneName();
280             if (zoneName.size() > 0)
281             {
282                 HashTable<label>::const_iterator iter = zoneMap.find(zoneName);
284                 label zoneI;
286                 if (iter == zoneMap.end())
287                 {
288                     zoneI = freeZoneI++;
290                     Info<< "    " << zoneI << '\t' << zoneName << endl;
292                     zoneMap.insert(zoneName, zoneI);
293                 }
294                 else
295                 {
296                     zoneI = iter();
297                 }
299                 forAll(blockCells, i)
300                 {
301                     zoneCells[zoneI].append(cellI++);
302                 }
303             }
304             else
305             {
306                 cellI += b.cells().size();
307             }
308         }
311         List<cellZone*> cz(zoneMap.size());
313         Info<< nl << "Writing cell zones as cellSets" << endl;
315         forAllConstIter(HashTable<label>, zoneMap, iter)
316         {
317             label zoneI = iter();
319             cz[zoneI]= new cellZone
320             (
321                 iter.key(),
322                 zoneCells[zoneI].shrink(),
323                 zoneI,
324                 mesh.cellZones()
325             );
327             // Write as cellSet for ease of processing
328             cellSet cset(mesh, iter.key(), zoneCells[zoneI].shrink());
329             cset.write();
330         }
332         mesh.pointZones().setSize(0);
333         mesh.faceZones().setSize(0);
334         mesh.cellZones().setSize(0);
335         mesh.addZones(List<pointZone*>(0), List<faceZone*>(0), cz);
336     }
338     // Set the precision of the points data to 10
339     IOstream::defaultPrecision(10);
341     Info << nl << "Writing polyMesh" << endl;
342     mesh.removeFiles(mesh.instance());
343     if (!mesh.write())
344     {
345         FatalErrorIn(args.executable())
346             << "Failed writing polyMesh."
347             << exit(FatalError);
348     }
350     Info<< nl << "end" << endl;
351     return 0;
355 // ************************************************************************* //