All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
refineHexMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 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  );
66  (
67  "noFields",
68  "do not update fields"
69  );
70 
71  #include "setRootCase.H"
72  #include "createTime.H"
74 
75  const bool overwrite = args.optionFound("overwrite");
76  const bool minSet = args.optionFound("minSet");
77  const bool fields = !args.optionFound("noFields");
78 
79  #include "createNamedMesh.H"
80  const word oldInstance = mesh.pointsInstance();
81 
82  word cellSetName(args.args()[1]);
83 
84  Info<< "Reading cells to refine from cellSet " << cellSetName
85  << nl << endl;
86 
87  cellSet cellsToRefine(mesh, cellSetName);
88 
89  Info<< "Read " << returnReduce(cellsToRefine.size(), sumOp<label>())
90  << " cells to refine from cellSet " << cellSetName << nl
91  << endl;
92 
93 
94  // Read objects in time directory
96 
97  if (fields) Info<< "Reading geometric fields" << nl << endl;
98 
99  #include "readVolFields.H"
100  #include "readSurfaceFields.H"
101  #include "readPointFields.H"
102 
103  Info<< endl;
104 
105 
106  // Construct refiner without unrefinement. Read existing point/cell level.
108 
109  // Some stats
110  Info<< "Read mesh:" << nl
111  << " cells:" << mesh.globalData().nTotalCells() << nl
112  << " faces:" << mesh.globalData().nTotalFaces() << nl
113  << " points:" << mesh.globalData().nTotalPoints() << nl
114  << " cellLevel :"
115  << " min:" << gMin(meshCutter.cellLevel())
116  << " max:" << gMax(meshCutter.cellLevel()) << nl
117  << " pointLevel :"
118  << " min:" << gMin(meshCutter.pointLevel())
119  << " max:" << gMax(meshCutter.pointLevel()) << nl
120  << endl;
121 
122 
123  // Maintain 2:1 ratio
124  labelList newCellsToRefine
125  (
126  meshCutter.consistentRefinement
127  (
128  cellsToRefine.toc(),
129  !minSet // extend set
130  )
131  );
132 
133  // Mesh changing engine.
134  polyTopoChange meshMod(mesh);
135 
136  // Play refinement commands into mesh changer.
137  meshCutter.setRefinement(newCellsToRefine, meshMod);
138 
139  if (!overwrite)
140  {
141  runTime++;
142  }
143 
144  // Create mesh, return map from old to new mesh.
145  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
146 
147  // Update fields
148  mesh.updateMesh(map);
149 
150  // Update numbering of cells/vertices.
151  meshCutter.updateMesh(map);
152 
153  // Optionally inflate mesh
154  if (map().hasMotionPoints())
155  {
156  mesh.movePoints(map().preMotionPoints());
157  }
158 
159  Info<< "Refined from " << returnReduce(map().nOldCells(), sumOp<label>())
160  << " to " << mesh.globalData().nTotalCells() << " cells." << nl << endl;
161 
162  if (overwrite)
163  {
164  mesh.setInstance(oldInstance);
165  meshCutter.setInstance(oldInstance);
166  }
167  Info<< "Writing mesh to " << runTime.timeName() << endl;
168 
169  mesh.write();
170  meshCutter.write();
171 
172  Info<< "End\n" << endl;
173 
174  return 0;
175 }
176 
177 
178 // ************************************************************************* //
Foam::surfaceFields.
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
void off()
Switch the function objects off.
Type gMin(const FieldField< Field, Type > &f)
label nTotalPoints() const
Return total number of points in decomposed mesh. Not.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
label nTotalCells() const
Return total number of cells in decomposed mesh.
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:153
Field reading functions for post-processing utilities.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:636
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1093
void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
Definition: meshCutter.C:998
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:802
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:795
label nTotalFaces() const
Return total number of faces in decomposed mesh. Not.
dynamicFvMesh & mesh
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
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1394
void setRefinement(const cellCuts &cuts, polyTopoChange &meshMod)
Do actual cutting with cut description. Inserts mesh changes.
Definition: meshCutter.C:525
static const char nl
Definition: Ostream.H:260
Type gMax(const FieldField< Field, Type > &f)
objects
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:716
const functionObjectList & functionObjects() const
Return the list of function objects.
Definition: Time.H:410
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:52
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
Foam::argList args(argc, argv)
Namespace for OpenFOAM.
const stringList & args() const
Return arguments.
Definition: argListI.H:72