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"
173  runTime.functionObjects().off();
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.headerOk())
197  {
198  Info<< "Refining according to "
199  << dictIO.path(typeGlobalFile<IOdictionary>()) << nl << endl;
200  refineDict = IOdictionary(dictIO);
201  }
202  else
203  {
205  << "Cannot open specified refinement dictionary "
206  << dictIO.path(typeGlobalFile<IOdictionary>())
207  << exit(FatalError);
208  }
209  }
210  else if (!refineAllCells)
211  {
212  if (dictIO.headerOk())
213  {
214  Info<< "Refining according to "
215  << dictIO.path(typeGlobalFile<IOdictionary>()) << nl << endl;
216  refineDict = IOdictionary(dictIO);
217  }
218  else
219  {
220  Info<< "Refinement dictionary "
221  << dictIO.path(typeGlobalFile<IOdictionary>()) << " not found"
222  << nl << endl;
223  }
224  }
225 
226  if (refineDict.size())
227  {
228  const word setName(refineDict.lookup("set"));
229 
230  cellSet cells(mesh, setName);
231 
232  Info<< "Read " << returnReduce(cells.size(), sumOp<label>())
233  << " cells from cellSet "
234  << cells.instance()/cells.local()/cells.name()
235  << endl << endl;
236 
237  refCells = cells.toc();
238  }
239  else
240  {
241  Info<< "Refining all cells" << nl << endl;
242 
243  // Select all cells
244  refCells = identity(mesh.nCells());
245 
246  if (mesh.nGeometricD() == 3)
247  {
248  Info<< "3D case; refining all directions" << nl << endl;
249 
250  wordList directions(3);
251  directions[0] = "e1";
252  directions[1] = "e2";
253  directions[2] = "e3";
254  refineDict.add("directions", directions);
255 
256  // Use hex cutter
257  refineDict.add("useHexTopology", "true");
258  }
259  else
260  {
261  const Vector<label> dirs(mesh.geometricD());
262 
263  wordList directions(2);
264 
265  if (dirs.x() == -1)
266  {
267  Info<< "2D case; refining in directions y,z\n" << endl;
268  directions[0] = "e2";
269  directions[1] = "e3";
270  }
271  else if (dirs.y() == -1)
272  {
273  Info<< "2D case; refining in directions x,z\n" << endl;
274  directions[0] = "e1";
275  directions[1] = "e3";
276  }
277  else
278  {
279  Info<< "2D case; refining in directions x,y\n" << endl;
280  directions[0] = "e1";
281  directions[1] = "e2";
282  }
283 
284  refineDict.add("directions", directions);
285 
286  // Use standard cutter
287  refineDict.add("useHexTopology", "false");
288  }
289 
290  refineDict.add("coordinateSystem", "global");
291 
292  dictionary coeffsDict;
293  coeffsDict.add("e1", vector(1, 0, 0));
294  coeffsDict.add("e2", vector(0, 1, 0));
295  refineDict.add("globalCoeffs", coeffsDict);
296 
297  refineDict.add("geometricCut", "false");
298  refineDict.add("writeMesh", "false");
299  }
300 
301 
302  string oldTimeName(runTime.timeName());
303 
304  if (!overwrite)
305  {
306  runTime++;
307  }
308 
309 
310  // Multi-directional refinement (does multiple iterations)
311  multiDirRefinement multiRef(mesh, refCells, refineDict);
312 
313 
314  // Write resulting mesh
315  if (overwrite)
316  {
317  mesh.setInstance(oldInstance);
318  }
319  mesh.write();
320 
321 
322  // Get list of cell splits.
323  // (is for every cell in old mesh the cells they have been split into)
324  const labelListList& oldToNew = multiRef.addedCells();
325 
326 
327  // Create cellSet with added cells for easy inspection
328  cellSet newCells(mesh, "refinedCells", refCells.size());
329 
330  forAll(oldToNew, oldCelli)
331  {
332  const labelList& added = oldToNew[oldCelli];
333 
334  forAll(added, i)
335  {
336  newCells.insert(added[i]);
337  }
338  }
339 
340  Info<< "Writing refined cells ("
341  << returnReduce(newCells.size(), sumOp<label>())
342  << ") to cellSet "
343  << newCells.instance()/newCells.local()/newCells.name()
344  << endl << endl;
345 
346  newCells.write();
347 
348 
349 
350  //
351  // Invert cell split to construct map from new to old
352  //
353 
354  labelIOList newToOld
355  (
356  IOobject
357  (
358  "cellMap",
359  runTime.timeName(),
361  mesh,
364  ),
365  mesh.nCells()
366  );
367  newToOld.note() =
368  "From cells in mesh at "
369  + runTime.timeName()
370  + " to cells in mesh at "
371  + oldTimeName;
372 
373 
374  forAll(oldToNew, oldCelli)
375  {
376  const labelList& added = oldToNew[oldCelli];
377 
378  if (added.size())
379  {
380  forAll(added, i)
381  {
382  newToOld[added[i]] = oldCelli;
383  }
384  }
385  else
386  {
387  // Unrefined cell
388  newToOld[oldCelli] = oldCelli;
389  }
390  }
391 
392  Info<< "Writing map from new to old cell to "
393  << newToOld.relativeObjectPath() << nl << endl;
394 
395  newToOld.write();
396 
397  printEdgeStats(mesh);
398 
399  Info<< "End\n" << endl;
400 
401  return 0;
402 }
403 
404 
405 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
virtual Ostream & write(const char)=0
Write character.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
Templated form of IOobject providing type information for file reading and header type checking...
Definition: IOobject.H:537
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:328
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
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
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
fvMesh & mesh
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:905
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:1153
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:1211
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:894
string & note()
Return non-constant access to the optional note.
Definition: IOobject.H:333
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:882
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
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
vector vec(const pointField &) const
Return the vector (end - start)
Definition: edgeI.H:175
static const char nl
Definition: Ostream.H:260
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:77
label nEdges() const
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
A bit-packed bool list.
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:76
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:98
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:864