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-2024 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 "volFields.H"
44 #include "pointFields.H"
45 #include "ReadFields.H"
46 
47 using namespace Foam;
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 int main(int argc, char *argv[])
52 {
53  #include "addOverwriteOption.H"
54  #include "addRegionOption.H"
55  argList::validArgs.append("cellSet");
57  (
58  "minSet",
59  "remove cells from input cellSet to keep to 2:1 ratio"
60  " (default is to extend set)"
61  );
63  (
64  "noFields",
65  "do not update fields"
66  );
67 
68  #include "setRootCase.H"
70 
71  const bool overwrite = args.optionFound("overwrite");
72  const bool minSet = args.optionFound("minSet");
73  const bool fields = !args.optionFound("noFields");
74 
76  const word oldInstance = mesh.pointsInstance();
77 
78  word cellSetName(args.args()[1]);
79 
80  Info<< "Reading cells to refine from cellSet " << cellSetName
81  << nl << endl;
82 
83  cellSet cellsToRefine(mesh, cellSetName);
84 
85  Info<< "Read " << returnReduce(cellsToRefine.size(), sumOp<label>())
86  << " cells to refine from cellSet " << cellSetName << nl
87  << endl;
88 
89 
90  // Read objects in time directory
91  IOobjectList objects(mesh, runTime.name());
92 
93  if (fields) Info<< "Reading geometric fields" << nl << endl;
94 
95  #include "readVolFields.H"
96  // Cannot reliable map surfaceFields to new internal faces
97  // #include "readSurfaceFields.H"
98  #include "readPointFields.H"
99 
100  Info<< endl;
101 
102 
103  // Construct refiner without unrefinement. Read existing point/cell level.
104  hexRef8 meshCutter(mesh);
105 
106  // Some stats
107  Info<< "Read mesh:" << nl
108  << " cells:" << mesh.globalData().nTotalCells() << nl
109  << " faces:" << mesh.globalData().nTotalFaces() << nl
110  << " points:" << mesh.globalData().nTotalPoints() << nl
111  << " cellLevel :"
112  << " min:" << gMin(meshCutter.cellLevel())
113  << " max:" << gMax(meshCutter.cellLevel()) << nl
114  << " pointLevel :"
115  << " min:" << gMin(meshCutter.pointLevel())
116  << " max:" << gMax(meshCutter.pointLevel()) << nl
117  << endl;
118 
119 
120  // Maintain 2:1 ratio
121  labelList newCellsToRefine
122  (
123  meshCutter.consistentRefinement
124  (
125  cellsToRefine.toc(),
126  !minSet // extend set
127  )
128  );
129 
130  // Mesh changing engine.
131  polyTopoChange meshMod(mesh);
132 
133  // Play refinement commands into mesh changer.
134  meshCutter.setRefinement(newCellsToRefine, meshMod);
135 
136  if (!overwrite)
137  {
138  runTime++;
139  }
140 
141  // Create mesh, return map from old to new mesh.
142  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh);
143 
144  // Update fields
145  mesh.topoChange(map);
146 
147  // Update numbering of cells/vertices.
148  meshCutter.topoChange(map);
149 
150  Info<< "Refined from " << returnReduce(map().nOldCells(), sumOp<label>())
151  << " to " << mesh.globalData().nTotalCells() << " cells." << nl << endl;
152 
153  if (overwrite)
154  {
155  mesh.setInstance(oldInstance);
156  meshCutter.setInstance(oldInstance);
157  }
158  Info<< "Writing mesh to " << runTime.name() << endl;
159 
160  mesh.write();
161  meshCutter.write();
162 
163  Info<< "End\n" << endl;
164 
165  return 0;
166 }
167 
168 
169 // ************************************************************************* //
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:131
void topoChange(const polyTopoChangeMap &)
Force recalculation of locally stored data on topological change.
Definition: meshCutter.C:916
void setRefinement(const cellCuts &cuts, polyTopoChange &meshMod)
Do actual cutting with cut description. Inserts mesh changes.
Definition: meshCutter.C:470
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:228
Namespace for OpenFOAM.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
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:266
Type gMax(const FieldField< Field, Type > &f)
objects
Foam::argList args(argc, argv)