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 | Copyright (C) 2011-2016 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  refinementLevel
26 
27 Description
28  Tries to figure out what the refinement level is on refined cartesian
29  meshes. Run BEFORE snapping.
30 
31  Writes
32  - volScalarField 'refinementLevel' with current refinement level.
33  - cellSet 'refCells' which are the cells that need to be refined to satisfy
34  2:1 refinement.
35 
36  Works by dividing cells into volume bins.
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #include "argList.H"
41 #include "Time.H"
42 #include "polyMesh.H"
43 #include "cellSet.H"
44 #include "SortableList.H"
45 #include "labelIOList.H"
46 #include "fvMesh.H"
47 #include "volFields.H"
48 
49 using namespace Foam;
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 // Return true if any cells had to be split to keep a difference between
54 // neighbouring refinement levels < limitDiff. Puts cells into refCells and
55 // update refLevel to account for refinement.
56 bool limitRefinementLevel
57 (
58  const primitiveMesh& mesh,
59  labelList& refLevel,
60  cellSet& refCells
61 )
62 {
63  const labelListList& cellCells = mesh.cellCells();
64 
65  label oldNCells = refCells.size();
66 
67  forAll(cellCells, celli)
68  {
69  const labelList& cCells = cellCells[celli];
70 
71  forAll(cCells, i)
72  {
73  if (refLevel[cCells[i]] > (refLevel[celli]+1))
74  {
75  // Found neighbour with >=2 difference in refLevel.
76  refCells.insert(celli);
77  refLevel[celli]++;
78  break;
79  }
80  }
81  }
82 
83  if (refCells.size() > oldNCells)
84  {
85  Info<< "Added an additional " << refCells.size() - oldNCells
86  << " cells to satisfy 1:2 refinement level"
87  << endl;
88 
89  return true;
90  }
91  else
92  {
93  return false;
94  }
95 }
96 
97 
98 
99 int main(int argc, char *argv[])
100 {
102  (
103  "readLevel",
104  "read level from refinementLevel file"
105  );
106 
107  #include "setRootCase.H"
108  #include "createTime.H"
109  #include "createPolyMesh.H"
110 
111  Info<< "Dividing cells into bins depending on cell volume.\nThis will"
112  << " correspond to refinement levels for a mesh with only 2x2x2"
113  << " refinement\n"
114  << "The upper range for every bin is always 1.1 times the lower range"
115  << " to allow for some truncation error."
116  << nl << endl;
117 
118  const bool readLevel = args.optionFound("readLevel");
119 
120  const scalarField& vols = mesh.cellVolumes();
121 
122  SortableList<scalar> sortedVols(vols);
123 
124  // All cell labels, sorted per bin.
126 
127  // Lower/upper limits
128  DynamicList<scalar> lowerLimits;
129  DynamicList<scalar> upperLimits;
130 
131  // Create bin0. Have upperlimit as factor times lowerlimit.
132  bins.append(DynamicList<label>());
133  lowerLimits.append(sortedVols[0]);
134  upperLimits.append(1.1 * lowerLimits.last());
135 
136  forAll(sortedVols, i)
137  {
138  if (sortedVols[i] > upperLimits.last())
139  {
140  // New value outside of current bin
141 
142  // Shrink old bin.
143  DynamicList<label>& bin = bins.last();
144 
145  bin.shrink();
146 
147  Info<< "Collected " << bin.size() << " elements in bin "
148  << lowerLimits.last() << " .. "
149  << upperLimits.last() << endl;
150 
151  // Create new bin.
152  bins.append(DynamicList<label>());
153  lowerLimits.append(sortedVols[i]);
154  upperLimits.append(1.1 * lowerLimits.last());
155 
156  Info<< "Creating new bin " << lowerLimits.last()
157  << " .. " << upperLimits.last()
158  << endl;
159  }
160 
161  // Append to current bin.
162  DynamicList<label>& bin = bins.last();
163 
164  bin.append(sortedVols.indices()[i]);
165  }
166  Info<< endl;
167 
168  bins.last().shrink();
169  bins.shrink();
170  lowerLimits.shrink();
171  upperLimits.shrink();
172 
173 
174  //
175  // Write to cellSets.
176  //
177 
178  Info<< "Volume bins:" << nl;
179  forAll(bins, binI)
180  {
181  const DynamicList<label>& bin = bins[binI];
182 
183  cellSet cells(mesh, "vol" + name(binI), bin.size());
184 
185  forAll(bin, i)
186  {
187  cells.insert(bin[i]);
188  }
189 
190  Info<< " " << lowerLimits[binI] << " .. " << upperLimits[binI]
191  << " : writing " << bin.size() << " cells to cellSet "
192  << cells.name() << endl;
193 
194  cells.write();
195  }
196 
197 
198 
199  //
200  // Convert bins into refinement level.
201  //
202 
203 
204  // Construct fvMesh to be able to construct volScalarField
205 
206  fvMesh fMesh
207  (
208  IOobject
209  (
211  runTime.timeName(),
212  runTime
213  ),
214  xferCopy(mesh.points()), // could we safely re-use the data?
215  xferCopy(mesh.faces()),
216  xferCopy(mesh.cells())
217  );
218 
219  // Add the boundary patches
220  const polyBoundaryMesh& patches = mesh.boundaryMesh();
221 
222  List<polyPatch*> p(patches.size());
223 
224  forAll(p, patchi)
225  {
226  p[patchi] = patches[patchi].clone(fMesh.boundaryMesh()).ptr();
227  }
228 
229  fMesh.addFvPatches(p);
230 
231 
232  // Refinement level
233  IOobject refHeader
234  (
235  "refinementLevel",
236  runTime.timeName(),
238  runTime
239  );
240 
241  if (!readLevel && refHeader.headerOk())
242  {
244  << "Detected " << refHeader.name() << " file in "
245  << polyMesh::defaultRegion << " directory. Please remove to"
246  << " recreate it or use the -readLevel option to use it"
247  << endl;
248  return 1;
249  }
250 
251 
252  labelIOList refLevel
253  (
254  IOobject
255  (
256  "refinementLevel",
257  runTime.timeName(),
258  mesh,
261  ),
262  labelList(mesh.nCells(), 0)
263  );
264 
265  if (readLevel)
266  {
267  refLevel = labelIOList(refHeader);
268  }
269 
270  // Construct volScalarField with same info for post processing
271  volScalarField postRefLevel
272  (
273  IOobject
274  (
275  "refinementLevel",
276  runTime.timeName(),
277  mesh,
280  ),
281  fMesh,
282  dimensionedScalar("zero", dimless/dimTime, 0)
283  );
284 
285  // Set cell values
286  forAll(bins, binI)
287  {
288  const DynamicList<label>& bin = bins[binI];
289 
290  forAll(bin, i)
291  {
292  refLevel[bin[i]] = bins.size() - binI - 1;
293  postRefLevel[bin[i]] = refLevel[bin[i]];
294  }
295  }
296 
297  volScalarField::Boundary& postRefLevelBf =
298  postRefLevel.boundaryFieldRef();
299 
300  // For volScalarField: set boundary values to same as cell.
301  // Note: could also put
302  // zeroGradient b.c. on postRefLevel and do evaluate.
303  forAll(postRefLevel.boundaryField(), patchi)
304  {
305  const polyPatch& pp = patches[patchi];
306 
307  fvPatchScalarField& bField = postRefLevelBf[patchi];
308 
309  Info<< "Setting field for patch "<< endl;
310 
311  forAll(bField, facei)
312  {
313  label own = mesh.faceOwner()[pp.start() + facei];
314 
315  bField[facei] = postRefLevel[own];
316  }
317  }
318 
319  Info<< "Determined current refinement level and writing to "
320  << postRefLevel.name() << " (as volScalarField; for post processing)"
321  << nl
322  << polyMesh::defaultRegion/refLevel.name()
323  << " (as labelIOList; for meshing)" << nl
324  << endl;
325 
326  refLevel.write();
327  postRefLevel.write();
328 
329 
330  // Find out cells to refine to keep to 2:1 refinement level restriction
331 
332  // Cells to refine
333  cellSet refCells(mesh, "refCells", 100);
334 
335  while
336  (
337  limitRefinementLevel
338  (
339  mesh,
340  refLevel, // current refinement level
341  refCells // cells to refine
342  )
343  )
344  {}
345 
346  if (refCells.size())
347  {
348  Info<< "Collected " << refCells.size() << " cells that need to be"
349  << " refined to get closer to overall 2:1 refinement level limit"
350  << nl
351  << "Written cells to be refined to cellSet " << refCells.name()
352  << nl << endl;
353 
354  refCells.write();
355 
356  Info<< "After refinement this tool can be run again to see if the 2:1"
357  << " limit is observed all over the mesh" << nl << endl;
358  }
359  else
360  {
361  Info<< "All cells in the mesh observe the 2:1 refinement level limit"
362  << nl << endl;
363  }
364 
365  Info<< "\nEnd\n" << endl;
366  return 0;
367 }
368 
369 
370 // ************************************************************************* //
const labelListList & cellCells() const
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
A list that is sorted upon construction or when explicitly requested with the sort() method...
Definition: List.H:68
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:74
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:306
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
virtual const pointField & points() const =0
Return mesh points.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Xfer< T > xferCopy(const T &)
Construct by copying the contents of the arg.
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
const cellList & cells() const
label nCells() const
dynamicFvMesh & mesh
const cellShapeList & cells
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:300
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
Foam::autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:239
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:292
List< label > labelList
A List of labels.
Definition: labelList.H:56
Foam::polyBoundaryMesh.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:240
static const char nl
Definition: Ostream.H:262
const scalarField & cellVolumes() const
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A collection of cell labels.
Definition: cellSet.H:48
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual const faceList & faces() const =0
Return faces.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual bool write() const
Write using setting from DB.
messageStream Info
virtual const labelList & faceOwner() const =0
Face face-owner addresing.
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:83
volScalarField & p
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
T & last()
Return the last element of the list.
Definition: UListI.H:128
Foam::argList args(argc, argv)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
const word & name() const
Return name.
Definition: IOobject.H:260
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:42
Namespace for OpenFOAM.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29