Forward compatibility: flex
[foam-extend-3.2.git] / applications / utilities / mesh / advanced / refineHexMesh / refineHexMesh.C
blob3fb3d26a46096741da8cebee2dd673d96f3778d9
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 Description
25     Refines a hex mesh by 2x2x2 cell splitting.
27 \*---------------------------------------------------------------------------*/
29 #include "fvMesh.H"
30 #include "pointMesh.H"
31 #include "argList.H"
32 #include "foamTime.H"
33 #include "hexRef8.H"
34 #include "cellSet.H"
35 #include "OFstream.H"
36 #include "meshTools.H"
37 #include "IFstream.H"
38 #include "directTopoChange.H"
39 #include "mapPolyMesh.H"
40 #include "volMesh.H"
41 #include "surfaceMesh.H"
42 #include "volFields.H"
43 #include "surfaceFields.H"
44 #include "pointFields.H"
45 #include "ReadFields.H"
47 using namespace Foam;
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 // Main program:
52 int main(int argc, char *argv[])
54     argList::validOptions.insert("overwrite", "");
55     argList::validArgs.append("cellSet");
56 #   include "setRootCase.H"
57 #   include "createTime.H"
58     runTime.functionObjects().off();
59 #   include "createMesh.H"
60     const word oldInstance = mesh.pointsInstance();
62     const pointMesh& pMesh = pointMesh::New(mesh);
64     word cellSetName(args.args()[1]);
65     bool overwrite = args.optionFound("overwrite");
67     Info<< "Reading cells to refine from cellSet " << cellSetName
68         << nl << endl;
70     cellSet cellsToRefine(mesh, cellSetName);
72     Info<< "Read " << returnReduce(cellsToRefine.size(), sumOp<label>())
73         << " cells to refine from cellSet " << cellSetName << nl
74         << endl;
77     // Read objects in time directory
78     IOobjectList objects(mesh, runTime.timeName());
80     // Read vol fields.
82     PtrList<volScalarField> vsFlds;
83     ReadFields(mesh, objects, vsFlds);
85     PtrList<volVectorField> vvFlds;
86     ReadFields(mesh, objects, vvFlds);
88     PtrList<volSphericalTensorField> vstFlds;
89     ReadFields(mesh, objects, vstFlds);
91     PtrList<volSymmTensorField> vsymtFlds;
92     ReadFields(mesh, objects, vsymtFlds);
94     PtrList<volTensorField> vtFlds;
95     ReadFields(mesh, objects, vtFlds);
97     // Read surface fields.
99     PtrList<surfaceScalarField> ssFlds;
100     ReadFields(mesh, objects, ssFlds);
102     PtrList<surfaceVectorField> svFlds;
103     ReadFields(mesh, objects, svFlds);
105     PtrList<surfaceSphericalTensorField> sstFlds;
106     ReadFields(mesh, objects, sstFlds);
108     PtrList<surfaceSymmTensorField> ssymtFlds;
109     ReadFields(mesh, objects, ssymtFlds);
111     PtrList<surfaceTensorField> stFlds;
112     ReadFields(mesh, objects, stFlds);
114     // Read point fields
115     PtrList<pointScalarField> psFlds;
116     ReadFields(pMesh, objects, psFlds);
118     PtrList<pointVectorField> pvFlds;
119     ReadFields(pMesh, objects, pvFlds);
123     // Construct refiner without unrefinement. Read existing point/cell level.
124     hexRef8 meshCutter(mesh);
126     // Some stats
127     Info<< "Read mesh:" << nl
128         << "    cells:" << mesh.globalData().nTotalCells() << nl
129         << "    faces:" << mesh.globalData().nTotalFaces() << nl
130         << "    points:" << mesh.globalData().nTotalPoints() << nl
131         << "    cellLevel :"
132         << " min:" << gMin(meshCutter.cellLevel())
133         << " max:" << gMax(meshCutter.cellLevel()) << nl
134         << "    pointLevel :"
135         << " min:" << gMin(meshCutter.pointLevel())
136         << " max:" << gMax(meshCutter.pointLevel()) << nl
137         << endl;
140     // Maintain 2:1 ratio
141     labelList newCellsToRefine
142     (
143         meshCutter.consistentRefinement
144         (
145             cellsToRefine.toc(),
146             true                  // extend set
147         )
148     );
150     // Mesh changing engine.
151     directTopoChange meshMod(mesh);
153     // Play refinement commands into mesh changer.
154     meshCutter.setRefinement(newCellsToRefine, meshMod);
156     if (!overwrite)
157     {
158         runTime++;
159     }
161     // Create mesh, return map from old to new mesh.
162     autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
164     // Update fields
165     mesh.updateMesh(map);
167     // Update numbering of cells/vertices.
168     meshCutter.updateMesh(map);
170     // Optionally inflate mesh
171     if (map().hasMotionPoints())
172     {
173         mesh.movePoints(map().preMotionPoints());
174     }
176     Pout<< "Refined from " << returnReduce(map().nOldCells(), sumOp<label>())
177         << " to " << mesh.globalData().nTotalCells() << " cells."
178         << nl << endl;
180     if (overwrite)
181     {
182         mesh.setInstance(oldInstance);
183         meshCutter.setInstance(oldInstance);
184     }
185     Info<< "Writing mesh to " << runTime.timeName() << endl;
187     mesh.write();
188     meshCutter.write();
190     Info<< "End\n" << endl;
192     return 0;
196 // ************************************************************************* //