Go to the documentation of this file.
52 for (
label faceI = 0; faceI <
mesh.nInternalFaces(); faceI++)
54 scalar ownDelta =
delta[
mesh.faceOwner()[faceI]];
56 scalar neiDelta =
delta[
mesh.faceNeighbour()[faceI]];
60 if (ownDelta > maxDeltaRatio_ * neiDelta)
62 changedFaces.
append(faceI);
65 else if (neiDelta > maxDeltaRatio_ * ownDelta)
67 changedFaces.
append(faceI);
84 scalar ownDelta =
delta[
mesh.faceOwner()[meshFaceI]];
86 changedFaces.
append(meshFaceI);
112 forAll(geometricDelta, cellI)
114 cellDeltaData[cellI] = geometricDelta[cellI];
129 mesh.globalData().nTotalCells()+1,
135 delta_[cellI] = cellDeltaData[cellI].delta();
177 geometricDelta_().read(coeffsDict);
178 coeffsDict.
lookup(
"maxDeltaRatio") >> maxDeltaRatio_;
185 geometricDelta_().correct();
187 if (turbulenceModel_.mesh().changing())
virtual void read(const dictionary &)
Read the LESdelta dictionary.
A class for handling words, derived from string.
#define forAll(list, i)
Loop across all elements in list.
Public member class used by mesh-wave to propagate the delta-ratio.
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
smoothDelta(const smoothDelta &)
Disallow default bitwise copy construct and assignment.
Mesh consisting of general polyhedral cells.
Smoothed delta which takes a given simple geometric delta and applies smoothing to it such that the r...
autoPtr< LESdelta > geometricDelta_
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
autoPtr< compressible::turbulenceModel > turbulence
A patch is a list of labels that address the faces in the global face list.
const turbulenceModel & turbulenceModel_
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Abstract base class for turbulence models (RAS, LES and laminar).
A list of keyword definitions, which are a keyword followed by any number of values (e....
Macros for easy insertion into run-time selection tables.
Mesh data needed to do the Finite Volume discretisation.
label start() const
Return start label of this patch in the polyMesh face list.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
void correctBoundaryConditions()
Correct boundary field.
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Abstract base class for LES deltas.
defineTypeNameAndDebug(cubeRootVolDelta, 0)
static autoPtr< LESdelta > New(const word &name, const turbulenceModel &turbulence, const dictionary &dict)
Return a reference to the selected LES delta.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
const fvMesh & mesh() const
addToRunTimeSelectionTable(LESfluidThermoCompressibleTurbulenceModel, SmagorinskyLESfluidThermoCompressibleTurbulenceModel, dictionary)
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Generic GeometricField class.
word name(const complex &)
Return a string representation of a complex.
stressControl lookup("compactNormalStress") >> compactNormalStress
void setChangedFaces(const polyMesh &mesh, const volScalarField &delta, DynamicList< label > &changedFaces, DynamicList< deltaData > &changedFacesInfo)
Fill changedFaces (with face labels) and changedFacesInfo.