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-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  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 "polyMesh.H"
42 #include "Time.H"
43 #include "cellSet.H"
44 #include "multiDirRefinement.H"
45 #include "labelIOList.H"
46 #include "IOdictionary.H"
47 #include "syncTools.H"
48 #include "systemDict.H"
49 
50 using namespace Foam;
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 // Max cos angle for edges to be considered aligned with axis.
55 static const scalar edgeTol = 1e-3;
56 
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 
80  PackedBoolList isMasterEdge(syncTools::getMasterEdges(mesh));
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 
141  Info<< "Mesh edge statistics:" << nl
142  << " x aligned : number:" << nX << "\tminLen:" << minX
143  << "\tmaxLen:" << maxX << nl
144  << " y aligned : number:" << nY << "\tminLen:" << minY
145  << "\tmaxLen:" << maxY << nl
146  << " z aligned : number:" << nZ << "\tminLen:" << minZ
147  << "\tmaxLen:" << maxZ << nl
148  << " other : number:" << nEdges - nX - nY - nZ
149  << "\tminLen:" << minOther
150  << "\tmaxLen:" << maxOther << nl << endl;
151 }
152 
153 
154 int main(int argc, char *argv[])
155 {
157  (
158  "refine cells in multiple directions"
159  );
160 
161  #include "addOverwriteOption.H"
162  #include "addRegionOption.H"
163  #include "addDictOption.H"
164 
166  (
167  "all",
168  "Refine all cells"
169  );
170 
171  #include "setRootCase.H"
173  #include "createNamedPolyMesh.H"
174  const word oldInstance = mesh.pointsInstance();
175 
176  printEdgeStats(mesh);
177 
178  //
179  // Read/construct control dictionary
180  //
181 
182  const bool readDict = args.optionFound("dict");
183  const bool refineAllCells = args.optionFound("all");
184  const bool overwrite = args.optionFound("overwrite");
185 
186  // List of cells to refine
187  labelList refCells;
188 
189  // Dictionary to control refinement
190  const word dictName("refineMeshDict");
191  typeIOobject<IOdictionary> dictIO(systemDictIO(dictName, args, runTime));
192  dictionary refineDict;
193  if (readDict)
194  {
195  if (dictIO.headerOk())
196  {
197  Info<< "Refining according to "
199  << nl << endl;
200  refineDict = IOdictionary(dictIO);
201  }
202  else
203  {
205  << "Cannot open specified refinement dictionary "
207  << exit(FatalError);
208  }
209  }
210  else if (!refineAllCells)
211  {
212  if (dictIO.headerOk())
213  {
214  Info<< "Refining according to "
216  << nl << endl;
217  refineDict = IOdictionary(dictIO);
218  }
219  else
220  {
221  Info<< "Refinement dictionary "
223  << " not found" << nl << endl;
224  }
225  }
226 
227  if (refineDict.size())
228  {
229  const word setName(refineDict.lookup("set"));
230 
231  cellSet cells(mesh, setName);
232 
233  Info<< "Read " << returnReduce(cells.size(), sumOp<label>())
234  << " cells from cellSet "
235  << cells.instance()/cells.local()/cells.name()
236  << endl << endl;
237 
238  refCells = cells.toc();
239  }
240  else
241  {
242  Info<< "Refining all cells" << nl << endl;
243 
244  // Select all cells
245  refCells = identityMap(mesh.nCells());
246 
247  if (mesh.nGeometricD() == 3)
248  {
249  Info<< "3D case; refining all directions" << nl << endl;
250 
251  wordList directions(3);
252  directions[0] = "e1";
253  directions[1] = "e2";
254  directions[2] = "e3";
255  refineDict.add("directions", directions);
256 
257  // Use hex cutter
258  refineDict.add("useHexTopology", "true");
259  }
260  else
261  {
262  const Vector<label> dirs(mesh.geometricD());
263 
264  wordList directions(2);
265 
266  if (dirs.x() == -1)
267  {
268  Info<< "2D case; refining in directions y,z\n" << endl;
269  directions[0] = "e2";
270  directions[1] = "e3";
271  }
272  else if (dirs.y() == -1)
273  {
274  Info<< "2D case; refining in directions x,z\n" << endl;
275  directions[0] = "e1";
276  directions[1] = "e3";
277  }
278  else
279  {
280  Info<< "2D case; refining in directions x,y\n" << endl;
281  directions[0] = "e1";
282  directions[1] = "e2";
283  }
284 
285  refineDict.add("directions", directions);
286 
287  // Use standard cutter
288  refineDict.add("useHexTopology", "false");
289  }
290 
291  refineDict.add("coordinateSystem", "global");
292 
293  dictionary coeffsDict;
294  coeffsDict.add("e1", vector(1, 0, 0));
295  coeffsDict.add("e2", vector(0, 1, 0));
296  refineDict.add("globalCoeffs", coeffsDict);
297 
298  refineDict.add("geometricCut", "false");
299  refineDict.add("writeMesh", "false");
300  }
301 
302 
303  string oldTimeName(runTime.name());
304 
305  if (!overwrite)
306  {
307  runTime++;
308  }
309 
310 
311  // Multi-directional refinement (does multiple iterations)
312  multiDirRefinement multiRef(mesh, refCells, refineDict);
313 
314 
315  // Write resulting mesh
316  if (overwrite)
317  {
318  mesh.setInstance(oldInstance);
319  }
320  mesh.write();
321 
322 
323  // Get list of cell splits.
324  // (is for every cell in old mesh the cells they have been split into)
325  const labelListList& oldToNew = multiRef.addedCells();
326 
327 
328  // Create cellSet with added cells for easy inspection
329  cellSet newCells(mesh, "refinedCells", refCells.size());
330 
331  forAll(oldToNew, oldCelli)
332  {
333  const labelList& added = oldToNew[oldCelli];
334 
335  forAll(added, i)
336  {
337  newCells.insert(added[i]);
338  }
339  }
340 
341  Info<< "Writing refined cells ("
342  << returnReduce(newCells.size(), sumOp<label>())
343  << ") to cellSet "
344  << newCells.instance()/newCells.local()/newCells.name()
345  << endl << endl;
346 
347  newCells.write();
348 
349 
350 
351  //
352  // Invert cell split to construct map from new to old
353  //
354 
355  labelIOList newToOld
356  (
357  IOobject
358  (
359  "cellMap",
360  runTime.name(),
362  mesh,
365  ),
366  mesh.nCells()
367  );
368  newToOld.note() =
369  "From cells in mesh at "
370  + runTime.name()
371  + " to cells in mesh at "
372  + oldTimeName;
373 
374 
375  forAll(oldToNew, oldCelli)
376  {
377  const labelList& added = oldToNew[oldCelli];
378 
379  if (added.size())
380  {
381  forAll(added, i)
382  {
383  newToOld[added[i]] = oldCelli;
384  }
385  }
386  else
387  {
388  // Unrefined cell
389  newToOld[oldCelli] = oldCelli;
390  }
391  }
392 
393  Info<< "Writing map from new to old cell to "
394  << newToOld.relativeObjectPath() << nl << endl;
395 
396  newToOld.write();
397 
398  printEdgeStats(mesh);
399 
400  Info<< "End\n" << endl;
401 
402  return 0;
403 }
404 
405 
406 // ************************************************************************* //
scalar y
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
virtual Ostream & write(const char)=0
Write character.
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
A collection of cell labels.
Definition: cellSet.H:51
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:1169
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
Does multiple pass refinement to refine cells in multiple directions.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:1054
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:1019
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1360
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:271
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:78
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:1043
label nEdges() const
label nCells() const
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
virtual bool write(const bool write=true) const
Write using setting from DB.
static PackedBoolList getMasterEdges(const polyMesh &)
Get per edge whether it is uncoupled or a master of a.
Definition: syncTools.C:109
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:531
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
int main(int argc, char *argv[])
Definition: financialFoam.C:44
const cellShapeList & cells
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const doubleScalar e
Definition: doubleScalar.H:105
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:251
IOobject systemDictIO(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion)
Definition: systemDict.C:33
messageStream Info
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)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
error FatalError
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
static const char nl
Definition: Ostream.H:260
Foam::argList args(argc, argv)
Trait for obtaining global write status.
Definition: IOobject.H:511