Loading...
Searching...
No Matches
LienLeschziner.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) 2019-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 "LienLeschziner.H"
30#include "wallDist.H"
31#include "bound.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
40namespace RASModels
41{
42
43// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44
47
48// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
49
51{
52 const volScalarField yStar(sqrt(k_)*y_/nu());
54 return
55 (scalar(1) - exp(-Anu_*yStar))
56 /((scalar(1) + SMALL) - exp(-Aeps_*yStar));
57}
58
59
63
64 return scalar(1) - 0.3*exp(-sqr(Rt));
65}
66
67
69{
70 const volScalarField yStar(sqrt(k_)*y_/nu());
71 const volScalarField le
72 (
73 kappa_*y_*((scalar(1) + SMALL) - exp(-Aeps_*yStar))
74 );
76 return
77 (Ceps2_*pow(Cmu_, 0.75))
78 *(f2*sqrt(k_)*epsilon_/le)*exp(-AE_*sqr(yStar));
79}
80
81
83{
85 nut_.correctBoundaryConditions();
86}
87
88
89// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
90
92(
95 const volVectorField& U,
96 const surfaceScalarField& alphaRhoPhi,
98 const transportModel& transport,
99 const word& propertiesName,
100 const word& type
101)
102:
104 (
105 type,
106 alpha,
107 rho,
108 U,
109 alphaRhoPhi,
110 phi,
111 transport,
112 propertiesName
113 ),
114
115 Ceps1_
116 (
117 dimensioned<scalar>::getOrAddToDict
118 (
119 "Ceps1",
120 coeffDict_,
121 1.44
122 )
123 ),
124 Ceps2_
125 (
126 dimensioned<scalar>::getOrAddToDict
127 (
128 "Ceps2",
129 coeffDict_,
130 1.92
131 )
132 ),
133 sigmak_
134 (
135 dimensioned<scalar>::getOrAddToDict
136 (
137 "sigmak",
138 coeffDict_,
139 1.0
140 )
141 ),
142 sigmaEps_
143 (
144 dimensioned<scalar>::getOrAddToDict
145 (
146 "sigmaEps",
147 coeffDict_,
148 1.3
149 )
150 ),
151 Cmu_
152 (
153 dimensioned<scalar>::getOrAddToDict
154 (
155 "Cmu",
156 coeffDict_,
157 0.09
158 )
159 ),
160 kappa_
161 (
162 dimensioned<scalar>::getOrAddToDict
163 (
164 "kappa",
165 coeffDict_,
166 0.41
167 )
168 ),
169 Anu_
170 (
171 dimensioned<scalar>::getOrAddToDict
172 (
173 "Anu",
174 coeffDict_,
175 0.016
176 )
177 ),
178 Aeps_
179 (
180 dimensioned<scalar>::getOrAddToDict
181 (
182 "Aeps",
183 coeffDict_,
184 0.263
185 )
186 ),
187 AE_
188 (
189 dimensioned<scalar>::getOrAddToDict
190 (
191 "AE",
192 coeffDict_,
193 0.00222
194 )
195 ),
196
197 k_
198 (
200 (
201 IOobject::groupName("k", alphaRhoPhi.group()),
202 runTime_.timeName(),
203 mesh_,
204 IOobject::MUST_READ,
205 IOobject::AUTO_WRITE
206 ),
207 mesh_
208 ),
209
210 epsilon_
211 (
213 (
214 IOobject::groupName("epsilon", alphaRhoPhi.group()),
215 runTime_.timeName(),
216 mesh_,
217 IOobject::MUST_READ,
218 IOobject::AUTO_WRITE
219 ),
220 mesh_
221 ),
222
223 y_(wallDist::New(mesh_).y())
224{
225 bound(k_, kMin_);
226 bound(epsilon_, epsilonMin_);
227
228 if (type == typeName)
229 {
230 printCoeffs(type);
231 }
232}
233
234
235// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
236
238{
240 {
241 Ceps1_.readIfPresent(coeffDict());
242 Ceps2_.readIfPresent(coeffDict());
243 sigmak_.readIfPresent(coeffDict());
244 sigmaEps_.readIfPresent(coeffDict());
245 Cmu_.readIfPresent(coeffDict());
246 kappa_.readIfPresent(coeffDict());
247 Anu_.readIfPresent(coeffDict());
248 Aeps_.readIfPresent(coeffDict());
249 AE_.readIfPresent(coeffDict());
250
251 return true;
252 }
253
254 return false;
255}
256
257
259{
260 if (!turbulence_)
261 {
262 return;
263 }
264
266
267 tmp<volTensorField> tgradU = fvc::grad(U_);
269 (
270 GName(),
271 nut_*(tgradU() && twoSymm(tgradU()))
272 );
273 tgradU.clear();
274
275 // Update epsilon and G at the wall
277 // Push any changed cell values to coupled neighbours
279
280 const volScalarField f2(this->f2());
281
282 // Dissipation equation
284 (
286 + fvm::div(phi_, epsilon_)
288 ==
291 + E(f2)
292 );
293
294 epsEqn.ref().relax();
295 epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
296 solve(epsEqn);
297 bound(epsilon_, epsilonMin_);
298
299
300 // Turbulent kinetic energy equation
302 (
303 fvm::ddt(k_)
304 + fvm::div(phi_, k_)
306 ==
307 G
309 );
310
311 kEqn.ref().relax();
312 solve(kEqn);
313 bound(k_, kMin_);
314
315 correctNut();
316}
317
318
319// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
320
321} // End namespace RASModels
322} // End namespace incompressible
323} // End namespace Foam
324
325// ************************************************************************* //
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.
Bound the given scalar field if it has gone unbounded.
void evaluateCoupled(const UPstream::commsTypes commsType=UPstream::defaultCommsType)
Evaluate boundary conditions on coupled patches of the given type, using specified or default comms.
void updateCoeffs()
Update the boundary condition coefficients.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
An abstract base class for patches that couple regions of the computational domain e....
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Generic dimensioned Type class.
bool readIfPresent(const dictionary &dict)
Update the value of dimensioned<Type> if found in the dictionary, lookup in dictionary with the name(...
Eddy viscosity turbulence model base class.
eddyViscosity(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
virtual bool read()=0
Re-read model coefficients if they have changed.
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
incompressible::RASModel::transportModel transportModel
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Lien and Leschziner low-Reynolds number k-epsilon turbulence model for incompressible flows.
LienLeschziner(const geometricOneField &alpha, const geometricOneField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
const volScalarField & y_
Wall distance.
tmp< volScalarField > E(const volScalarField &f2) const
virtual bool read()
Re-read model coefficients if they have changed.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
A class for managing temporary objects.
Definition tmp.H:75
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition tmpI.H:289
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
Base-class for all transport models used by the incompressible turbulence models.
Interface to run-time selectable methods to calculate the distance-to-wall and normal-to-wall fields.
Definition wallDist.H:74
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
word timeName
Definition getTimeIndex.H:3
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvmDiv.C:41
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvmDdt.C:41
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
RASModel< turbulenceModel > RASModel
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition bound.C:29
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensionedScalar sqrt(const dimensionedScalar &ds)
volScalarField & nu
volScalarField & alpha
CEqn solve()