All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2021 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"
172  #include "createTime.H"
174  #include "createNamedPolyMesh.H"
175  const word oldInstance = mesh.pointsInstance();
176 
177  printEdgeStats(mesh);
178 
179  //
180  // Read/construct control dictionary
181  //
182 
183  const bool readDict = args.optionFound("dict");
184  const bool refineAllCells = args.optionFound("all");
185  const bool overwrite = args.optionFound("overwrite");
186 
187  // List of cells to refine
188  labelList refCells;
189 
190  // Dictionary to control refinement
191  const word dictName("refineMeshDict");
193  dictionary refineDict;
194  if (readDict)
195  {
196  if (dictIO.typeHeaderOk<IOdictionary>(true))
197  {
198  Info<< "Refining according to " << dictIO.path() << nl << endl;
199  refineDict = IOdictionary(dictIO);
200  }
201  else
202  {
204  << "Cannot open specified refinement dictionary "
205  << dictIO.path() << exit(FatalError);
206  }
207  }
208  else if (!refineAllCells)
209  {
210  if (dictIO.typeHeaderOk<IOdictionary>(true))
211  {
212  Info<< "Refining according to " << dictIO.path() << nl << endl;
213  refineDict = IOdictionary(dictIO);
214  }
215  else
216  {
217  Info<< "Refinement dictionary " << dictIO.path() << " not found"
218  << nl << endl;
219  }
220  }
221 
222  if (refineDict.size())
223  {
224  const word setName(refineDict.lookup("set"));
225 
226  cellSet cells(mesh, setName);
227 
228  Info<< "Read " << returnReduce(cells.size(), sumOp<label>())
229  << " cells from cellSet "
230  << cells.instance()/cells.local()/cells.name()
231  << endl << endl;
232 
233  refCells = cells.toc();
234  }
235  else
236  {
237  Info<< "Refining all cells" << nl << endl;
238 
239  // Select all cells
240  refCells = identity(mesh.nCells());
241 
242  if (mesh.nGeometricD() == 3)
243  {
244  Info<< "3D case; refining all directions" << nl << endl;
245 
246  wordList directions(3);
247  directions[0] = "e1";
248  directions[1] = "e2";
249  directions[2] = "e3";
250  refineDict.add("directions", directions);
251 
252  // Use hex cutter
253  refineDict.add("useHexTopology", "true");
254  }
255  else
256  {
257  const Vector<label> dirs(mesh.geometricD());
258 
259  wordList directions(2);
260 
261  if (dirs.x() == -1)
262  {
263  Info<< "2D case; refining in directions y,z\n" << endl;
264  directions[0] = "e2";
265  directions[1] = "e3";
266  }
267  else if (dirs.y() == -1)
268  {
269  Info<< "2D case; refining in directions x,z\n" << endl;
270  directions[0] = "e1";
271  directions[1] = "e3";
272  }
273  else
274  {
275  Info<< "2D case; refining in directions x,y\n" << endl;
276  directions[0] = "e1";
277  directions[1] = "e2";
278  }
279 
280  refineDict.add("directions", directions);
281 
282  // Use standard cutter
283  refineDict.add("useHexTopology", "false");
284  }
285 
286  refineDict.add("coordinateSystem", "global");
287 
288  dictionary coeffsDict;
289  coeffsDict.add("e1", vector(1, 0, 0));
290  coeffsDict.add("e2", vector(0, 1, 0));
291  refineDict.add("globalCoeffs", coeffsDict);
292 
293  refineDict.add("geometricCut", "false");
294  refineDict.add("writeMesh", "false");
295  }
296 
297 
298  string oldTimeName(runTime.timeName());
299 
300  if (!overwrite)
301  {
302  runTime++;
303  }
304 
305 
306  // Multi-directional refinement (does multiple iterations)
307  multiDirRefinement multiRef(mesh, refCells, refineDict);
308 
309 
310  // Write resulting mesh
311  if (overwrite)
312  {
313  mesh.setInstance(oldInstance);
314  }
315  mesh.write();
316 
317 
318  // Get list of cell splits.
319  // (is for every cell in old mesh the cells they have been split into)
320  const labelListList& oldToNew = multiRef.addedCells();
321 
322 
323  // Create cellSet with added cells for easy inspection
324  cellSet newCells(mesh, "refinedCells", refCells.size());
325 
326  forAll(oldToNew, oldCelli)
327  {
328  const labelList& added = oldToNew[oldCelli];
329 
330  forAll(added, i)
331  {
332  newCells.insert(added[i]);
333  }
334  }
335 
336  Info<< "Writing refined cells ("
337  << returnReduce(newCells.size(), sumOp<label>())
338  << ") to cellSet "
339  << newCells.instance()/newCells.local()/newCells.name()
340  << endl << endl;
341 
342  newCells.write();
343 
344 
345 
346  //
347  // Invert cell split to construct map from new to old
348  //
349 
350  labelIOList newToOld
351  (
352  IOobject
353  (
354  "cellMap",
355  runTime.timeName(),
357  mesh,
360  ),
361  mesh.nCells()
362  );
363  newToOld.note() =
364  "From cells in mesh at "
365  + runTime.timeName()
366  + " to cells in mesh at "
367  + oldTimeName;
368 
369 
370  forAll(oldToNew, oldCelli)
371  {
372  const labelList& added = oldToNew[oldCelli];
373 
374  if (added.size())
375  {
376  forAll(added, i)
377  {
378  newToOld[added[i]] = oldCelli;
379  }
380  }
381  else
382  {
383  // Unrefined cell
384  newToOld[oldCelli] = oldCelli;
385  }
386  }
387 
388  Info<< "Writing map from new to old cell to "
389  << newToOld.localObjectPath() << nl << endl;
390 
391  newToOld.write();
392 
393  printEdgeStats(mesh);
394 
395  Info<< "End\n" << endl;
396 
397  return 0;
398 }
399 
400 
401 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual Ostream & write(const char)=0
Write character.
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
void off()
Switch the function objects off.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:312
Does multiple pass refinement to refine cells in multiple directions.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Set of directions for each cell in the mesh. Either uniform and size=1 or one set of directions per c...
Definition: directions.H:64
label nCells() const
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
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:636
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:825
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:1133
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1131
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:814
string & note()
Return non-constant access to the optional note.
Definition: IOobject.H:321
static PackedBoolList getMasterEdges(const polyMesh &)
Get per edge whether it is uncoupled or a master of a.
Definition: syncTools.C:109
scalar y
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:802
dynamicFvMesh & mesh
const cellShapeList & cells
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
A class for handling words, derived from string.
Definition: word.H:59
vector vec(const pointField &) const
Return the vector (end - start)
Definition: edgeI.H:175
static const char nl
Definition: Ostream.H:260
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:32
label nEdges() const
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
A bit-packed bool list.
const functionObjectList & functionObjects() const
Return the list of function objects.
Definition: Time.H:410
A collection of cell labels.
Definition: cellSet.H:48
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
IOobject systemDictIO(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion)
Definition: systemDict.C:33
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
virtual bool write(const bool write=true) const
Write using setting from DB.
Foam::argList args(argc, argv)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
const word dictName("noiseDict")
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844