Go to the documentation of this file.
62 bool limitRefinementLevel
71 label oldNCells = refCells.size();
75 const labelList& cCells = cellCells[celli];
79 if (refLevel[cCells[i]] > (refLevel[celli]+1))
89 if (refCells.size() > oldNCells)
91 Info<<
"Added an additional " << refCells.size() - oldNCells
92 <<
" cells to satisfy 1:2 refinement level"
102 int main(
int argc,
char *argv[])
106 "Attempt to determine refinement levels of a refined cartesian mesh.\n"
107 "Run BEFORE snapping!"
113 "Read level from refinementLevel file"
120 Info<<
"Dividing cells into bins depending on cell volume.\nThis will"
121 <<
" correspond to refinement levels for a mesh with only 2x2x2"
123 <<
"The upper range for every bin is always 1.1 times the lower range"
124 <<
" to allow for some truncation error."
127 const bool readLevel =
args.
found(
"readLevel");
142 lowerLimits.
append(sortedVols[0]);
143 upperLimits.
append(1.1 * lowerLimits.last());
147 if (sortedVols[i] > upperLimits.last())
156 Info<<
"Collected " << bin.size() <<
" elements in bin "
157 << lowerLimits.last() <<
" .. "
158 << upperLimits.last() <<
endl;
162 lowerLimits.
append(sortedVols[i]);
163 upperLimits.
append(1.1 * lowerLimits.last());
165 Info<<
"Creating new bin " << lowerLimits.last()
166 <<
" .. " << upperLimits.last()
173 bin.
append(sortedVols.indices()[i]);
187 Info<<
"Volume bins:" <<
nl;
195 Info<<
" " << lowerLimits[binI] <<
" .. " << upperLimits[binI]
196 <<
" : writing " << bin.size() <<
" cells to cellSet "
233 p[patchi] =
patches[patchi].
clone(fMesh.boundaryMesh()).ptr();
236 fMesh.addFvPatches(
p);
248 if (!readLevel && refHeader.typeHeaderOk<
labelIOList>(
true))
251 <<
"Detected " << refHeader.name() <<
" file in "
253 <<
" recreate it or use the -readLevel option to use it"
299 refLevel[bin[i]] = bins.size() - binI - 1;
300 postRefLevel[bin[i]] = refLevel[bin[i]];
305 postRefLevel.boundaryFieldRef();
310 forAll(postRefLevel.boundaryField(), patchi)
316 Info<<
"Setting field for patch "<<
endl;
322 bField[facei] = postRefLevel[own];
326 Info<<
"Determined current refinement level and writing to "
327 << postRefLevel.name() <<
" (as volScalarField; for post processing)"
330 <<
" (as labelIOList; for meshing)" <<
nl
334 postRefLevel.write();
355 Info<<
"Collected " << refCells.size() <<
" cells that need to be"
356 <<
" refined to get closer to overall 2:1 refinement level limit"
358 <<
"Written cells to be refined to cellSet " << refCells.
name()
363 Info<<
"After refinement this tool can be run again to see if the 2:1"
364 <<
" limit is observed all over the mesh" <<
nl <<
endl;
368 Info<<
"All cells in the mesh observe the 2:1 refinement level limit"
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual const pointField & points() const
Defines the attributes of an object for which implicit objectRegistry management is supported,...
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
static word defaultRegion
static constexpr const zero Zero
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
static word timeName(const scalar t, const int precision=precision_)
const cellList & cells() const
static void addNote(const string ¬e)
const polyBoundaryMesh & boundaryMesh() const
Ostream & endl(Ostream &os)
DynamicList< T, SizeMin > & shrink()
IOList< label > labelIOList
Label container classes.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
label nCells() const noexcept
Generic templated field type.
virtual bool write(const bool valid=true) const
A patch is a list of labels that address the faces in the global face list.
DynamicList< T, SizeMin > & append(const T &val)
virtual const labelList & faceOwner() const
const labelListList & cellCells() const
PtrList< T > clone(Args &&... args) const
A list that is sorted upon construction or when explicitly requested with the sort() method.
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
const scalarField & cellVolumes() const
A collection of cell labels.
const word & name() const noexcept
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
virtual const faceList & faces() const
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
bool insert(const Key &key)
const polyBoundaryMesh & patches
word name(const expressions::valueTypeCode typeCode)
Generic GeometricField class.
Foam::argList args(argc, argv)
#define WarningInFunction
const dimensionSet dimless
bool found(const word &optName) const
Cell-face mesh analysis engine.