Go to the documentation of this file.
47 void Foam::LESModels::smoothDelta::setChangedFaces
51 DynamicList<label>& changedFaces,
52 DynamicList<deltaData>& changedFacesInfo
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);
90 changedFacesInfo.append(deltaData(ownDelta));
95 changedFaces.shrink();
96 changedFacesInfo.shrink();
100 void 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();
148 Foam::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())
static autoPtr< LESdelta > New(const word &name, const turbulenceModel &turbulence, const dictionary &dict, const word &lookupName="delta")
Defines the attributes of an object for which implicit objectRegistry management is supported,...
A class for handling words, derived from Foam::string.
Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
constexpr const char *const group
Smoothed delta which takes a given simple geometric delta and applies smoothing to it such that the r...
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Abstract base class for turbulence models (RAS, LES and laminar).
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Macros for easy insertion into run-time selection tables.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Abstract base class for LES deltas.
defineTypeNameAndDebug(cubeRootVolDelta, 0)
fileName::Type type(const fileName &name, const bool followLink=true)
virtual void read(const dictionary &dict)
addToRunTimeSelectionTable(LESfluidThermoCompressibleTurbulenceModel, SmagorinskyLESfluidThermoCompressibleTurbulenceModel, dictionary)
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
word name(const expressions::valueTypeCode typeCode)
static word groupName(StringType base, const word &group)