47void Foam::LESModels::smoothDelta::setChangedFaces
55 for (label facei = 0; facei <
mesh.nInternalFaces(); facei++)
57 scalar ownDelta =
delta[
mesh.faceOwner()[facei]];
59 scalar neiDelta =
delta[
mesh.faceNeighbour()[facei]];
63 if (ownDelta > maxDeltaRatio_ * neiDelta)
65 changedFaces.
append(facei);
66 changedFacesInfo.
append(deltaData(ownDelta));
68 else if (neiDelta > maxDeltaRatio_ * ownDelta)
70 changedFaces.
append(facei);
71 changedFacesInfo.
append(deltaData(neiDelta));
79 const polyPatch&
patch =
mesh.boundaryMesh()[patchi];
85 label meshFacei =
patch.start() + patchFacei;
87 scalar ownDelta =
delta[
mesh.faceOwner()[meshFacei]];
89 changedFaces.
append(meshFacei);
100void Foam::LESModels::smoothDelta::calcDelta()
102 const fvMesh&
mesh = turbulenceModel_.mesh();
107 DynamicList<label> changedFaces(
mesh.nFaces()/100 + 100);
108 DynamicList<deltaData> changedFacesInfo(changedFaces.
size());
110 setChangedFaces(
mesh, geometricDelta, changedFaces, changedFacesInfo);
113 List<deltaData> cellDeltaData(
mesh.nCells());
115 forAll(geometricDelta, celli)
117 cellDeltaData[celli] = geometricDelta[celli];
121 List<deltaData> faceDeltaData(
mesh.nFaces());
125 FaceCellWave<deltaData, scalar> deltaCalc
132 mesh.globalData().nTotalCells()+1,
138 delta_[celli] = cellDeltaData[celli].delta();
142 delta_.correctBoundaryConditions();
148Foam::LESModels::smoothDelta::smoothDelta
162 dict.optionalSubDict(
type() +
"Coeffs")
167 dict.optionalSubDict(
type() +
"Coeffs").get<scalar>(
"maxDeltaRatio")
180 geometricDelta_().read(coeffsDict);
181 coeffsDict.readEntry(
"maxDeltaRatio", maxDeltaRatio_);
188 geometricDelta_().correct();
190 if (turbulenceModel_.mesh().changing())
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
void append(const T &val)
Copy append an element to the end of this list.
DynamicList< T, SizeMin > & shrink()
Calls shrink_to_fit() and returns a reference to the DynamicList.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Public class used by mesh-wave to propagate the delta-ratio.
virtual void read(const dictionary &dict)
Read the LESdelta dictionary.
Abstract base class for LES deltas.
LESdelta(const LESdelta &)=delete
No copy construct.
static autoPtr< LESdelta > New(const word &name, const turbulenceModel &turbulence, const dictionary &dict, const word &lookupName="delta")
Return a reference to the selected LES delta.
const turbulenceModel & turbulence() const
Return turbulenceModel reference.
const turbulenceModel & turbulenceModel_
void size(const label n)
Older name for setAddressableSize.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect,...
Mesh consisting of general polyhedral cells.
Abstract base class for turbulence models (RAS, LES and laminar).
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Namespace for LES SGS models.
const std::string patch
OpenFOAM patch number as a std::string.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
#define forAll(list, i)
Loop across all elements in list.