refinementLevel.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-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
27 Application
28  refinementLevel
29 
30 Group
31  grpMeshAdvancedUtilities
32 
33 Description
34  Attempt to determine the refinement levels of a refined cartesian mesh.
35  Run BEFORE snapping.
36 
37  Writes
38  - volScalarField 'refinementLevel' with current refinement level.
39  - cellSet 'refCells' which are the cells that need to be refined to satisfy
40  2:1 refinement.
41 
42  Works by dividing cells into volume bins.
43 
44 \*---------------------------------------------------------------------------*/
45 
46 #include "argList.H"
47 #include "Time.H"
48 #include "polyMesh.H"
49 #include "cellSet.H"
50 #include "SortList.H"
51 #include "labelIOList.H"
52 #include "fvMesh.H"
53 #include "volFields.H"
54 
55 using namespace Foam;
56 
57 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 
59 // Return true if any cells had to be split to keep a difference between
60 // neighbouring refinement levels < limitDiff. Puts cells into refCells and
61 // update refLevel to account for refinement.
62 bool limitRefinementLevel
63 (
64  const primitiveMesh& mesh,
65  labelList& refLevel,
66  cellSet& refCells
67 )
68 {
69  const labelListList& cellCells = mesh.cellCells();
70 
71  label oldNCells = refCells.size();
72 
73  forAll(cellCells, celli)
74  {
75  const labelList& cCells = cellCells[celli];
76 
77  forAll(cCells, i)
78  {
79  if (refLevel[cCells[i]] > (refLevel[celli]+1))
80  {
81  // Found neighbour with >=2 difference in refLevel.
82  refCells.insert(celli);
83  refLevel[celli]++;
84  break;
85  }
86  }
87  }
88 
89  if (refCells.size() > oldNCells)
90  {
91  Info<< "Added an additional " << refCells.size() - oldNCells
92  << " cells to satisfy 1:2 refinement level"
93  << endl;
94 
95  return true;
96  }
97 
98  return false;
99 }
100 
101 
102 int main(int argc, char *argv[])
103 {
105  (
106  "Attempt to determine refinement levels of a refined cartesian mesh.\n"
107  "Run BEFORE snapping!"
108  );
109 
111  (
112  "readLevel",
113  "Read level from refinementLevel file"
114  );
115 
116  #include "setRootCase.H"
117  #include "createTime.H"
118  #include "createPolyMesh.H"
119 
120  Info<< "Dividing cells into bins depending on cell volume.\nThis will"
121  << " correspond to refinement levels for a mesh with only 2x2x2"
122  << " refinement\n"
123  << "The upper range for every bin is always 1.1 times the lower range"
124  << " to allow for some truncation error."
125  << nl << endl;
126 
127  const bool readLevel = args.found("readLevel");
128 
129  const scalarField& vols = mesh.cellVolumes();
130 
131  SortList<scalar> sortedVols(vols);
132 
133  // All cell labels, sorted per bin.
135 
136  // Lower/upper limits
138 
139  // Create bin0. Have upperlimit as factor times lowerlimit.
140  bins.emplace_back();
141  limits.emplace_back(sortedVols[0], 1.1*sortedVols[0]);
142 
143  forAll(sortedVols, i)
144  {
145  if (sortedVols[i] > limits.back().max())
146  {
147  // New value outside of current bin
148  Info<< "Collected " << bins.back() << " elements in bin "
149  << limits.back().min() << " .. "
150  << limits.back().max() << endl;
151 
152  // Create new bin.
153  bins.emplace_back();
154  limits.emplace_back(sortedVols[i], 1.1 * sortedVols[i]);
155 
156  Info<< "Creating new bin "
157  << limits.back().min() << " .. "
158  << limits.back().max() << endl;
159  }
160 
161  // Add to current bin.
162  bins.back().push_back(sortedVols.indices()[i]);
163  }
164  Info<< endl;
165 
166  //
167  // Write to cellSets.
168  //
169 
170  Info<< "Volume bins:" << nl;
171  forAll(bins, bini)
172  {
173  const auto& bin = bins[bini];
174 
175  cellSet cells(mesh, "vol" + Foam::name(bini), bin.size());
176  cells.insert(bin);
177 
178  Info<< " " << limits[bini].min() << " .. " << limits[bini].max()
179  << " : writing " << bin.size() << " cells to cellSet "
180  << cells.name() << endl;
181 
182  cells.write();
183  }
184 
185 
186 
187  //
188  // Convert bins into refinement level.
189  //
190 
191 
192  // Construct fvMesh to be able to construct volScalarField
193 
194  fvMesh fMesh
195  (
196  IOobject
197  (
199  runTime.timeName(),
200  runTime,
203  ),
204  pointField(mesh.points()), // Could we safely re-use the data?
205  faceList(mesh.faces()),
206  cellList(mesh.cells())
207  );
208 
209  // Add the boundary patches
211 
212  polyPatchList newPatches(patches.size());
213 
214  forAll(newPatches, patchi)
215  {
216  newPatches.set
217  (
218  patchi,
219  patches[patchi].clone(fMesh.boundaryMesh())
220  );
221  }
222 
223  fMesh.addFvPatches(newPatches);
224 
225 
226  // Refinement level
227  IOobject refHeader
228  (
229  "refinementLevel",
230  runTime.timeName(),
232  runTime
233  );
234 
235  if (!readLevel && refHeader.typeHeaderOk<labelIOList>(true))
236  {
238  << "Detected " << refHeader.name() << " file in "
239  << polyMesh::defaultRegion << " directory. Please remove to"
240  << " recreate it or use the -readLevel option to use it"
241  << endl;
242  return 1;
243  }
244 
245 
246  labelIOList refLevel
247  (
248  IOobject
249  (
250  "refinementLevel",
251  runTime.timeName(),
252  mesh,
255  ),
257  );
258 
259  if (readLevel)
260  {
261  refLevel = labelIOList(refHeader);
262  }
263 
264  // Construct volScalarField with same info for post processing
265  volScalarField postRefLevel
266  (
267  IOobject
268  (
269  "refinementLevel",
270  runTime.timeName(),
271  mesh,
274  ),
275  fMesh,
277  );
278 
279  // Set cell values
280  forAll(bins, bini)
281  {
282  const auto& bin = bins[bini];
283 
284  forAll(bin, i)
285  {
286  refLevel[bin[i]] = bins.size() - bini - 1;
287  postRefLevel[bin[i]] = refLevel[bin[i]];
288  }
289  }
290 
291  volScalarField::Boundary& postRefLevelBf =
292  postRefLevel.boundaryFieldRef();
293 
294  // For volScalarField: set boundary values to same as cell.
295  // Note: could also put
296  // zeroGradient b.c. on postRefLevel and do evaluate.
297  forAll(postRefLevel.boundaryField(), patchi)
298  {
299  const polyPatch& pp = patches[patchi];
300 
301  fvPatchScalarField& bField = postRefLevelBf[patchi];
302 
303  Info<< "Setting field for patch "<< endl;
304 
305  forAll(bField, facei)
306  {
307  label own = mesh.faceOwner()[pp.start() + facei];
308 
309  bField[facei] = postRefLevel[own];
310  }
311  }
312 
313  Info<< "Determined current refinement level and writing to "
314  << postRefLevel.name() << " (as volScalarField; for post processing)"
315  << nl
316  << polyMesh::defaultRegion/refLevel.name()
317  << " (as labelIOList; for meshing)" << nl
318  << endl;
319 
320  refLevel.write();
321  postRefLevel.write();
322 
323 
324  // Find out cells to refine to keep to 2:1 refinement level restriction
325 
326  // Cells to refine
327  cellSet refCells(mesh, "refCells", 100);
328 
329  while
330  (
331  limitRefinementLevel
332  (
333  mesh,
334  refLevel, // current refinement level
335  refCells // cells to refine
336  )
337  )
338  {}
339 
340  if (refCells.size())
341  {
342  Info<< "Collected " << refCells.size() << " cells that need to be"
343  << " refined to get closer to overall 2:1 refinement level limit"
344  << nl
345  << "Written cells to be refined to cellSet " << refCells.name()
346  << nl << endl;
347 
348  refCells.write();
349 
350  Info<< "After refinement this tool can be run again to see if the 2:1"
351  << " limit is observed all over the mesh" << nl << endl;
352  }
353  else
354  {
355  Info<< "All cells in the mesh observe the 2:1 refinement level limit"
356  << nl << endl;
357  }
358 
359  Info<< "\nEnd\n" << endl;
360  return 0;
361 }
362 
363 
364 // ************************************************************************* //
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:462
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
Required Variables.
List< cell > cellList
List of cell.
Definition: cellListFwd.H:39
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:75
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:374
const cellList & cells() const
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:232
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
label size() const noexcept
The number of elements in table.
Definition: HashTable.H:358
const T & min() const noexcept
The min value (first)
Definition: MinMaxI.H:107
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
An indirect list with addressing based on sorting. The list is sorted upon construction or when expli...
Definition: SortList.H:50
IOList< label > labelIOList
IO for a List of label.
Definition: labelIOList.H:32
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const cellShapeList & cells
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:51
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:609
T & emplace_back(Args &&... args)
Construct an element at the end of the list, return reference to the new list element.
Definition: DynamicListI.H:538
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:406
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1116
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label nCells() const noexcept
Number of mesh cells.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
A collection of cell labels.
Definition: cellSet.H:47
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const polyBoundaryMesh & patches
Nothing to be read.
Automatically write from objectRegistry::writeObject()
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const T & max() const noexcept
The max value (second)
Definition: MinMaxI.H:121
T & back()
Access last element of the list, position [size()-1].
Definition: UListI.H:244
messageStream Info
Information stream (stdout output on master, null elsewhere)
List< label > labelList
A List of labels.
Definition: List.H:62
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
Foam::argList args(argc, argv)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
const labelListList & cellCells() const
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
const scalarField & cellVolumes() const
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127