1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 3.2
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
9 This file is part of foam-extend.
11 foam-extend is free software: you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by the
13 Free Software Foundation, either version 3 of the License, or (at your
14 option) any later version.
16 foam-extend is distributed in the hope that it will be useful, but
17 WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 General Public License for more details.
21 You should have received a copy of the GNU General Public License
22 along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
25 Tries to figure out what the refinement level is on refined cartesian
26 meshes. Run BEFORE snapping.
29 - volScalarField 'refinementLevel' with current refinement level.
30 - cellSet 'refCells' which are the cells that need to be refined to satisfy
33 Works by dividing cells into volume bins.
35 \*---------------------------------------------------------------------------*/
38 #include "objectRegistry.H"
42 #include "SortableList.H"
43 #include "labelIOList.H"
45 #include "volFields.H"
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 // Return true if any cells had to be split to keep a difference between
52 // neighbouring refinement levels < limitDiff. Puts cells into refCells and
53 // update refLevel to account for refinement.
54 bool limitRefinementLevel
56 const primitiveMesh& mesh,
61 const labelListList& cellCells = mesh.cellCells();
63 label oldNCells = refCells.size();
65 forAll(cellCells, cellI)
67 const labelList& cCells = cellCells[cellI];
71 if (refLevel[cCells[i]] > (refLevel[cellI]+1))
73 // Found neighbour with >=2 difference in refLevel.
74 refCells.insert(cellI);
81 if (refCells.size() > oldNCells)
83 Info<< "Added an additional " << refCells.size() - oldNCells
84 << " cells to satisfy 1:2 refinement level"
98 int main(int argc, char *argv[])
100 argList::validOptions.insert("readLevel", "");
102 # include "setRootCase.H"
103 # include "createTime.H"
104 # include "createPolyMesh.H"
106 Info<< "Dividing cells into bins depending on cell volume.\nThis will"
107 << " correspond to refinement levels for a mesh with only 2x2x2"
109 << "The upper range for every bin is always 1.1 times the lower range"
110 << " to allow for some truncation error."
113 bool readLevel = args.optionFound("readLevel");
115 const scalarField& vols = mesh.cellVolumes();
117 SortableList<scalar> sortedVols(vols);
119 // All cell labels, sorted per bin.
120 DynamicList<DynamicList<label> > bins;
122 // Lower/upper limits
123 DynamicList<scalar> lowerLimits;
124 DynamicList<scalar> upperLimits;
126 // Create bin0. Have upperlimit as factor times lowerlimit.
127 bins.append(DynamicList<label>());
128 lowerLimits.append(sortedVols[0]);
129 upperLimits.append(1.1*lowerLimits[lowerLimits.size()-1]);
131 forAll(sortedVols, i)
133 if (sortedVols[i] > upperLimits[upperLimits.size()-1])
135 // New value outside of current bin
138 DynamicList<label>& bin = bins[bins.size()-1];
142 Info<< "Collected " << bin.size() << " elements in bin "
143 << lowerLimits[lowerLimits.size()-1] << " .. "
144 << upperLimits[upperLimits.size()-1] << endl;
147 bins.append(DynamicList<label>());
148 lowerLimits.append(sortedVols[i]);
149 upperLimits.append(1.1*lowerLimits[lowerLimits.size()-1]);
151 Info<< "Creating new bin " << lowerLimits[lowerLimits.size()-1]
152 << " .. " << upperLimits[upperLimits.size()-1]
156 // Append to current bin.
157 DynamicList<label>& bin = bins[bins.size()-1];
159 bin.append(sortedVols.indices()[i]);
163 bins[bins.size()-1].shrink();
165 lowerLimits.shrink();
166 upperLimits.shrink();
170 // Write to cellSets.
173 Info<< "Volume bins:" << nl;
176 const DynamicList<label>& bin = bins[binI];
178 cellSet cells(mesh, "vol" + name(binI), bin.size());
182 cells.insert(bin[i]);
185 Info<< " " << lowerLimits[binI] << " .. " << upperLimits[binI]
186 << " : writing " << bin.size() << " cells to cellSet "
187 << cells.name() << endl;
195 // Convert bins into refinement level.
199 // Construct fvMesh to be able to construct volScalarField
205 fvMesh::defaultRegion,
209 xferCopy(mesh.points()), // could we safely re-use the data?
210 xferCopy(mesh.faces()),
211 xferCopy(mesh.cells())
214 // Add the boundary patches
215 const polyBoundaryMesh& patches = mesh.boundaryMesh();
217 List<polyPatch*> p(patches.size());
221 p[patchI] = patches[patchI].clone(fMesh.boundaryMesh()).ptr();
224 fMesh.addFvPatches(p);
232 polyMesh::defaultRegion,
236 if (!readLevel && refHeader.headerOk())
238 WarningIn(args.executable())
239 << "Detected " << refHeader.name() << " file in "
240 << polyMesh::defaultRegion << " directory. Please remove to"
241 << " recreate it or use the -readLevel option to use it"
257 labelList(mesh.nCells(), 0)
262 refLevel = labelIOList(refHeader);
265 // Construct volScalarField with same info for post processing
266 volScalarField postRefLevel
277 dimensionedScalar("zero", dimless/dimTime, 0)
283 const DynamicList<label>& bin = bins[binI];
287 refLevel[bin[i]] = bins.size() - binI - 1;
288 postRefLevel[bin[i]] = refLevel[bin[i]];
293 // For volScalarField: set boundary values to same as cell.
294 // Note: could also put
295 // zeroGradient b.c. on postRefLevel and do evaluate.
296 forAll(postRefLevel.boundaryField(), patchI)
298 const polyPatch& pp = patches[patchI];
300 fvPatchScalarField& bField = postRefLevel.boundaryField()[patchI];
302 Info<< "Setting field for patch "<< endl;
304 forAll(bField, faceI)
306 label own = mesh.faceOwner()[pp.start() + faceI];
308 bField[faceI] = postRefLevel[own];
312 Info<< "Determined current refinement level and writing to "
313 << postRefLevel.name() << " (as volScalarField; for post processing)"
315 << polyMesh::defaultRegion/refLevel.name()
316 << " (as labelIOList; for meshing)" << nl
320 postRefLevel.write();
323 // Find out cells to refine to keep to 2:1 refinement level restriction
326 cellSet refCells(mesh, "refCells", 100);
333 refLevel, // current refinement level
334 refCells // cells to refine
341 Info<< "Collected " << refCells.size() << " cells that need to be"
342 << " refined to get closer to overall 2:1 refinement level limit"
344 << "Written cells to be refined to cellSet " << refCells.name()
349 Info<< "After refinement this tool can be run again to see if the 2:1"
350 << " limit is observed all over the mesh" << nl << endl;
354 Info<< "All cells in the mesh observe the 2:1 refinement level limit"
358 Info << nl << "End" << endl;
364 // ************************************************************************* //