File: refinementLevel.C

package info (click to toggle)
openfoam 4.1%2Bdfsg1-1
  • links: PTS, VCS
  • area: main
  • in suites: stretch
  • size: 163,028 kB
  • ctags: 58,990
  • sloc: cpp: 830,760; sh: 10,227; ansic: 8,215; xml: 745; lex: 437; awk: 194; sed: 91; makefile: 77; python: 18
file content (370 lines) | stat: -rw-r--r-- 9,646 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
    \\  /    A nd           | Copyright (C) 2011-2016 OpenFOAM Foundation
     \\/     M anipulation  |
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    You should have received a copy of the GNU General Public License
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.

Application
    refinementLevel

Description
    Tries to figure out what the refinement level is on refined cartesian
    meshes. Run BEFORE snapping.

    Writes
    - volScalarField 'refinementLevel' with current refinement level.
    - cellSet 'refCells' which are the cells that need to be refined to satisfy
      2:1 refinement.

    Works by dividing cells into volume bins.

\*---------------------------------------------------------------------------*/

#include "argList.H"
#include "Time.H"
#include "polyMesh.H"
#include "cellSet.H"
#include "SortableList.H"
#include "labelIOList.H"
#include "fvMesh.H"
#include "volFields.H"

using namespace Foam;

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

// Return true if any cells had to be split to keep a difference between
// neighbouring refinement levels < limitDiff. Puts cells into refCells and
// update refLevel to account for refinement.
bool limitRefinementLevel
(
    const primitiveMesh& mesh,
    labelList& refLevel,
    cellSet& refCells
)
{
    const labelListList& cellCells = mesh.cellCells();

    label oldNCells = refCells.size();

    forAll(cellCells, celli)
    {
        const labelList& cCells = cellCells[celli];

        forAll(cCells, i)
        {
            if (refLevel[cCells[i]] > (refLevel[celli]+1))
            {
                // Found neighbour with >=2 difference in refLevel.
                refCells.insert(celli);
                refLevel[celli]++;
                break;
            }
        }
    }

    if (refCells.size() > oldNCells)
    {
        Info<< "Added an additional " << refCells.size() - oldNCells
            << " cells to satisfy 1:2 refinement level"
            << endl;

        return true;
    }
    else
    {
        return false;
    }
}



int main(int argc, char *argv[])
{
    argList::addBoolOption
    (
        "readLevel",
        "read level from refinementLevel file"
    );

    #include "setRootCase.H"
    #include "createTime.H"
    #include "createPolyMesh.H"

    Info<< "Dividing cells into bins depending on cell volume.\nThis will"
        << " correspond to refinement levels for a mesh with only 2x2x2"
        << " refinement\n"
        << "The upper range for every bin is always 1.1 times the lower range"
        << " to allow for some truncation error."
        << nl << endl;

    const bool readLevel = args.optionFound("readLevel");

    const scalarField& vols = mesh.cellVolumes();

    SortableList<scalar> sortedVols(vols);

    // All cell labels, sorted per bin.
    DynamicList<DynamicList<label>> bins;

    // Lower/upper limits
    DynamicList<scalar> lowerLimits;
    DynamicList<scalar> upperLimits;

    // Create bin0. Have upperlimit as factor times lowerlimit.
    bins.append(DynamicList<label>());
    lowerLimits.append(sortedVols[0]);
    upperLimits.append(1.1 * lowerLimits.last());

    forAll(sortedVols, i)
    {
        if (sortedVols[i] > upperLimits.last())
        {
            // New value outside of current bin

            // Shrink old bin.
            DynamicList<label>& bin = bins.last();

            bin.shrink();

            Info<< "Collected " << bin.size() << " elements in bin "
                << lowerLimits.last() << " .. "
                << upperLimits.last() << endl;

            // Create new bin.
            bins.append(DynamicList<label>());
            lowerLimits.append(sortedVols[i]);
            upperLimits.append(1.1 * lowerLimits.last());

            Info<< "Creating new bin " << lowerLimits.last()
                << " .. " << upperLimits.last()
                << endl;
        }

        // Append to current bin.
        DynamicList<label>& bin = bins.last();

        bin.append(sortedVols.indices()[i]);
    }
    Info<< endl;

    bins.last().shrink();
    bins.shrink();
    lowerLimits.shrink();
    upperLimits.shrink();


    //
    // Write to cellSets.
    //

    Info<< "Volume bins:" << nl;
    forAll(bins, binI)
    {
        const DynamicList<label>& bin = bins[binI];

        cellSet cells(mesh, "vol" + name(binI), bin.size());

        forAll(bin, i)
        {
            cells.insert(bin[i]);
        }

        Info<< "    " << lowerLimits[binI] << " .. " << upperLimits[binI]
            << "  : writing " << bin.size() << " cells to cellSet "
            << cells.name() << endl;

        cells.write();
    }



    //
    // Convert bins into refinement level.
    //


    // Construct fvMesh to be able to construct volScalarField

    fvMesh fMesh
    (
        IOobject
        (
            fvMesh::defaultRegion,
            runTime.timeName(),
            runTime
        ),
        xferCopy(mesh.points()),   // could we safely re-use the data?
        xferCopy(mesh.faces()),
        xferCopy(mesh.cells())
    );

    // Add the boundary patches
    const polyBoundaryMesh& patches = mesh.boundaryMesh();

    List<polyPatch*> p(patches.size());

    forAll(p, patchi)
    {
        p[patchi] = patches[patchi].clone(fMesh.boundaryMesh()).ptr();
    }

    fMesh.addFvPatches(p);


    // Refinement level
    IOobject refHeader
    (
        "refinementLevel",
        runTime.timeName(),
        polyMesh::defaultRegion,
        runTime
    );

    if (!readLevel && refHeader.headerOk())
    {
        WarningInFunction
            << "Detected " << refHeader.name() << " file in "
            << polyMesh::defaultRegion <<  " directory. Please remove to"
            << " recreate it or use the -readLevel option to use it"
            << endl;
        return 1;
    }


    labelIOList refLevel
    (
        IOobject
        (
            "refinementLevel",
            runTime.timeName(),
            mesh,
            IOobject::NO_READ,
            IOobject::AUTO_WRITE
        ),
        labelList(mesh.nCells(), 0)
    );

    if (readLevel)
    {
        refLevel = labelIOList(refHeader);
    }

    // Construct volScalarField with same info for post processing
    volScalarField postRefLevel
    (
        IOobject
        (
            "refinementLevel",
            runTime.timeName(),
            mesh,
            IOobject::NO_READ,
            IOobject::NO_WRITE
        ),
        fMesh,
        dimensionedScalar("zero", dimless/dimTime, 0)
    );

    // Set cell values
    forAll(bins, binI)
    {
        const DynamicList<label>& bin = bins[binI];

        forAll(bin, i)
        {
            refLevel[bin[i]] = bins.size() - binI - 1;
            postRefLevel[bin[i]] = refLevel[bin[i]];
        }
    }

    volScalarField::Boundary& postRefLevelBf =
        postRefLevel.boundaryFieldRef();

    // For volScalarField: set boundary values to same as cell.
    // Note: could also put
    // zeroGradient b.c. on postRefLevel and do evaluate.
    forAll(postRefLevel.boundaryField(), patchi)
    {
        const polyPatch& pp = patches[patchi];

        fvPatchScalarField& bField = postRefLevelBf[patchi];

        Info<< "Setting field for patch "<< endl;

        forAll(bField, facei)
        {
            label own = mesh.faceOwner()[pp.start() + facei];

            bField[facei] = postRefLevel[own];
        }
    }

    Info<< "Determined current refinement level and writing to "
        << postRefLevel.name() << " (as volScalarField; for post processing)"
        << nl
        << polyMesh::defaultRegion/refLevel.name()
        << " (as labelIOList; for meshing)" << nl
        << endl;

    refLevel.write();
    postRefLevel.write();


    // Find out cells to refine to keep to 2:1 refinement level restriction

    // Cells to refine
    cellSet refCells(mesh, "refCells", 100);

    while
    (
        limitRefinementLevel
        (
            mesh,
            refLevel,       // current refinement level
            refCells        // cells to refine
        )
    )
    {}

    if (refCells.size())
    {
        Info<< "Collected " << refCells.size() << " cells that need to be"
            << " refined to get closer to overall 2:1 refinement level limit"
            << nl
            << "Written cells to be refined to cellSet " << refCells.name()
            << nl << endl;

        refCells.write();

        Info<< "After refinement this tool can be run again to see if the 2:1"
            << " limit is observed all over the mesh" << nl << endl;
    }
    else
    {
        Info<< "All cells in the mesh observe the 2:1 refinement level limit"
            << nl << endl;
    }

    Info<< "\nEnd\n" << endl;
    return 0;
}


// ************************************************************************* //