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 | Copyright (C) 2011-2014 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  refineHexMesh
26 
27 Description
28  Refines a hex mesh by 2x2x2 cell splitting.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "fvMesh.H"
33 #include "pointMesh.H"
34 #include "argList.H"
35 #include "Time.H"
36 #include "hexRef8.H"
37 #include "cellSet.H"
38 #include "OFstream.H"
39 #include "meshTools.H"
40 #include "IFstream.H"
41 #include "polyTopoChange.H"
42 #include "mapPolyMesh.H"
43 #include "volMesh.H"
44 #include "surfaceMesh.H"
45 #include "volFields.H"
46 #include "surfaceFields.H"
47 #include "pointFields.H"
48 #include "ReadFields.H"
49 
50 using namespace Foam;
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 int main(int argc, char *argv[])
55 {
56  #include "addOverwriteOption.H"
57  #include "addRegionOption.H"
58  argList::validArgs.append("cellSet");
60  (
61  "minSet",
62  "remove cells from input cellSet to keep to 2:1 ratio"
63  " (default is to extend set)"
64  );
65 
66  #include "setRootCase.H"
67  #include "createTime.H"
68  runTime.functionObjects().off();
69  #include "createNamedMesh.H"
70  const word oldInstance = mesh.pointsInstance();
71 
72  word cellSetName(args.args()[1]);
73  const bool overwrite = args.optionFound("overwrite");
74 
75  const bool minSet = args.optionFound("minSet");
76 
77  Info<< "Reading cells to refine from cellSet " << cellSetName
78  << nl << endl;
79 
80  cellSet cellsToRefine(mesh, cellSetName);
81 
82  Info<< "Read " << returnReduce(cellsToRefine.size(), sumOp<label>())
83  << " cells to refine from cellSet " << cellSetName << nl
84  << endl;
85 
86 
87  // Read objects in time directory
88  IOobjectList objects(mesh, runTime.timeName());
89 
90  // Read vol fields.
91 
93  ReadFields(mesh, objects, vsFlds);
94 
96  ReadFields(mesh, objects, vvFlds);
97 
99  ReadFields(mesh, objects, vstFlds);
100 
101  PtrList<volSymmTensorField> vsymtFlds;
102  ReadFields(mesh, objects, vsymtFlds);
103 
105  ReadFields(mesh, objects, vtFlds);
106 
107  // Read surface fields.
108 
110  ReadFields(mesh, objects, ssFlds);
111 
113  ReadFields(mesh, objects, svFlds);
114 
116  ReadFields(mesh, objects, sstFlds);
117 
119  ReadFields(mesh, objects, ssymtFlds);
120 
122  ReadFields(mesh, objects, stFlds);
123 
124  // Read point fields
126  ReadFields(pointMesh::New(mesh), objects, psFlds);
127 
129  ReadFields(pointMesh::New(mesh), objects, pvFlds);
130 
131 
132  // Construct refiner without unrefinement. Read existing point/cell level.
134 
135  // Some stats
136  Info<< "Read mesh:" << nl
137  << " cells:" << mesh.globalData().nTotalCells() << nl
138  << " faces:" << mesh.globalData().nTotalFaces() << nl
139  << " points:" << mesh.globalData().nTotalPoints() << nl
140  << " cellLevel :"
141  << " min:" << gMin(meshCutter.cellLevel())
142  << " max:" << gMax(meshCutter.cellLevel()) << nl
143  << " pointLevel :"
144  << " min:" << gMin(meshCutter.pointLevel())
145  << " max:" << gMax(meshCutter.pointLevel()) << nl
146  << endl;
147 
148 
149  // Maintain 2:1 ratio
150  labelList newCellsToRefine
151  (
152  meshCutter.consistentRefinement
153  (
154  cellsToRefine.toc(),
155  !minSet // extend set
156  )
157  );
158 
159  // Mesh changing engine.
160  polyTopoChange meshMod(mesh);
161 
162  // Play refinement commands into mesh changer.
163  meshCutter.setRefinement(newCellsToRefine, meshMod);
164 
165  if (!overwrite)
166  {
167  runTime++;
168  }
169 
170  // Create mesh, return map from old to new mesh.
171  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
172 
173  // Update fields
174  mesh.updateMesh(map);
175 
176  // Update numbering of cells/vertices.
177  meshCutter.updateMesh(map);
178 
179  // Optionally inflate mesh
180  if (map().hasMotionPoints())
181  {
182  mesh.movePoints(map().preMotionPoints());
183  }
184 
185  Info<< "Refined from " << returnReduce(map().nOldCells(), sumOp<label>())
186  << " to " << mesh.globalData().nTotalCells() << " cells." << nl << endl;
187 
188  if (overwrite)
189  {
190  mesh.setInstance(oldInstance);
191  meshCutter.setInstance(oldInstance);
192  }
193  Info<< "Writing mesh to " << runTime.timeName() << endl;
194 
195  mesh.write();
196  meshCutter.write();
197 
198  Info<< "End\n" << endl;
199 
200  return 0;
201 }
202 
203 
204 // ************************************************************************* //
Foam::surfaceFields.
wordList ReadFields(const Mesh &mesh, const IOobjectList &objects, PtrList< GeoField > &fields, const bool syncPar=true)
Read all fields of the specified type.
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
Type gMin(const FieldField< Field, Type > &f)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:154
Field reading functions for post-processing utilities.
static const pointMesh & New(const polyMesh &mesh)
void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
Definition: meshCutter.C:972
const stringList & args() const
Return arguments.
Definition: argListI.H:66
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:800
dynamicFvMesh & mesh
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
Refinement of (split) hexes using polyTopoChange.
Definition: hexRef8.H:64
A class for handling words, derived from string.
Definition: word.H:59
Cuts (splits) cells.
Definition: meshCutter.H:134
label nTotalFaces() const
Return total number of faces in decomposed mesh. Not.
label nTotalPoints() const
Return total number of points in decomposed mesh. Not.
void setRefinement(const cellCuts &cuts, polyTopoChange &meshMod)
Do actual cutting with cut description. Inserts mesh changes.
Definition: meshCutter.C:524
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1140
static const char nl
Definition: Ostream.H:262
Type gMax(const FieldField< Field, Type > &f)
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:32
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:721
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:62
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:772
A collection of cell labels.
Definition: cellSet.H:48
Direct mesh changes based on v1.3 polyTopoChange syntax.
messageStream Info
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:83
Foam::argList args(argc, argv)
label nTotalCells() const
Return total number of cells in decomposed mesh.
virtual bool write() const
Write mesh using IO settings from time.
Definition: fvMesh.C:870
Namespace for OpenFOAM.