refineHexMesh.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2014 OpenFOAM Foundation
9 Copyright (C) 2016 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Application
28 refineHexMesh
29
30Group
31 grpMeshAdvancedUtilities
32
33Description
34 Refine a hex mesh by 2x2x2 cell splitting for the specified cellSet.
35
36\*---------------------------------------------------------------------------*/
37
38#include "fvMesh.H"
39#include "pointMesh.H"
40#include "argList.H"
41#include "Time.H"
42#include "hexRef8.H"
43#include "cellSet.H"
44#include "Fstream.H"
45#include "meshTools.H"
46#include "polyTopoChange.H"
47#include "mapPolyMesh.H"
48#include "volMesh.H"
49#include "surfaceMesh.H"
50#include "volFields.H"
51#include "surfaceFields.H"
52#include "pointFields.H"
53#include "ReadFields.H"
54#include "processorMeshes.H"
55
56using namespace Foam;
57
58// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59
60int main(int argc, char *argv[])
61{
62 argList::addNote
63 (
64 "Refine a hex mesh by 2x2x2 cell splitting for the specified cellSet"
65 );
66 #include "addOverwriteOption.H"
67 #include "addRegionOption.H"
68 argList::addArgument("cellSet");
69 argList::addBoolOption
70 (
71 "minSet",
72 "Remove cells from input cellSet to keep to 2:1 ratio"
73 " (default is to extend set)"
74 );
75
76 argList::noFunctionObjects(); // Never use function objects
77
78 #include "setRootCase.H"
79 #include "createTime.H"
80 #include "createNamedMesh.H"
81
82 const word oldInstance = mesh.pointsInstance();
83
84 word cellSetName(args[1]);
85 const bool overwrite = args.found("overwrite");
86
87 const bool minSet = args.found("minSet");
88
89 Info<< "Reading cells to refine from cellSet " << cellSetName
90 << nl << endl;
91
92 cellSet cellsToRefine(mesh, cellSetName);
93
94 Info<< "Read " << returnReduce(cellsToRefine.size(), sumOp<label>())
95 << " cells to refine from cellSet " << cellSetName << nl
96 << endl;
97
98
99 // Read objects in time directory
100 IOobjectList objects(mesh, runTime.timeName());
101
102 // Read vol fields.
103
105 ReadFields(mesh, objects, vsFlds);
106
108 ReadFields(mesh, objects, vvFlds);
109
111 ReadFields(mesh, objects, vstFlds);
112
114 ReadFields(mesh, objects, vsymtFlds);
115
117 ReadFields(mesh, objects, vtFlds);
118
119 // Read surface fields.
120
122 ReadFields(mesh, objects, ssFlds);
123
125 ReadFields(mesh, objects, svFlds);
126
128 ReadFields(mesh, objects, sstFlds);
129
131 ReadFields(mesh, objects, ssymtFlds);
132
134 ReadFields(mesh, objects, stFlds);
135
136 // Read point fields
138 ReadFields(pointMesh::New(mesh), objects, psFlds);
139
141 ReadFields(pointMesh::New(mesh), objects, pvFlds);
142
143
144 // Construct refiner without unrefinement. Read existing point/cell level.
146
147 // Some stats
148 Info<< "Read mesh:" << nl
149 << " cells:" << mesh.globalData().nTotalCells() << nl
150 << " faces:" << mesh.globalData().nTotalFaces() << nl
151 << " points:" << mesh.globalData().nTotalPoints() << nl
152 << " cellLevel :"
153 << " min:" << gMin(meshCutter.cellLevel())
154 << " max:" << gMax(meshCutter.cellLevel()) << nl
155 << " pointLevel :"
156 << " min:" << gMin(meshCutter.pointLevel())
157 << " max:" << gMax(meshCutter.pointLevel()) << nl
158 << endl;
159
160
161 // Maintain 2:1 ratio
162 labelList newCellsToRefine
163 (
164 meshCutter.consistentRefinement
165 (
166 cellsToRefine.toc(),
167 !minSet // extend set
168 )
169 );
170
171 // Mesh changing engine.
172 polyTopoChange meshMod(mesh);
173
174 // Play refinement commands into mesh changer.
175 meshCutter.setRefinement(newCellsToRefine, meshMod);
176
177 if (!overwrite)
178 {
179 ++runTime;
180 }
181
182 // Create mesh, return map from old to new mesh.
183 autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
184
185 // Update fields
186 mesh.updateMesh(map());
187
188 // Update numbering of cells/vertices.
189 meshCutter.updateMesh(map());
190
191 // Optionally inflate mesh
192 if (map().hasMotionPoints())
193 {
194 mesh.movePoints(map().preMotionPoints());
195 }
196
197 Info<< "Refined from " << returnReduce(map().nOldCells(), sumOp<label>())
198 << " to " << mesh.globalData().nTotalCells() << " cells." << nl << endl;
199
200 if (overwrite)
201 {
202 mesh.setInstance(oldInstance);
203 meshCutter.setInstance(oldInstance);
204 }
205 Info<< "Writing mesh to " << runTime.timeName() << endl;
206
207 mesh.write();
208 meshCutter.write();
209 topoSet::removeFiles(mesh);
210 processorMeshes::removeFiles(mesh);
211
212 Info<< "End\n" << endl;
213
214 return 0;
215}
216
217
218// ************************************************************************* //
Field reading functions for post-processing utilities.
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:59
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:178
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A collection of cell labels.
Definition: cellSet.H:54
Refinement of (split) hexes using polyTopoChange.
Definition: hexRef8.H:68
Cuts (splits) cells.
Definition: meshCutter.H:141
void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
Definition: meshCutter.C:995
void setRefinement(const cellCuts &cuts, polyTopoChange &meshMod)
Do actual cutting with cut description. Inserts mesh changes.
Definition: meshCutter.C:522
Direct mesh changes based on v1.3 polyTopoChange syntax.
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
engineTime & runTime
Required Variables.
Namespace for OpenFOAM.
wordList ReadFields(const typename GeoMesh::Mesh &mesh, const IOobjectList &objects, PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields, const bool syncPar=true, const bool readOldTime=false)
Read Geometric fields of templated type.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
Type gMin(const FieldField< Field, Type > &f)
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
Type gMax(const FieldField< Field, Type > &f)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
Foam::argList args(argc, argv)
Foam::surfaceFields.