Forward compatibility: flex
[foam-extend-3.2.git] / applications / utilities / mesh / generation / blockMesh / blockMeshApp.C
blob87d65a94f9a82c1ff4c7e4d374860d26055aecef
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 Application
25     blockMesh
27 Description
28     A multi-block mesh generator.
30     Uses the block mesh description found in
31     \a constant/polyMesh/blockMeshDict
32     (or \a constant/\<region\>/polyMesh/blockMeshDict).
34 Usage
36     - blockMesh [OPTION]
38     \param -blockTopology \n
39     Write the topology as a set of edges in OBJ format.
41     \param -region \<name\> \n
42     Specify an alternative mesh region.
44     \param -dict \<filename\> \n
45     Specify alternative dictionary for the block mesh description.
47 \*---------------------------------------------------------------------------*/
49 #include "objectRegistry.H"
50 #include "foamTime.H"
51 #include "IOdictionary.H"
52 #include "IOPtrList.H"
54 #include "blockMesh.H"
55 #include "preservePatchTypes.H"
56 #include "emptyPolyPatch.H"
57 #include "cellSet.H"
59 #include "argList.H"
60 #include "OSspecific.H"
61 #include "OFstream.H"
63 #include "Pair.H"
64 #include "mapPolyMesh.H"
65 #include "polyTopoChanger.H"
66 #include "slidingInterface.H"
68 using namespace Foam;
70 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71 // Main program:
73 int main(int argc, char *argv[])
75     argList::noParallel();
76     argList::validOptions.insert("blockTopology", "");
77     argList::validOptions.insert("dict", "dictionary");
78 #   include "addRegionOption.H"
79 #   include "setRootCase.H"
80 #   include "createTime.H"
82     const word dictName("blockMeshDict");
84     word regionName;
85     fileName polyMeshDir;
87     if (args.optionFound("region"))
88     {
89         // constant/<region>/polyMesh/blockMeshDict
90         regionName  = args.option("region");
91         polyMeshDir = regionName/polyMesh::meshSubDir;
93         Info<< nl << "Generating mesh for region " << regionName << endl;
94     }
95     else
96     {
97         // constant/polyMesh/blockMeshDict
98         regionName  = polyMesh::defaultRegion;
99         polyMeshDir = polyMesh::meshSubDir;
100     }
102     autoPtr<IOobject> meshDictIoPtr;
104     if (args.optionFound("dict"))
105     {
106         fileName dictPath(args.option("dict"));
108         meshDictIoPtr.set
109         (
110             new IOobject
111             (
112                 (
113                     isDir(dictPath)
114                   ? dictPath/dictName
115                   : dictPath
116                 ),
117                 runTime,
118                 IOobject::MUST_READ,
119                 IOobject::NO_WRITE,
120                 false
121             )
122         );
123     }
124     else
125     {
126         meshDictIoPtr.set
127         (
128             new IOobject
129             (
130                 dictName,
131                 runTime.constant(),
132                 polyMeshDir,
133                 runTime,
134                 IOobject::MUST_READ,
135                 IOobject::NO_WRITE,
136                 false
137             )
138         );
139     }
141     if (!meshDictIoPtr->headerOk())
142     {
143         FatalErrorIn(args.executable())
144             << "Cannot open mesh description file\n    "
145             << meshDictIoPtr->objectPath()
146             << nl
147             << exit(FatalError);
148     }
150     Info<< nl << "Creating block mesh from\n    "
151         << meshDictIoPtr->objectPath() << nl << endl;
153     blockMesh::verbose(true);
155     IOdictionary meshDict(meshDictIoPtr());
156     blockMesh blocks(meshDict, regionName);
159     if (args.optionFound("blockTopology"))
160     {
161         // Write mesh as edges.
162         {
163             fileName objMeshFile("blockTopology.obj");
165             OFstream str(runTime.path()/objMeshFile);
167             Info<< nl << "Dumping block structure as Lightwave obj format"
168                 << " to " << objMeshFile << endl;
170             blocks.writeTopology(str);
171         }
173         // Write centres of blocks
174         {
175             fileName objCcFile("blockCentres.obj");
177             OFstream str(runTime.path()/objCcFile);
179             Info<< nl << "Dumping block centres as Lightwave obj format"
180                 << " to " << objCcFile << endl;
182             const polyMesh& topo = blocks.topology();
184             const pointField& cellCentres = topo.cellCentres();
186             forAll(cellCentres, cellI)
187             {
188                 //point cc = b.blockShape().centre(b.points());
189                 const point& cc = cellCentres[cellI];
191                 str << "v " << cc.x() << ' ' << cc.y() << ' ' << cc.z() << nl;
192             }
193         }
195         Info<< nl << "end" << endl;
197         return 0;
198     }
201     Info<< nl << "Creating polyMesh from blockMesh" << endl;
203     word defaultFacesName = "defaultFaces";
204     word defaultFacesType = emptyPolyPatch::typeName;
205     polyMesh mesh
206     (
207         IOobject
208         (
209             regionName,
210             runTime.constant(),
211             runTime
212         ),
213         xferCopy(blocks.points()),           // could we re-use space?
214         blocks.cells(),
215         blocks.patches(),
216         blocks.patchNames(),
217         blocks.patchDicts(),
218         defaultFacesName,
219         defaultFacesType
220     );
223     // Read in a list of dictionaries for the merge patch pairs
224     if (meshDict.found("mergePatchPairs"))
225     {
226         List<Pair<word> > mergePatchPairs
227         (
228             meshDict.lookup("mergePatchPairs")
229         );
231 #       include "mergePatchPairs.H"
232     }
233     else
234     {
235         Info<< nl << "There are no merge patch pairs edges" << endl;
236     }
239     // Set any cellZones (note: cell labelling unaffected by above
240     // mergePatchPairs)
242     label nZones = blocks.numZonedBlocks();
244     if (nZones > 0)
245     {
246         Info<< nl << "Adding cell zones" << endl;
248         // Map from zoneName to cellZone index
249         HashTable<label> zoneMap(nZones);
251         // Cells per zone.
252         List<DynamicList<label> > zoneCells(nZones);
254         // Running cell counter
255         label cellI = 0;
257         // Largest zone so far
258         label freeZoneI = 0;
260         forAll(blocks, blockI)
261         {
262             const block& b = blocks[blockI];
263             const labelListList& blockCells = b.cells();
264             const word& zoneName = b.blockDef().zoneName();
266             if (zoneName.size())
267             {
268                 HashTable<label>::const_iterator iter = zoneMap.find(zoneName);
270                 label zoneI;
272                 if (iter == zoneMap.end())
273                 {
274                     zoneI = freeZoneI++;
276                     Info<< "    " << zoneI << '\t' << zoneName << endl;
278                     zoneMap.insert(zoneName, zoneI);
279                 }
280                 else
281                 {
282                     zoneI = iter();
283                 }
285                 forAll(blockCells, i)
286                 {
287                     zoneCells[zoneI].append(cellI++);
288                 }
289             }
290             else
291             {
292                 cellI += b.cells().size();
293             }
294         }
297         List<cellZone*> cz(zoneMap.size());
299         Info<< nl << "Writing cell zones as cellSets" << endl;
301         forAllConstIter(HashTable<label>, zoneMap, iter)
302         {
303             label zoneI = iter();
305             cz[zoneI] = new cellZone
306             (
307                 iter.key(),
308                 zoneCells[zoneI].shrink(),
309                 zoneI,
310                 mesh.cellZones()
311             );
313             // Write as cellSet for ease of processing
314             cellSet cset
315             (
316                 mesh,
317                 iter.key(),
318                 labelHashSet(zoneCells[zoneI].shrink())
319             );
320             cset.write();
321         }
323         mesh.pointZones().setSize(0);
324         mesh.faceZones().setSize(0);
325         mesh.cellZones().setSize(0);
326         mesh.addZones(List<pointZone*>(0), List<faceZone*>(0), cz);
327     }
329     // Set the precision of the points data to 10
330     IOstream::defaultPrecision(max(10u, IOstream::defaultPrecision()));
332     Info<< nl << "Writing polyMesh" << endl;
333     mesh.removeFiles();
334     if (!mesh.write())
335     {
336         FatalErrorIn(args.executable())
337             << "Failed writing polyMesh."
338             << exit(FatalError);
339     }
342     //
343     // write some information
344     //
345     {
346         const polyPatchList& patches = mesh.boundaryMesh();
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;
357         Info<< "----------------" << nl
358             << "Patches" << nl
359             << "----------------" << nl;
361         forAll(patches, patchI)
362         {
363             const polyPatch& p = patches[patchI];
365             Info<< "  " << "patch " << patchI
366                 << " (start: " << p.start()
367                 << " size: " << p.size()
368                 << ") name: " << p.name()
369                 << nl;
370         }
371     }
373     Info<< "\nEnd\n" << endl;
375     return 0;
379 // ************************************************************************* //