FreeFOAM The Cross-Platform CFD Toolkit
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) 1991-2010 OpenCFD Ltd.
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 Usage
39 
40  - refinementLevel [OPTIONS]
41 
42  @param -readLevel \n
43  Read reference refinementLevel file.
44 
45  @param -case <dir>\n
46  Case directory.
47 
48  @param -parallel \n
49  Run in parallel.
50 
51  @param -help \n
52  Display help message.
53 
54  @param -doc \n
55  Display Doxygen API documentation page for this application.
56 
57  @param -srcDoc \n
58  Display Doxygen source documentation page for this application.
59 
60 \*---------------------------------------------------------------------------*/
61 
62 #include <OpenFOAM/argList.H>
63 #include <OpenFOAM/Time.H>
64 #include <OpenFOAM/polyMesh.H>
65 #include <meshTools/cellSet.H>
66 #include <OpenFOAM/SortableList.H>
67 #include <OpenFOAM/labelIOList.H>
68 #include <finiteVolume/fvMesh.H>
69 #include <finiteVolume/volFields.H>
70 
71 using namespace Foam;
72 
73 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
74 
75 // Return true if any cells had to be split to keep a difference between
76 // neighbouring refinement levels < limitDiff. Puts cells into refCells and
77 // update refLevel to account for refinement.
78 bool limitRefinementLevel
79 (
80  const primitiveMesh& mesh,
81  labelList& refLevel,
82  cellSet& refCells
83 )
84 {
85  const labelListList& cellCells = mesh.cellCells();
86 
87  label oldNCells = refCells.size();
88 
89  forAll(cellCells, cellI)
90  {
91  const labelList& cCells = cellCells[cellI];
92 
93  forAll(cCells, i)
94  {
95  if (refLevel[cCells[i]] > (refLevel[cellI]+1))
96  {
97  // Found neighbour with >=2 difference in refLevel.
98  refCells.insert(cellI);
99  refLevel[cellI]++;
100  break;
101  }
102  }
103  }
104 
105  if (refCells.size() > oldNCells)
106  {
107  Info<< "Added an additional " << refCells.size() - oldNCells
108  << " cells to satisfy 1:2 refinement level"
109  << endl;
110 
111  return true;
112  }
113  else
114  {
115  return false;
116  }
117 }
118 
119 
120 // Main program:
121 
122 int main(int argc, char *argv[])
123 {
124  argList::validOptions.insert("readLevel", "");
125 
126 # include <OpenFOAM/setRootCase.H>
127 # include <OpenFOAM/createTime.H>
128 # include <OpenFOAM/createPolyMesh.H>
129 
130  Info<< "Dividing cells into bins depending on cell volume.\nThis will"
131  << " correspond to refinement levels for a mesh with only 2x2x2"
132  << " refinement\n"
133  << "The upper range for every bin is always 1.1 times the lower range"
134  << " to allow for some truncation error."
135  << nl << endl;
136 
137  bool readLevel = args.optionFound("readLevel");
138 
139  const scalarField& vols = mesh.cellVolumes();
140 
141  SortableList<scalar> sortedVols(vols);
142 
143  // All cell labels, sorted per bin.
144  // Use autoPtr to help GCC compilers (< 4.3) which suffer from core
145  // language defect 391, i.e. require a copy-ctor to be present when passing
146  // a temporary as an rvalue.
148 
149  // Lower/upper limits
150  DynamicList<scalar> lowerLimits;
151  DynamicList<scalar> upperLimits;
152 
153  // Create bin0. Have upperlimit as factor times lowerlimit.
155  lowerLimits.append(sortedVols[0]);
156  upperLimits.append(1.1*lowerLimits[lowerLimits.size()-1]);
157 
158  forAll(sortedVols, i)
159  {
160  if (sortedVols[i] > upperLimits[upperLimits.size()-1])
161  {
162  // New value outside of current bin
163 
164  // Shrink old bin.
165  DynamicList<label>& bin = bins[bins.size()-1]();
166 
167  bin.shrink();
168 
169  Info<< "Collected " << bin.size() << " elements in bin "
170  << lowerLimits[lowerLimits.size()-1] << " .. "
171  << upperLimits[upperLimits.size()-1] << endl;
172 
173  // Create new bin.
175  lowerLimits.append(sortedVols[i]);
176  upperLimits.append(1.1*lowerLimits[lowerLimits.size()-1]);
177 
178  Info<< "Creating new bin " << lowerLimits[lowerLimits.size()-1]
179  << " .. " << upperLimits[upperLimits.size()-1]
180  << endl;
181  }
182 
183  // Append to current bin.
184  DynamicList<label>& bin = bins[bins.size()-1]();
185 
186  bin.append(sortedVols.indices()[i]);
187  }
188  Info<< endl;
189 
190  bins[bins.size()-1]->shrink();
191  bins.shrink();
192  lowerLimits.shrink();
193  upperLimits.shrink();
194 
195 
196  //
197  // Write to cellSets.
198  //
199 
200  Info<< "Volume bins:" << nl;
201  forAll(bins, binI)
202  {
203  const DynamicList<label>& bin = bins[binI]();
204 
205  cellSet cells(mesh, "vol" + name(binI), bin.size());
206 
207  forAll(bin, i)
208  {
209  cells.insert(bin[i]);
210  }
211 
212  Info<< " " << lowerLimits[binI] << " .. " << upperLimits[binI]
213  << " : writing " << bin.size() << " cells to cellSet "
214  << cells.name() << endl;
215 
216  cells.write();
217  }
218 
219 
220 
221  //
222  // Convert bins into refinement level.
223  //
224 
225 
226  // Construct fvMesh to be able to construct volScalarField
227 
228  fvMesh fMesh
229  (
230  IOobject
231  (
233  runTime.timeName(),
234  runTime
235  ),
236  xferCopy(mesh.points()), // could we safely re-use the data?
237  xferCopy(mesh.faces()),
238  xferCopy(mesh.cells())
239  );
240 
241  // Add the boundary patches
242  const polyBoundaryMesh& patches = mesh.boundaryMesh();
243 
244  List<polyPatch*> p(patches.size());
245 
246  forAll (p, patchI)
247  {
248  p[patchI] = patches[patchI].clone(fMesh.boundaryMesh()).ptr();
249  }
250 
251  fMesh.addFvPatches(p);
252 
253 
254  // Refinement level
255  IOobject refHeader
256  (
257  "refinementLevel",
258  runTime.timeName(),
260  runTime
261  );
262 
263  if (!readLevel && refHeader.headerOk())
264  {
266  << "Detected " << refHeader.name() << " file in "
267  << polyMesh::defaultRegion << " directory. Please remove to"
268  << " recreate it or use the -readLevel option to use it"
269  << endl;
270  return 1;
271  }
272 
273 
274  labelIOList refLevel
275  (
276  IOobject
277  (
278  "refinementLevel",
279  runTime.timeName(),
280  mesh,
283  ),
284  labelList(mesh.nCells(), 0)
285  );
286 
287  if (readLevel)
288  {
289  refLevel = labelIOList(refHeader);
290  }
291 
292  // Construct volScalarField with same info for post processing
293  volScalarField postRefLevel
294  (
295  IOobject
296  (
297  "refinementLevel",
298  runTime.timeName(),
299  mesh,
302  ),
303  fMesh,
304  dimensionedScalar("zero", dimless/dimTime, 0)
305  );
306 
307  // Set cell values
308  forAll(bins, binI)
309  {
310  const DynamicList<label>& bin = bins[binI]();
311 
312  forAll(bin, i)
313  {
314  refLevel[bin[i]] = bins.size() - binI - 1;
315  postRefLevel[bin[i]] = refLevel[bin[i]];
316  }
317  }
318 
319 
320  // For volScalarField: set boundary values to same as cell.
321  // Note: could also put
322  // zeroGradient b.c. on postRefLevel and do evaluate.
323  forAll(postRefLevel.boundaryField(), patchI)
324  {
325  const polyPatch& pp = patches[patchI];
326 
327  fvPatchScalarField& bField = postRefLevel.boundaryField()[patchI];
328 
329  Info<< "Setting field for patch "<< endl;
330 
331  forAll(bField, faceI)
332  {
333  label own = mesh.faceOwner()[pp.start() + faceI];
334 
335  bField[faceI] = postRefLevel[own];
336  }
337  }
338 
339  Info<< "Determined current refinement level and writing to "
340  << postRefLevel.name() << " (as volScalarField; for post processing)"
341  << nl
342  << polyMesh::defaultRegion/refLevel.name()
343  << " (as labelIOList; for meshing)" << nl
344  << endl;
345 
346  refLevel.write();
347  postRefLevel.write();
348 
349 
350  // Find out cells to refine to keep to 2:1 refinement level restriction
351 
352  // Cells to refine
353  cellSet refCells(mesh, "refCells", 100);
354 
355  while
356  (
357  limitRefinementLevel
358  (
359  mesh,
360  refLevel, // current refinement level
361  refCells // cells to refine
362  )
363  )
364  {}
365 
366  if (refCells.size())
367  {
368  Info<< "Collected " << refCells.size() << " cells that need to be"
369  << " refined to get closer to overall 2:1 refinement level limit"
370  << nl
371  << "Written cells to be refined to cellSet " << refCells.name()
372  << nl << endl;
373 
374  refCells.write();
375 
376  Info<< "After refinement this tool can be run again to see if the 2:1"
377  << " limit is observed all over the mesh" << nl << endl;
378  }
379  else
380  {
381  Info<< "All cells in the mesh observe the 2:1 refinement level limit"
382  << nl << endl;
383  }
384 
385  Info << nl << "End" << endl;
386 
387  return 0;
388 }
389 
390 
391 // ************************ vim: set sw=4 sts=4 et: ************************ //