Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062 #include <OpenFOAM/argList.H>
00063 #include <OpenFOAM/Time.H>
00064 #include <OpenFOAM/polyMesh.H>
00065 #include <meshTools/cellSet.H>
00066 #include <OpenFOAM/SortableList.H>
00067 #include <OpenFOAM/labelIOList.H>
00068 #include <finiteVolume/fvMesh.H>
00069 #include <finiteVolume/volFields.H>
00070
00071 using namespace Foam;
00072
00073
00074
00075
00076
00077
00078 bool limitRefinementLevel
00079 (
00080 const primitiveMesh& mesh,
00081 labelList& refLevel,
00082 cellSet& refCells
00083 )
00084 {
00085 const labelListList& cellCells = mesh.cellCells();
00086
00087 label oldNCells = refCells.size();
00088
00089 forAll(cellCells, cellI)
00090 {
00091 const labelList& cCells = cellCells[cellI];
00092
00093 forAll(cCells, i)
00094 {
00095 if (refLevel[cCells[i]] > (refLevel[cellI]+1))
00096 {
00097
00098 refCells.insert(cellI);
00099 refLevel[cellI]++;
00100 break;
00101 }
00102 }
00103 }
00104
00105 if (refCells.size() > oldNCells)
00106 {
00107 Info<< "Added an additional " << refCells.size() - oldNCells
00108 << " cells to satisfy 1:2 refinement level"
00109 << endl;
00110
00111 return true;
00112 }
00113 else
00114 {
00115 return false;
00116 }
00117 }
00118
00119
00120
00121
00122 int main(int argc, char *argv[])
00123 {
00124 argList::validOptions.insert("readLevel", "");
00125
00126 # include <OpenFOAM/setRootCase.H>
00127 # include <OpenFOAM/createTime.H>
00128 # include <OpenFOAM/createPolyMesh.H>
00129
00130 Info<< "Dividing cells into bins depending on cell volume.\nThis will"
00131 << " correspond to refinement levels for a mesh with only 2x2x2"
00132 << " refinement\n"
00133 << "The upper range for every bin is always 1.1 times the lower range"
00134 << " to allow for some truncation error."
00135 << nl << endl;
00136
00137 bool readLevel = args.optionFound("readLevel");
00138
00139 const scalarField& vols = mesh.cellVolumes();
00140
00141 SortableList<scalar> sortedVols(vols);
00142
00143
00144
00145
00146
00147 DynamicList<autoPtr<DynamicList<label> > > bins;
00148
00149
00150 DynamicList<scalar> lowerLimits;
00151 DynamicList<scalar> upperLimits;
00152
00153
00154 bins.append(autoPtr<DynamicList<label> >(new DynamicList<label>()));
00155 lowerLimits.append(sortedVols[0]);
00156 upperLimits.append(1.1*lowerLimits[lowerLimits.size()-1]);
00157
00158 forAll(sortedVols, i)
00159 {
00160 if (sortedVols[i] > upperLimits[upperLimits.size()-1])
00161 {
00162
00163
00164
00165 DynamicList<label>& bin = bins[bins.size()-1]();
00166
00167 bin.shrink();
00168
00169 Info<< "Collected " << bin.size() << " elements in bin "
00170 << lowerLimits[lowerLimits.size()-1] << " .. "
00171 << upperLimits[upperLimits.size()-1] << endl;
00172
00173
00174 bins.append(autoPtr<DynamicList<label> >(new DynamicList<label>()));
00175 lowerLimits.append(sortedVols[i]);
00176 upperLimits.append(1.1*lowerLimits[lowerLimits.size()-1]);
00177
00178 Info<< "Creating new bin " << lowerLimits[lowerLimits.size()-1]
00179 << " .. " << upperLimits[upperLimits.size()-1]
00180 << endl;
00181 }
00182
00183
00184 DynamicList<label>& bin = bins[bins.size()-1]();
00185
00186 bin.append(sortedVols.indices()[i]);
00187 }
00188 Info<< endl;
00189
00190 bins[bins.size()-1]->shrink();
00191 bins.shrink();
00192 lowerLimits.shrink();
00193 upperLimits.shrink();
00194
00195
00196
00197
00198
00199
00200 Info<< "Volume bins:" << nl;
00201 forAll(bins, binI)
00202 {
00203 const DynamicList<label>& bin = bins[binI]();
00204
00205 cellSet cells(mesh, "vol" + name(binI), bin.size());
00206
00207 forAll(bin, i)
00208 {
00209 cells.insert(bin[i]);
00210 }
00211
00212 Info<< " " << lowerLimits[binI] << " .. " << upperLimits[binI]
00213 << " : writing " << bin.size() << " cells to cellSet "
00214 << cells.name() << endl;
00215
00216 cells.write();
00217 }
00218
00219
00220
00221
00222
00223
00224
00225
00226
00227
00228 fvMesh fMesh
00229 (
00230 IOobject
00231 (
00232 fvMesh::defaultRegion,
00233 runTime.timeName(),
00234 runTime
00235 ),
00236 xferCopy(mesh.points()),
00237 xferCopy(mesh.faces()),
00238 xferCopy(mesh.cells())
00239 );
00240
00241
00242 const polyBoundaryMesh& patches = mesh.boundaryMesh();
00243
00244 List<polyPatch*> p(patches.size());
00245
00246 forAll (p, patchI)
00247 {
00248 p[patchI] = patches[patchI].clone(fMesh.boundaryMesh()).ptr();
00249 }
00250
00251 fMesh.addFvPatches(p);
00252
00253
00254
00255 IOobject refHeader
00256 (
00257 "refinementLevel",
00258 runTime.timeName(),
00259 polyMesh::defaultRegion,
00260 runTime
00261 );
00262
00263 if (!readLevel && refHeader.headerOk())
00264 {
00265 WarningIn(args.executable())
00266 << "Detected " << refHeader.name() << " file in "
00267 << polyMesh::defaultRegion << " directory. Please remove to"
00268 << " recreate it or use the -readLevel option to use it"
00269 << endl;
00270 return 1;
00271 }
00272
00273
00274 labelIOList refLevel
00275 (
00276 IOobject
00277 (
00278 "refinementLevel",
00279 runTime.timeName(),
00280 mesh,
00281 IOobject::NO_READ,
00282 IOobject::AUTO_WRITE
00283 ),
00284 labelList(mesh.nCells(), 0)
00285 );
00286
00287 if (readLevel)
00288 {
00289 refLevel = labelIOList(refHeader);
00290 }
00291
00292
00293 volScalarField postRefLevel
00294 (
00295 IOobject
00296 (
00297 "refinementLevel",
00298 runTime.timeName(),
00299 mesh,
00300 IOobject::NO_READ,
00301 IOobject::NO_WRITE
00302 ),
00303 fMesh,
00304 dimensionedScalar("zero", dimless/dimTime, 0)
00305 );
00306
00307
00308 forAll(bins, binI)
00309 {
00310 const DynamicList<label>& bin = bins[binI]();
00311
00312 forAll(bin, i)
00313 {
00314 refLevel[bin[i]] = bins.size() - binI - 1;
00315 postRefLevel[bin[i]] = refLevel[bin[i]];
00316 }
00317 }
00318
00319
00320
00321
00322
00323 forAll(postRefLevel.boundaryField(), patchI)
00324 {
00325 const polyPatch& pp = patches[patchI];
00326
00327 fvPatchScalarField& bField = postRefLevel.boundaryField()[patchI];
00328
00329 Info<< "Setting field for patch "<< endl;
00330
00331 forAll(bField, faceI)
00332 {
00333 label own = mesh.faceOwner()[pp.start() + faceI];
00334
00335 bField[faceI] = postRefLevel[own];
00336 }
00337 }
00338
00339 Info<< "Determined current refinement level and writing to "
00340 << postRefLevel.name() << " (as volScalarField; for post processing)"
00341 << nl
00342 << polyMesh::defaultRegion/refLevel.name()
00343 << " (as labelIOList; for meshing)" << nl
00344 << endl;
00345
00346 refLevel.write();
00347 postRefLevel.write();
00348
00349
00350
00351
00352
00353 cellSet refCells(mesh, "refCells", 100);
00354
00355 while
00356 (
00357 limitRefinementLevel
00358 (
00359 mesh,
00360 refLevel,
00361 refCells
00362 )
00363 )
00364 {}
00365
00366 if (refCells.size())
00367 {
00368 Info<< "Collected " << refCells.size() << " cells that need to be"
00369 << " refined to get closer to overall 2:1 refinement level limit"
00370 << nl
00371 << "Written cells to be refined to cellSet " << refCells.name()
00372 << nl << endl;
00373
00374 refCells.write();
00375
00376 Info<< "After refinement this tool can be run again to see if the 2:1"
00377 << " limit is observed all over the mesh" << nl << endl;
00378 }
00379 else
00380 {
00381 Info<< "All cells in the mesh observe the 2:1 refinement level limit"
00382 << nl << endl;
00383 }
00384
00385 Info << nl << "End" << endl;
00386
00387 return 0;
00388 }
00389
00390
00391