62bool 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"
102int 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");
141 limits.emplace_back(sortedVols[0], 1.1*sortedVols[0]);
145 if (sortedVols[i] >
limits.back().max())
148 Info<<
"Collected " << bins.
back() <<
" elements in bin "
149 <<
limits.back().min() <<
" .. "
154 limits.emplace_back(sortedVols[i], 1.1 * sortedVols[i]);
156 Info<<
"Creating new bin "
157 <<
limits.back().min() <<
" .. "
162 bins.
back().push_back(sortedVols.indices()[i]);
170 Info<<
"Volume bins:" <<
nl;
173 const auto& bin = bins[bini];
179 <<
" : writing " << bin.size() <<
" cells to cellSet "
214 forAll(newPatches, patchi)
219 patches[patchi].clone(fMesh.boundaryMesh())
223 fMesh.addFvPatches(newPatches);
235 if (!readLevel && refHeader.typeHeaderOk<
labelIOList>(
true))
238 <<
"Detected " << refHeader.name() <<
" file in "
240 <<
" recreate it or use the -readLevel option to use it"
282 const auto& bin = bins[bini];
286 refLevel[bin[i]] = bins.
size() - bini - 1;
287 postRefLevel[bin[i]] = refLevel[bin[i]];
292 postRefLevel.boundaryFieldRef();
297 forAll(postRefLevel.boundaryField(), patchi)
303 Info<<
"Setting field for patch "<<
endl;
307 label own =
mesh.faceOwner()[
pp.start() + facei];
309 bField[facei] = postRefLevel[own];
313 Info<<
"Determined current refinement level and writing to "
314 << postRefLevel.name() <<
" (as volScalarField; for post processing)"
317 <<
" (as labelIOList; for meshing)" <<
nl
321 postRefLevel.write();
342 Info<<
"Collected " << refCells.
size() <<
" cells that need to be"
343 <<
" refined to get closer to overall 2:1 refinement level limit"
345 <<
"Written cells to be refined to cellSet " << refCells.
name()
350 Info<<
"After refinement this tool can be run again to see if the 2:1"
351 <<
" limit is observed all over the mesh" <<
nl <<
endl;
355 Info<<
"All cells in the mesh observe the 2:1 refinement level limit"
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
T & emplace_back(Args &&... args)
Construct an element at the end of the list, return reference to the new list element.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
label size() const noexcept
The number of elements in table.
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const word & name() const noexcept
Return the object name.
An indirect list with addressing based on sorting. The list is sorted upon construction or when expli...
T & back()
Access last element of the list, position [size()-1].
void size(const label n)
Older name for setAddressableSize.
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
static void addNote(const string ¬e)
Add extra notes for the usage information.
A collection of cell labels.
Mesh data needed to do the Finite Volume discretisation.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
static word defaultRegion
Return the default region name.
A patch is a list of labels that address the faces in the global face list.
Cell-face mesh analysis engine.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
const polyBoundaryMesh & patches
#define WarningInFunction
Report a warning using Foam::Warning.
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
const dimensionSet dimless
Dimensionless.
List< labelList > labelListList
List of labelList.
List< label > labelList
A List of labels.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
IOList< label > labelIOList
IO for a List of label.
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
List< cell > cellList
List of cell.
static constexpr const zero Zero
Global zero (0).
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fvPatchField< scalar > fvPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
Foam::argList args(argc, argv)
#define forAll(list, i)
Loop across all elements in list.