Loading...
Searching...
No Matches
vanDriestDelta.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2020-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "vanDriestDelta.H"
30#include "wallFvPatch.H"
32#include "wallDistAddressing.H"
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace LESModels
39{
42}
43}
44
45// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46
47void Foam::LESModels::vanDriestDelta::calcDelta()
48{
50
53 const volScalarField& nu = tnu();
55
56 volScalarField ystar
57 (
59 (
60 "ystar",
61 mesh.time().constant(),
62 mesh.thisDb(),
64 ),
65 mesh,
66 dimensionedScalar("ystar", dimLength, GREAT)
67 );
68
69 const fvPatchList& patches = mesh.boundary();
70 volScalarField::Boundary& ystarBf = ystar.boundaryFieldRef();
71
72 forAll(patches, patchi)
73 {
74 if (isA<wallFvPatch>(patches[patchi]))
75 {
76 const fvPatchVectorField& Uw = U.boundaryField()[patchi];
77 const scalarField& nuw = nu.boundaryField()[patchi];
78 const scalarField& nuSgsw = nuSgs().boundaryField()[patchi];
79
80 ystarBf[patchi] =
81 nuw/sqrt((nuw + nuSgsw)*mag(Uw.snGrad()) + VSMALL);
82 }
83 }
84
85 // Construct wall transporter
86 const auto& wDist = wallDistAddressing::New(mesh);
87
88 // Get distance to nearest wall
89 const auto& y = wDist.y();
90
91 // Get ystar from nearest wall
92 wDist.map(ystar, mapDistribute::transform());
93
94 // Calculate y/ystar (stored in ystar!) and do the clipping
95 constexpr scalar yPlusCutOff = 500;
96 // Allow for some precision loss from transformation/interpolation of GREAT
97 // (= unvisited value)(though ystar is scalar so should not be transformed)
98 constexpr scalar fuzzyGREAT = 0.5*GREAT;
99
100 ystar.dimensions().reset(y.dimensions()/ystar.dimensions());
101 forAll(y, celli)
102 {
103 const scalar yPlus = y[celli]/ystar[celli];
104 if (y[celli] > fuzzyGREAT || (yPlus > yPlusCutOff))
105 {
106 // unvisited : y is GREAT, ystar is 1.0
107 ystar[celli] = GREAT;
108 }
109 else
110 {
111 ystar[celli] = yPlus;
112 }
113 }
114
115 forAll(y.boundaryField(), patchi)
116 {
117 const auto& yp = y.boundaryField()[patchi];
118 auto& ystarp = ystar.boundaryFieldRef()[patchi];
119
120 forAll(yp, i)
121 {
122 const scalar yPlus = yp[i]/ystarp[i];
123 if (yp[i] > fuzzyGREAT || (yPlus > yPlusCutOff))
124 {
125 ystarp[i] = GREAT;
126 }
127 else
128 {
129 ystarp[i] = yPlus;
130 }
131 }
132 }
133 ystar.correctBoundaryConditions();
134
135 // Note: y/ystar stored in ystar
136 delta_ = min
137 (
138 static_cast<const volScalarField&>(geometricDelta_()),
139 //(kappa_/Cdelta_)*((scalar(1) + SMALL) - exp(-y/ystar/Aplus_))*y
140 (kappa_/Cdelta_)*((scalar(1) + SMALL) - exp(-ystar/Aplus_))*y
141 );
142}
143
144
145// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
146
147Foam::LESModels::vanDriestDelta::vanDriestDelta
148(
149 const word& name,
151 const dictionary& dict
152)
153:
155 geometricDelta_
156 (
158 (
159 IOobject::groupName("geometricDelta", turbulence.U().group()),
161 // Note: cannot use optionalSubDict - if no *Coeffs dict the
162 // code will get stuck in a loop attempting to read the delta entry
163 // - consider looking up "geometricDelta" instead of "delta"?
164 dict.subDict(type() + "Coeffs")
165 )
166 ),
167 kappa_(dict.getOrDefault<scalar>("kappa", 0.41)),
168 Aplus_
169 (
170 dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
171 (
172 "Aplus",
173 26.0
174 )
175 ),
176 Cdelta_
177 (
178 dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
179 (
180 "Cdelta",
181 0.158
182 )
183 )
184{
185 calcInterval_ = 1;
186 const dictionary& coeffsDict = dict.optionalSubDict(type() + "Coeffs");
187 if (!coeffsDict.readIfPresent("calcInterval", calcInterval_))
188 {
189 coeffsDict.readIfPresent("updateInterval", calcInterval_);
190 }
192 delta_ = geometricDelta_();
193}
194
195
196// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
197
199{
200 const dictionary& coeffsDict = dict.optionalSubDict(type() + "Coeffs");
201
202 geometricDelta_().read(coeffsDict);
203 dict.readIfPresent<scalar>("kappa", kappa_);
204 coeffsDict.readIfPresent<scalar>("Aplus", Aplus_);
205 coeffsDict.readIfPresent<scalar>("Cdelta", Cdelta_);
206 calcInterval_ = 1;
207 if (!coeffsDict.readIfPresent<label>("calcInterval", calcInterval_))
208 {
209 coeffsDict.readIfPresent("updateInterval", calcInterval_);
210 }
211
212 calcDelta();
213}
214
215
217{
218 if
219 (
220 (calcInterval_ > 0)
221 && (turbulenceModel_.mesh().time().timeIndex() % calcInterval_) == 0
222 )
223 {
224 geometricDelta_().correct();
225 calcDelta();
226 }
227}
228
229
230// ************************************************************************* //
scalar y
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
@ NO_REGISTER
Do not request registration (bool: false).
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
virtual void read(const dictionary &)
Read the LESdelta dictionary.
Abstract base class for LES deltas.
Definition LESdelta.H:50
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.
Definition LESdelta.C:61
const turbulenceModel & turbulence() const
Return turbulenceModel reference.
Definition LESdelta.H:146
volScalarField delta_
Definition LESdelta.H:57
const turbulenceModel & turbulenceModel_
Definition LESdelta.H:55
static FOAM_NO_DANGLING_REFERENCE const wallDistAddressing & New(const fvMesh &mesh, Args &&... args)
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
A class for managing temporary objects.
Definition tmp.H:75
Abstract base class for turbulence models (RAS, LES and laminar).
const volVectorField & U() const
Access function to velocity field.
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
const fvMesh & mesh() const
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
scalar yPlus
Namespace for LES SGS models.
Namespace for OpenFOAM.
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.
dimensionedScalar exp(const dimensionedScalar &ds)
GeometricField< vector, fvPatchField, volMesh > volVectorField
PtrList< fvPatch > fvPatchList
Store lists of fvPatch as a PtrList.
Definition fvPatch.H:64
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar sqrt(const dimensionedScalar &ds)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fvPatchField< vector > fvPatchVectorField
volScalarField & nu
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299