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-2023 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 "polyTopoChangeMap.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"
73 
74  const bool overwrite = args.optionFound("overwrite");
75  const bool minSet = args.optionFound("minSet");
76  const bool fields = !args.optionFound("noFields");
77 
78  #include "createNamedMesh.H"
79  const word oldInstance = mesh.pointsInstance();
80 
81  word cellSetName(args.args()[1]);
82 
83  Info<< "Reading cells to refine from cellSet " << cellSetName
84  << nl << endl;
85 
86  cellSet cellsToRefine(mesh, cellSetName);
87 
88  Info<< "Read " << returnReduce(cellsToRefine.size(), sumOp<label>())
89  << " cells to refine from cellSet " << cellSetName << nl
90  << endl;
91 
92 
93  // Read objects in time directory
94  IOobjectList objects(mesh, runTime.name());
95 
96  if (fields) Info<< "Reading geometric fields" << nl << endl;
97 
98  #include "readVolFields.H"
99  #include "readSurfaceFields.H"
100  #include "readPointFields.H"
101 
102  Info<< endl;
103 
104 
105  // Construct refiner without unrefinement. Read existing point/cell level.
106  hexRef8 meshCutter(mesh);
107 
108  // Some stats
109  Info<< "Read mesh:" << nl
110  << " cells:" << mesh.globalData().nTotalCells() << nl
111  << " faces:" << mesh.globalData().nTotalFaces() << nl
112  << " points:" << mesh.globalData().nTotalPoints() << nl
113  << " cellLevel :"
114  << " min:" << gMin(meshCutter.cellLevel())
115  << " max:" << gMax(meshCutter.cellLevel()) << nl
116  << " pointLevel :"
117  << " min:" << gMin(meshCutter.pointLevel())
118  << " max:" << gMax(meshCutter.pointLevel()) << nl
119  << endl;
120 
121 
122  // Maintain 2:1 ratio
123  labelList newCellsToRefine
124  (
125  meshCutter.consistentRefinement
126  (
127  cellsToRefine.toc(),
128  !minSet // extend set
129  )
130  );
131 
132  // Mesh changing engine.
133  polyTopoChange meshMod(mesh);
134 
135  // Play refinement commands into mesh changer.
136  meshCutter.setRefinement(newCellsToRefine, meshMod);
137 
138  if (!overwrite)
139  {
140  runTime++;
141  }
142 
143  // Create mesh, return map from old to new mesh.
144  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh, false);
145 
146  // Update fields
147  mesh.topoChange(map);
148 
149  // Update numbering of cells/vertices.
150  meshCutter.topoChange(map);
151 
152  // Optionally inflate mesh
153  if (map().hasMotionPoints())
154  {
155  mesh.setPoints(map().preMotionPoints());
156  }
157 
158  Info<< "Refined from " << returnReduce(map().nOldCells(), sumOp<label>())
159  << " to " << mesh.globalData().nTotalCells() << " cells." << nl << endl;
160 
161  if (overwrite)
162  {
163  mesh.setInstance(oldInstance);
164  meshCutter.setInstance(oldInstance);
165  }
166  Info<< "Writing mesh to " << runTime.name() << endl;
167 
168  mesh.write();
169  meshCutter.write();
170 
171  Info<< "End\n" << endl;
172 
173  return 0;
174 }
175 
176 
177 // ************************************************************************* //
Field reading functions for post-processing utilities.
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:53
virtual Ostream & write(const char)=0
Write character.
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
const stringList & args() const
Return arguments.
Definition: argListI.H:72
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:153
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A collection of cell labels.
Definition: cellSet.H:51
Refinement of (split) hexes using polyTopoChange.
Definition: hexRef8.H:65
Cuts (splits) cells.
Definition: meshCutter.H:137
void topoChange(const polyTopoChangeMap &)
Force recalculation of locally stored data on topological change.
Definition: meshCutter.C:998
void setRefinement(const cellCuts &cuts, polyTopoChange &meshMod)
Do actual cutting with cut description. Inserts mesh changes.
Definition: meshCutter.C:525
Direct mesh changes based on v1.3 polyTopoChange syntax.
A class for handling words, derived from string.
Definition: word.H:62
int main(int argc, char *argv[])
Definition: financialFoam.C:44
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:230
Namespace for OpenFOAM.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Type gMin(const FieldField< Field, Type > &f)
static const char nl
Definition: Ostream.H:260
Type gMax(const FieldField< Field, Type > &f)
objects
Foam::argList args(argc, argv)
Foam::surfaceFields.