Loading...
Searching...
No Matches
refineMesh.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2018-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Application
28 refineMesh
29
30Group
31 grpMeshManipulationUtilities
32
33Description
34 Utility to refine cells in multiple directions.
35
36 Command-line option handling:
37 - If -all specified or no refineMeshDict exists or, refine all cells
38 - If -dict <file> specified refine according to <file>
39 - If refineMeshDict exists refine according to refineMeshDict
40
41 When the refinement or all cells is selected apply 3D refinement for 3D
42 cases and 2D refinement for 2D cases.
43
44\*---------------------------------------------------------------------------*/
45
46#include "argList.H"
47#include "polyMesh.H"
48#include "Time.H"
49#include "cellSet.H"
50#include "multiDirRefinement.H"
51#include "labelIOList.H"
52#include "IOdictionary.H"
53#include "syncTools.H"
54
55using namespace Foam;
56
57// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58
59// Max cos angle for edges to be considered aligned with axis.
60static const scalar edgeTol = 1e-3;
61
62
63// Print edge statistics on mesh.
64void printEdgeStats(const polyMesh& mesh)
65{
66 label nAny(0);
67 label nX(0);
68 label nY(0);
69 label nZ(0);
70
71 scalarMinMax limitsAny(GREAT, -GREAT);
72 scalarMinMax limitsX(limitsAny);
73 scalarMinMax limitsY(limitsAny);
74 scalarMinMax limitsZ(limitsAny);
75
76 bitSet isMasterEdge(syncTools::getMasterEdges(mesh));
77
78 const edgeList& edges = mesh.edges();
79
80 for (const label edgei : isMasterEdge)
81 {
82 const edge& e = edges[edgei];
83
84 vector eVec(e.vec(mesh.points()));
85
86 scalar eMag = mag(eVec);
87
88 eVec /= eMag;
89
90 if (mag(eVec.x()) > 1-edgeTol)
91 {
92 limitsX.add(eMag);
93 nX++;
94 }
95 else if (mag(eVec.y()) > 1-edgeTol)
96 {
97 limitsY.add(eMag);
98 nY++;
99 }
100 else if (mag(eVec.z()) > 1-edgeTol)
101 {
102 limitsZ.add(eMag);
103 nZ++;
104 }
105 else
106 {
107 limitsAny.add(eMag);
108 nAny++;
109 }
110 }
111
112 reduce(nX, sumOp<label>());
113 reduce(nY, sumOp<label>());
114 reduce(nZ, sumOp<label>());
115 reduce(nAny, sumOp<label>());
116
117 reduce(limitsX, sumOp<scalarMinMax>());
118 reduce(limitsY, sumOp<scalarMinMax>());
119 reduce(limitsZ, sumOp<scalarMinMax>());
120 reduce(limitsAny, sumOp<scalarMinMax>());
121
122 Info<< "Mesh edge statistics:" << nl
123 << " x aligned : number:" << nX
124 << "\tminLen:" << limitsX.min() << "\tmaxLen:" << limitsX.max() << nl
125 << " y aligned : number:" << nY
126 << "\tminLen:" << limitsY.min() << "\tmaxLen:" << limitsY.max() << nl
127 << " z aligned : number:" << nZ
128 << "\tminLen:" << limitsZ.min() << "\tmaxLen:" << limitsZ.max() << nl
129 << " other : number:" << nAny
130 << "\tminLen:" << limitsAny.min()
131 << "\tmaxLen:" << limitsAny.max() << nl << endl;
132}
133
134
135int main(int argc, char *argv[])
136{
138 (
139 "Refine cells in multiple directions"
140 );
141
142 #include "addOverwriteOption.H"
143 #include "addRegionOption.H"
144
145 argList::addOption("dict", "file", "Alternative refineMeshDict");
146
148 (
149 "all",
150 "Refine all cells"
151 );
152
153 argList::noFunctionObjects(); // Never use function objects
154
155 #include "setRootCase.H"
156 #include "createTime.H"
157 #include "createNamedPolyMesh.H"
158
159 const word oldInstance = mesh.pointsInstance();
160
161 printEdgeStats(mesh);
162
163 //
164 // Read/construct control dictionary
165 //
166
167 const bool refineAllCells = args.found("all");
168 const bool overwrite = args.found("overwrite");
169
170 // List of cells to refine
171 labelList refCells;
172
173 // Dictionary to control refinement
174 dictionary refineDict;
175
176 // The -all option has precedence over -dict, or anything else
177 if (!refineAllCells)
178 {
179 const word dictName("refineMeshDict");
180
181 // Obtain dictPath here for messages
182 fileName dictPath = args.getOrDefault<fileName>("dict", "");
183
185 (
187 (
188 dictName,
189 runTime.system(),
190 mesh,
192 ),
194 );
195
196 // The reported dictPath will not be full resolved for the output
197 // (it will just be the -dict value) but this is purely cosmetic.
198
199 if (dictIO.typeHeaderOk<IOdictionary>(true))
200 {
201 refineDict = IOdictionary(dictIO);
202
203 Info<< "Refining according to ";
204
205 if (dictPath.empty())
206 {
207 Info<< dictName;
208 }
209 else
210 {
211 Info<< dictPath;
212 }
213 Info<< nl << endl;
214 }
215 else if (dictPath.empty())
216 {
217 Info<< "Refinement dictionary " << dictName << " not found" << nl;
218 }
219 else
220 {
222 << "Cannot open specified refinement dictionary "
223 << dictPath
224 << exit(FatalError);
225 }
226 }
227
228 if (refineDict.size())
229 {
230 const word setName(refineDict.get<word>("set"));
231
232 cellSet cells(mesh, setName);
233
234 Info<< "Read " << returnReduce(cells.size(), sumOp<label>())
235 << " cells from cellSet "
236 << cells.instance()/cells.local()/cells.name()
237 << endl << endl;
238
239 refCells = cells.toc();
240 }
241 else
242 {
243 Info<< "Refining all cells" << nl << endl;
244
245 // Select all cells
246 refCells = identity(mesh.nCells());
247
248 if (mesh.nGeometricD() == 3)
249 {
250 Info<< "3D case; refining all directions" << nl << endl;
251
253 directions[0] = "tan1";
254 directions[1] = "tan2";
255 directions[2] = "normal";
256 refineDict.add("directions", directions);
257
258 // Use hex cutter
259 refineDict.add("useHexTopology", "true");
260 }
261 else
262 {
263 const Vector<label> dirs(mesh.geometricD());
264
266
267 if (dirs.x() == -1)
268 {
269 Info<< "2D case; refining in directions y,z\n" << endl;
270 directions[0] = "tan2";
271 directions[1] = "normal";
272 }
273 else if (dirs.y() == -1)
274 {
275 Info<< "2D case; refining in directions x,z\n" << endl;
276 directions[0] = "tan1";
277 directions[1] = "normal";
278 }
279 else
280 {
281 Info<< "2D case; refining in directions x,y\n" << endl;
282 directions[0] = "tan1";
283 directions[1] = "tan2";
284 }
285
286 refineDict.add("directions", directions);
287
288 // Use standard cutter
289 refineDict.add("useHexTopology", "false");
290 }
291
292 refineDict.add("coordinateSystem", "global");
293
294 dictionary coeffsDict;
295 coeffsDict.add("tan1", vector(1, 0, 0));
296 coeffsDict.add("tan2", vector(0, 1, 0));
297 refineDict.add("globalCoeffs", coeffsDict);
298
299 refineDict.add("geometricCut", "false");
300 refineDict.add("writeMesh", "false");
301 }
302
303
304 string oldTimeName(runTime.timeName());
305
306 if (!overwrite)
307 {
308 ++runTime;
309 }
310
311
312 // Multi-directional refinement (does multiple iterations)
313 multiDirRefinement multiRef(mesh, refCells, refineDict);
314
315
316 // Write resulting mesh
317 if (overwrite)
318 {
319 mesh.setInstance(oldInstance);
320 }
321 mesh.write();
322
323
324 // Get list of cell splits.
325 // (is for every cell in old mesh the cells they have been split into)
326 const labelListList& oldToNew = multiRef.addedCells();
327
328
329 // Create cellSet with added cells for easy inspection
330 cellSet newCells(mesh, "refinedCells", refCells.size());
331
332 for (const labelList& added : oldToNew)
333 {
334 newCells.insert(added);
335 }
336
337 Info<< "Writing refined cells ("
338 << returnReduce(newCells.size(), sumOp<label>())
339 << ") to cellSet "
340 << newCells.instance()/newCells.local()/newCells.name()
341 << endl << endl;
342
343 newCells.write();
344
345
346 //
347 // Invert cell split to construct map from new to old
348 //
349
350 labelIOList newToOld
351 (
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 for (const label celli : added)
377 {
378 newToOld[celli] = 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.objectPath() << nl << endl;
390
391 newToOld.write();
392
393 printEdgeStats(mesh);
394
395 Info<< "End\n" << endl;
396
397 return 0;
398}
399
400
401// ************************************************************************* //
Required Classes.
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static IOobject selectIO(const IOobject &io, const fileName &altFile, const word &ioName="")
Return the IOobject, but also consider an alternative file name.
Definition IOobject.C:256
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition Vector.H:61
static void noFunctionObjects(bool addWithOption=false)
Remove '-noFunctionObjects' option and ignore any occurrences.
Definition argList.C:562
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition argList.C:389
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition argList.C:400
static void addNote(const string &note)
Add extra notes for the usage information.
Definition argList.C:477
A collection of cell labels.
Definition cellSet.H:50
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition dictionary.C:625
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 vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
A class for handling file names.
Definition fileName.H:75
Does multiple pass refinement to refine cells in multiple directions.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh").
Definition polyMesh.H:411
static bitSet getMasterEdges(const polyMesh &mesh)
Get per edge whether it is uncoupled or a master of a coupled set of edges.
Definition syncTools.C:90
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
dynamicFvMesh & mesh
engineTime & runTime
Required Classes.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const word dictName("faMeshDefinition")
fileName dictPath(runTime.system()/regionDir/faMesh::prefix()/areaRegionDir/dictName)
const cellShapeList & cells
return returnReduce(nRefine-oldNRefine, sumOp< label >())
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
List< word > wordList
List of word.
Definition fileName.H:60
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
IOList< label > labelIOList
IO for a List of label.
Definition labelIOList.H:32
messageStream Info
Information stream (stdout output on master, null elsewhere).
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition MinMax.H:97
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
IOobject dictIO
Foam::argList args(argc, argv)
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299