44 void Foam::LESModels::smoothDelta::setChangedFaces
48 DynamicList<labelPair>& changedFaces,
49 DynamicList<deltaData>& changedFacesInfo
52 forAll(mesh.owner(), facei)
54 const scalar ownDelta = delta[mesh.owner()[facei]];
55 const scalar neiDelta = delta[mesh.neighbour()[facei]];
58 if (ownDelta > maxDeltaRatio_*neiDelta)
60 changedFaces.append(
labelPair(-1, facei));
61 changedFacesInfo.append(deltaData(ownDelta));
63 else if (neiDelta > maxDeltaRatio_*ownDelta)
65 changedFaces.append(
labelPair(-1, facei));
66 changedFacesInfo.append(deltaData(neiDelta));
74 const fvPatch& patch = mesh.boundary()[
patchi];
80 const scalar ownDelta = delta[patch.faceCells()[patchFacei]];
83 changedFacesInfo.append(deltaData(ownDelta));
88 changedFaces.shrink();
89 changedFacesInfo.shrink();
93 void Foam::LESModels::smoothDelta::calcDelta()
95 const fvMesh& mesh = momentumTransportModel_.mesh();
100 DynamicList<labelPair> changedFaces(mesh.nFaces()/100 + 100);
101 DynamicList<deltaData> changedFacesInfo(changedFaces.size());
102 setChangedFaces(mesh, geometricDelta, changedFaces, changedFacesInfo);
105 List<deltaData> cellDeltaData(mesh.nCells());
106 forAll(geometricDelta, celli)
108 cellDeltaData[celli] = geometricDelta[celli];
112 List<deltaData> internalFaceDeltaData(mesh.nInternalFaces());
113 List<List<deltaData>> patchFaceDeltaData
115 FvFaceCellWave<deltaData>::template
116 sizesListList<List<List<deltaData>>>
118 FvFaceCellWave<deltaData>::template
119 listListSizes(mesh.boundary()),
125 FvFaceCellWave<deltaData, scalar> deltaCalc
130 internalFaceDeltaData,
133 mesh.globalData().nTotalCells() + 1,
139 delta_[celli] = cellDeltaData[celli].delta();
178 geometricDelta_().
read(coeffsDict);
179 coeffsDict.lookup(
"maxDeltaRatio") >> maxDeltaRatio_;
186 geometricDelta_().correct();
188 if (momentumTransportModel_.mesh().changing())
#define forAll(list, i)
Loop across all elements in list.
static autoPtr< LESdelta > New(const word &name, const momentumTransportModel &turbulence, const dictionary &dict)
Return a reference to the selected LES delta.
virtual void read(const dictionary &)
Read the LESdelta dictionary.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Abstract base class for LES deltas.
addToRunTimeSelectionTable(LESdelta, cubeRootVolDelta, dictionary)
bool read(Istream &, const bool keepHeader=false)
Read dictionary from Istream, optionally keeping the header.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
A class for handling words, derived from string.
Pair< label > labelPair
Label pair.
defineTypeNameAndDebug(cubeRootVolDelta, 0)
Abstract base class for turbulence models (RAS, LES and laminar).
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
smoothDelta(const word &name, const momentumTransportModel &turbulence, const dictionary &)
Construct from name, momentumTransportModel and dictionary.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.