refineMesh.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-2025 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  refineMesh
26 
27 Description
28  Utility to refine cells in multiple directions.
29 
30  Command-line option handling:
31  - If -all specified or no refineMeshDict exists or, refine all cells
32  - If -dict <file> specified refine according to <file>
33  - If refineMeshDict exists refine according to refineMeshDict
34 
35  When the refinement or all cells is selected apply 3D refinement for 3D
36  cases and 2D refinement for 2D cases.
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #include "argList.H"
41 #include "Time.H"
42 #include "polyMesh.H"
43 #include "zoneGenerator.H"
44 #include "cellSet.H"
45 #include "hexRef8.H"
46 #include "multiDirRefinement.H"
47 #include "polyTopoChange.H"
48 #include "syncTools.H"
49 #include "systemDict.H"
50 
51 using namespace Foam;
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 // Max cos angle for edges to be considered aligned with axis.
56 static const scalar edgeTol = 1e-3;
57 
58 // Print edge statistics on mesh.
59 void printEdgeStats(const polyMesh& mesh)
60 {
61  label nX = 0;
62  label nY = 0;
63  label nZ = 0;
64 
65  scalar minX = great;
66  scalar maxX = -great;
67  static const vector x(1, 0, 0);
68 
69  scalar minY = great;
70  scalar maxY = -great;
71  static const vector y(0, 1, 0);
72 
73  scalar minZ = great;
74  scalar maxZ = -great;
75  static const vector z(0, 0, 1);
76 
77  scalar minOther = great;
78  scalar maxOther = -great;
79 
81 
82  const edgeList& edges = mesh.edges();
83 
84  forAll(edges, edgeI)
85  {
86  if (isMasterEdge[edgeI])
87  {
88  const edge& e = edges[edgeI];
89 
90  vector eVec(e.vec(mesh.points()));
91 
92  scalar eMag = mag(eVec);
93 
94  eVec /= eMag;
95 
96  if (mag(eVec & x) > 1-edgeTol)
97  {
98  minX = min(minX, eMag);
99  maxX = max(maxX, eMag);
100  nX++;
101  }
102  else if (mag(eVec & y) > 1-edgeTol)
103  {
104  minY = min(minY, eMag);
105  maxY = max(maxY, eMag);
106  nY++;
107  }
108  else if (mag(eVec & z) > 1-edgeTol)
109  {
110  minZ = min(minZ, eMag);
111  maxZ = max(maxZ, eMag);
112  nZ++;
113  }
114  else
115  {
116  minOther = min(minOther, eMag);
117  maxOther = max(maxOther, eMag);
118  }
119  }
120  }
121 
122  label nEdges = mesh.nEdges();
123  reduce(nEdges, sumOp<label>());
124  reduce(nX, sumOp<label>());
125  reduce(nY, sumOp<label>());
126  reduce(nZ, sumOp<label>());
127 
128  reduce(minX, minOp<scalar>());
129  reduce(maxX, maxOp<scalar>());
130 
131  reduce(minY, minOp<scalar>());
132  reduce(maxY, maxOp<scalar>());
133 
134  reduce(minZ, minOp<scalar>());
135  reduce(maxZ, maxOp<scalar>());
136 
137  reduce(minOther, minOp<scalar>());
138  reduce(maxOther, maxOp<scalar>());
139 
140  Info<< "Mesh edge statistics:" << nl;
141  if (minX > great || maxX > -great)
142  {
143  Info<< " x aligned : number:" << nX << "\tminLen:" << minX
144  << "\tmaxLen:" << maxX << nl;
145  }
146  if (minY > great || maxY > -great)
147  {
148  Info<< " y aligned : number:" << nY << "\tminLen:" << minY
149  << "\tmayLen:" << maxY << nl;
150  }
151  if (minZ > great || maxZ > -great)
152  {
153  Info<< " z aligned : number:" << nZ << "\tminLen:" << minZ
154  << "\tmazLen:" << maxZ << nl;
155  }
156  if (minOther > great || maxOther > -great)
157  {
158  Info<< " other : number:" << nEdges - nX - nY - nZ
159  << "\tminLen:" << minOther
160  << "\tmaxLen:" << maxOther << nl;
161  }
162 
163  Info<< endl;
164 }
165 
166 
167 void refineZone
168 (
169  polyMesh& mesh,
170  const labelList& refCells,
171  autoPtr<hexRef8>& meshCutterPtr,
172  const dictionary& refineDict,
173  const dictionary& zoneDict
174 )
175 {
176  const label nCells0 = mesh.globalData().nTotalCells();
177 
178  if (meshCutterPtr.valid())
179  {
180  hexRef8& meshCutter = meshCutterPtr();
181 
182  // Maintain 2:1 ratio
183  const labelList newCellsToRefine
184  (
185  meshCutter.consistentRefinement
186  (
187  refCells,
188  true // extend set
189  )
190  );
191 
192  // Mesh changing engine.
193  polyTopoChange meshMod(mesh);
194 
195  // Play refinement commands into mesh changer.
196  meshCutter.setRefinement(newCellsToRefine, meshMod);
197 
198  // Create mesh, return map from old to new mesh.
199  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh);
200 
201  // Update mesh objects
202  mesh.topoChange(map);
203 
204  // Update hexRef8 cells/vertices
205  meshCutter.topoChange(map);
206  }
207  else
208  {
210  (
211  mesh,
212  refCells,
213  refineDict,
214  zoneDict.isDict("coordinates")
215  ? zoneDict.subDict("coordinates")
216  : refineDict.subDict("coordinates")
217  );
218  }
219 
220  Info<< " total number of cell increased from " << nCells0
221  << " to " << mesh.globalData().nTotalCells()
222  << nl << endl;
223 }
224 
225 
226 int main(int argc, char *argv[])
227 {
228  argList::addNote("Refine cells in multiple directions");
229 
230  #include "addNoOverwriteOption.H"
231  #include "addMeshOption.H"
232  #include "addRegionOption.H"
233  #include "addDictOption.H"
234 
236  (
237  "all",
238  "Refine all cells"
239  );
240 
241  #include "setRootCase.H"
243  #include "createSpecifiedPolyMesh.H"
244  const word oldInstance = mesh.pointsInstance();
245 
246  // Some stats
247  Info<< "Read mesh:" << nl
248  << " cells:" << mesh.globalData().nTotalCells() << nl
249  << " faces:" << mesh.globalData().nTotalFaces() << nl
250  << " points:" << mesh.globalData().nTotalPoints() << nl << endl;
251  printEdgeStats(mesh);
252 
253  const bool refineAllCells = args.optionFound("all");
254  #include "setNoOverwrite.H"
255 
256  autoPtr<hexRef8> meshCutterPtr;
257 
258  if (refineAllCells)
259  {
260  Info<< "Refining all cells" << nl << endl;
261 
262  // Select all cells
263  labelList refCells(identityMap(mesh.nCells()));
264 
265  dictionary refineDict;
266 
267  if (mesh.nGeometricD() == 3)
268  {
269  Info<< "3D case; refining all directions" << nl << endl;
270 
271  wordList directions(3);
272  directions[0] = "e1";
273  directions[1] = "e2";
274  directions[2] = "e3";
275  refineDict.add("directions", directions);
276 
277  // Use hex cutter
278  refineDict.add("useHexTopology", "true");
279  }
280  else
281  {
282  const Vector<label> dirs(mesh.geometricD());
283 
284  wordList directions(2);
285 
286  if (dirs.x() == -1)
287  {
288  Info<< "2D case; refining in directions y,z\n" << endl;
289  directions[0] = "e2";
290  directions[1] = "e3";
291  }
292  else if (dirs.y() == -1)
293  {
294  Info<< "2D case; refining in directions x,z\n" << endl;
295  directions[0] = "e1";
296  directions[1] = "e3";
297  }
298  else
299  {
300  Info<< "2D case; refining in directions x,y\n" << endl;
301  directions[0] = "e1";
302  directions[1] = "e2";
303  }
304 
305  refineDict.add("directions", directions);
306 
307  // Use standard cutter
308  refineDict.add("useHexTopology", "false");
309  }
310 
311  refineDict.add("coordinateSystem", "global");
312 
313  dictionary coeffDict;
314  coeffDict.add("e1", vector(1, 0, 0));
315  coeffDict.add("e2", vector(0, 1, 0));
316  refineDict.add("globalCoeffs", coeffDict);
317 
318  refineDict.add("geometricCut", "false");
319  refineDict.add("writeMesh", "false");
320 
321  multiDirRefinement multiRef(mesh, refCells, refineDict, refineDict);
322  }
323  else
324  {
325  const dictionary refineDict
326  (
327  systemDict("refineMeshDict", args, mesh)
328  );
329 
330  const bool hexRef8Refine(refineDict.lookupOrDefault("hexRef8", false));
331 
332  // Construct refiner without unrefinement.
333  // Read existing point/cell level.
334  if (hexRef8Refine)
335  {
336  meshCutterPtr = new hexRef8(mesh);
337  hexRef8& meshCutter = meshCutterPtr();
338 
339  Info<< "Refining using hexRef8" << nl
340  << " cellLevel :"
341  << " min:" << gMin(meshCutter.cellLevel())
342  << " max:" << gMax(meshCutter.cellLevel()) << nl
343  << " pointLevel :"
344  << " min:" << gMin(meshCutter.pointLevel())
345  << " max:" << gMax(meshCutter.pointLevel()) << nl
346  << endl;
347  }
348 
349  if (refineDict.found("set"))
350  {
351  const word setName(refineDict.lookup("set"));
352  const cellSet cells(mesh, setName);
353 
354  Info<< "Refining "
356  << " cells in set "
357  << cells.instance()/cells.local()/cells.name()
358  << endl;
359 
360  refineZone
361  (
362  mesh,
363  cells.toc(),
364  meshCutterPtr,
365  refineDict,
366  refineDict
367  );
368  }
369  else if (refineDict.found("zone"))
370  {
371  labelList refCells;
372 
373  if (refineDict.isDict("zone"))
374  {
376  (
378  (
379  "zone",
381  mesh,
382  refineDict.subDict("zone")
383  )
384  );
385 
386  refCells = zg->generate().cZone();
387 
388  Info<< "Refining "
389  << returnReduce(refCells.size(), sumOp<label>())
390  << " cells in zone " << zg->zoneName()
391  << " of type " << zg->type() << endl;
392  }
393  else
394  {
395  const word cellZoneName(refineDict.lookup("zone"));
396  refCells = mesh.cellZones()[cellZoneName];
397 
398  Info<< "Refining "
399  << returnReduce(refCells.size(), sumOp<label>())
400  << " cells in zone " << cellZoneName << endl;
401  }
402 
403  refineZone
404  (
405  mesh,
406  refCells,
407  meshCutterPtr,
408  refineDict,
409  refineDict.isDict("zone")
410  ? refineDict.subDict("zone")
411  : refineDict
412  );
413  }
414  else if (refineDict.found("zones"))
415  {
416  const dictionary& zones = refineDict.subDict("zones");
417 
418  forAllConstIter(dictionary, zones, iter)
419  {
420  const word& name = iter().keyword();
421  const dictionary& zoneDict = iter().dict();
422 
424  (
426  (
427  name,
429  mesh,
430  zoneDict
431  )
432  );
433 
434  const labelList refCells(zg->generate().cZone());
435 
436  Info<< "Refining "
437  << returnReduce(refCells.size(), sumOp<label>())
438  << " cells in zone " << zg->zoneName()
439  << " of type " << zg->type() << endl;
440 
441  refineZone
442  (
443  mesh,
444  refCells,
445  meshCutterPtr,
446  refineDict,
447  zoneDict
448  );
449  }
450  }
451  }
452 
453  printEdgeStats(mesh);
454 
455  if (overwrite)
456  {
457  mesh.setInstance(oldInstance);
458 
459  if (meshCutterPtr.valid())
460  {
461  meshCutterPtr->setInstance(oldInstance);
462  }
463  }
464  else
465  {
466  runTime++;
467  }
468 
469  Info<< "Writing mesh to ";
470  mesh.write();
472 
473  if (meshCutterPtr.valid())
474  {
475  Info<< "Writing hexRef8 refinement level files to "
476  << mesh.facesInstance()/mesh.meshDir() << endl;
477  meshCutterPtr->write();
478  }
479 
480  Info<< "\nEnd\n" << endl;
481 
482  return 0;
483 }
484 
485 
486 // ************************************************************************* //
scalar y
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
A bit-packed bool list.
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
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
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Definition: autoPtrI.H:83
A collection of cell labels.
Definition: cellSet.H:51
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:740
bool isDict(const word &) const
Check if entry is a sub-dictionary.
Definition: dictionary.C:803
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:849
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:1020
T lookupOrDefault(const word &, const T &, const bool writeDefault=writeOptionalEntries > 0) const
Find and return a T, if not found return the given default.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:539
Set of directions for each cell in the mesh. Either uniform and size=1 or one set of directions per c...
Definition: directions.H:67
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1785
virtual void topoChange(const polyTopoChangeMap &map)
Update mesh corresponding to the given map.
Definition: fvMesh.C:1369
label nTotalFaces() const
Return total number of faces in decomposed mesh. Not.
label nTotalPoints() const
Return total number of points in decomposed mesh. Not.
label nTotalCells() const
Return total number of cells in decomposed mesh.
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
Does multiple pass refinement to refine cells in multiple directions.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:994
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:1023
const cellZoneList & cellZones() const
Return cell zones.
Definition: polyMesh.H:449
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: polyMesh.C:982
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1521
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:988
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1331
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:91
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:1012
Direct mesh changes based on v1.3 polyTopoChange syntax.
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
label nCells() const
label nEdges() const
static PackedBoolList getMasterEdges(const polyMesh &)
Get per edge whether it is uncoupled or a master of a.
Definition: syncTools.C:109
A class for handling words, derived from string.
Definition: word.H:62
static autoPtr< zoneGenerator > New(const word &name, const polyMesh &mesh, const dictionary &dict)
Select constructed from name, mesh and dictionary.
Definition: zoneGenerator.C:66
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
int main(int argc, char *argv[])
Definition: financialFoam.C:44
const cellShapeList & cells
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:106
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
messageStream Info
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
Type gMin(const FieldField< Field, Type > &f)
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
IOdictionary systemDict(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion, const fileName &path=fileName::null)
Definition: systemDict.C:93
static const char nl
Definition: Ostream.H:267
Type gMax(const FieldField< Field, Type > &f)
const bool overwrite
Definition: setNoOverwrite.H:1
Foam::argList args(argc, argv)