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-2018 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 
49 using namespace Foam;
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 // Max cos angle for edges to be considered aligned with axis.
54 static const scalar edgeTol = 1e-3;
55 
56 
57 // Print edge statistics on mesh.
58 void printEdgeStats(const polyMesh& mesh)
59 {
60  label nX = 0;
61  label nY = 0;
62  label nZ = 0;
63 
64  scalar minX = great;
65  scalar maxX = -great;
66  static const vector x(1, 0, 0);
67 
68  scalar minY = great;
69  scalar maxY = -great;
70  static const vector y(0, 1, 0);
71 
72  scalar minZ = great;
73  scalar maxZ = -great;
74  static const vector z(0, 0, 1);
75 
76  scalar minOther = great;
77  scalar maxOther = -great;
78 
79  PackedBoolList isMasterEdge(syncTools::getMasterEdges(mesh));
80 
81  const edgeList& edges = mesh.edges();
82 
83  forAll(edges, edgeI)
84  {
85  if (isMasterEdge[edgeI])
86  {
87  const edge& e = edges[edgeI];
88 
89  vector eVec(e.vec(mesh.points()));
90 
91  scalar eMag = mag(eVec);
92 
93  eVec /= eMag;
94 
95  if (mag(eVec & x) > 1-edgeTol)
96  {
97  minX = min(minX, eMag);
98  maxX = max(maxX, eMag);
99  nX++;
100  }
101  else if (mag(eVec & y) > 1-edgeTol)
102  {
103  minY = min(minY, eMag);
104  maxY = max(maxY, eMag);
105  nY++;
106  }
107  else if (mag(eVec & z) > 1-edgeTol)
108  {
109  minZ = min(minZ, eMag);
110  maxZ = max(maxZ, eMag);
111  nZ++;
112  }
113  else
114  {
115  minOther = min(minOther, eMag);
116  maxOther = max(maxOther, eMag);
117  }
118  }
119  }
120 
121  label nEdges = mesh.nEdges();
122  reduce(nEdges, sumOp<label>());
123  reduce(nX, sumOp<label>());
124  reduce(nY, sumOp<label>());
125  reduce(nZ, sumOp<label>());
126 
127  reduce(minX, minOp<scalar>());
128  reduce(maxX, maxOp<scalar>());
129 
130  reduce(minY, minOp<scalar>());
131  reduce(maxY, maxOp<scalar>());
132 
133  reduce(minZ, minOp<scalar>());
134  reduce(maxZ, maxOp<scalar>());
135 
136  reduce(minOther, minOp<scalar>());
137  reduce(maxOther, maxOp<scalar>());
138 
139 
140  Info<< "Mesh edge statistics:" << nl
141  << " x aligned : number:" << nX << "\tminLen:" << minX
142  << "\tmaxLen:" << maxX << nl
143  << " y aligned : number:" << nY << "\tminLen:" << minY
144  << "\tmaxLen:" << maxY << nl
145  << " z aligned : number:" << nZ << "\tminLen:" << minZ
146  << "\tmaxLen:" << maxZ << nl
147  << " other : number:" << nEdges - nX - nY - nZ
148  << "\tminLen:" << minOther
149  << "\tmaxLen:" << maxOther << nl << endl;
150 }
151 
152 
153 int main(int argc, char *argv[])
154 {
156  (
157  "refine cells in multiple directions"
158  );
159 
160  #include "addOverwriteOption.H"
161  #include "addRegionOption.H"
162  #include "addDictOption.H"
163 
165  (
166  "all",
167  "Refine all cells"
168  );
169 
170  #include "setRootCase.H"
171  #include "createTime.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  dictionary refineDict;
191  const word dictName("refineMeshDict");
192 
193  if (readDict)
194  {
195  fileName dictPath = args["dict"];
196  if (isDir(dictPath))
197  {
198  dictPath = dictPath/dictName;
199  }
200 
202  (
203  dictPath,
204  mesh,
206  );
207 
208  if (!dictIO.typeHeaderOk<IOdictionary>(true))
209  {
211  << "Cannot open specified refinement dictionary "
212  << dictPath
213  << exit(FatalError);
214  }
215 
216  Info<< "Refining according to " << dictPath << nl << endl;
217 
218  refineDict = IOdictionary(dictIO);
219  }
220  else if (!refineAllCells)
221  {
223  (
224  dictName,
225  runTime.system(),
226  mesh,
228  );
229 
230  if (dictIO.typeHeaderOk<IOdictionary>(true))
231  {
232  Info<< "Refining according to " << dictName << nl << endl;
233 
234  refineDict = IOdictionary(dictIO);
235  }
236  else
237  {
238  Info<< "Refinement dictionary " << dictName << " not found" << endl;
239  }
240  }
241 
242  if (refineDict.size())
243  {
244  const word setName(refineDict.lookup("set"));
245 
246  cellSet cells(mesh, setName);
247 
248  Info<< "Read " << returnReduce(cells.size(), sumOp<label>())
249  << " cells from cellSet "
250  << cells.instance()/cells.local()/cells.name()
251  << endl << endl;
252 
253  refCells = cells.toc();
254  }
255  else
256  {
257  Info<< "Refining all cells" << nl << endl;
258 
259  // Select all cells
260  refCells = identity(mesh.nCells());
261 
262  if (mesh.nGeometricD() == 3)
263  {
264  Info<< "3D case; refining all directions" << nl << endl;
265 
266  wordList directions(3);
267  directions[0] = "tan1";
268  directions[1] = "tan2";
269  directions[2] = "normal";
270  refineDict.add("directions", directions);
271 
272  // Use hex cutter
273  refineDict.add("useHexTopology", "true");
274  }
275  else
276  {
277  const Vector<label> dirs(mesh.geometricD());
278 
279  wordList directions(2);
280 
281  if (dirs.x() == -1)
282  {
283  Info<< "2D case; refining in directions y,z\n" << endl;
284  directions[0] = "tan2";
285  directions[1] = "normal";
286  }
287  else if (dirs.y() == -1)
288  {
289  Info<< "2D case; refining in directions x,z\n" << endl;
290  directions[0] = "tan1";
291  directions[1] = "normal";
292  }
293  else
294  {
295  Info<< "2D case; refining in directions x,y\n" << endl;
296  directions[0] = "tan1";
297  directions[1] = "tan2";
298  }
299 
300  refineDict.add("directions", directions);
301 
302  // Use standard cutter
303  refineDict.add("useHexTopology", "false");
304  }
305 
306  refineDict.add("coordinateSystem", "global");
307 
308  dictionary coeffsDict;
309  coeffsDict.add("tan1", vector(1, 0, 0));
310  coeffsDict.add("tan2", vector(0, 1, 0));
311  refineDict.add("globalCoeffs", coeffsDict);
312 
313  refineDict.add("geometricCut", "false");
314  refineDict.add("writeMesh", "false");
315  }
316 
317 
318  string oldTimeName(runTime.timeName());
319 
320  if (!overwrite)
321  {
322  runTime++;
323  }
324 
325 
326  // Multi-directional refinement (does multiple iterations)
327  multiDirRefinement multiRef(mesh, refCells, refineDict);
328 
329 
330  // Write resulting mesh
331  if (overwrite)
332  {
333  mesh.setInstance(oldInstance);
334  }
335  mesh.write();
336 
337 
338  // Get list of cell splits.
339  // (is for every cell in old mesh the cells they have been split into)
340  const labelListList& oldToNew = multiRef.addedCells();
341 
342 
343  // Create cellSet with added cells for easy inspection
344  cellSet newCells(mesh, "refinedCells", refCells.size());
345 
346  forAll(oldToNew, oldCelli)
347  {
348  const labelList& added = oldToNew[oldCelli];
349 
350  forAll(added, i)
351  {
352  newCells.insert(added[i]);
353  }
354  }
355 
356  Info<< "Writing refined cells ("
357  << returnReduce(newCells.size(), sumOp<label>())
358  << ") to cellSet "
359  << newCells.instance()/newCells.local()/newCells.name()
360  << endl << endl;
361 
362  newCells.write();
363 
364 
365 
366  //
367  // Invert cell split to construct map from new to old
368  //
369 
370  labelIOList newToOld
371  (
372  IOobject
373  (
374  "cellMap",
375  runTime.timeName(),
377  mesh,
380  ),
381  mesh.nCells()
382  );
383  newToOld.note() =
384  "From cells in mesh at "
385  + runTime.timeName()
386  + " to cells in mesh at "
387  + oldTimeName;
388 
389 
390  forAll(oldToNew, oldCelli)
391  {
392  const labelList& added = oldToNew[oldCelli];
393 
394  if (added.size())
395  {
396  forAll(added, i)
397  {
398  newToOld[added[i]] = oldCelli;
399  }
400  }
401  else
402  {
403  // Unrefined cell
404  newToOld[oldCelli] = oldCelli;
405  }
406  }
407 
408  Info<< "Writing map from new to old cell to "
409  << newToOld.objectPath() << nl << endl;
410 
411  newToOld.write();
412 
413  printEdgeStats(mesh);
414 
415  Info<< "End\n" << endl;
416 
417  return 0;
418 }
419 
420 
421 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
bool typeHeaderOk(const bool checkType=true)
Read header (uses typeFilePath to find file) and check header.
A class for handling file names.
Definition: fileName.H:79
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:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
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:163
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:256
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
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:626
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:821
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
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:313
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
bool isDir(const fileName &, const bool followLink=true)
Does the name exist as a directory in the file system?
Definition: POSIX.C:539
A class for handling words, derived from string.
Definition: word.H:59
const word dictName("particleTrackDict")
const word & system() const
Return system name.
Definition: TimePaths.H:114
vector vec(const pointField &) const
Return the vector (end - start)
Definition: edgeI.H:175
static const char nl
Definition: Ostream.H:265
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)
IOobject dictIO(dictName, runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE)
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 > &)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
virtual Ostream & write(const token &)=0
Write next token to stream.
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:117
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:158
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
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583