1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 3.2
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
9 This file is part of foam-extend.
11 foam-extend is free software: you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by the
13 Free Software Foundation, either version 3 of the License, or (at your
14 option) any later version.
16 foam-extend is distributed in the hope that it will be useful, but
17 WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 General Public License for more details.
21 You should have received a copy of the GNU General Public License
22 along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "LienLeschzinerLowRe.H"
27 #include "wallFvPatch.H"
28 #include "addToRunTimeSelectionTable.H"
30 #include "backwardsCompatibilityWallFunctions.H"
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 namespace incompressible
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 defineTypeNameAndDebug(LienLeschzinerLowRe, 0);
44 addToRunTimeSelectionTable(RASModel, LienLeschzinerLowRe, dictionary);
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 LienLeschzinerLowRe::LienLeschzinerLowRe
50 const volVectorField& U,
51 const surfaceScalarField& phi,
52 transportModel& lamTransportModel
55 RASModel(typeName, U, phi, lamTransportModel),
59 dimensionedScalar::lookupOrAddToDict
68 dimensionedScalar::lookupOrAddToDict
77 dimensionedScalar::lookupOrAddToDict
86 dimensionedScalar::lookupOrAddToDict
95 dimensionedScalar::lookupOrAddToDict
104 dimensionedScalar::lookupOrAddToDict
113 dimensionedScalar::lookupOrAddToDict
122 dimensionedScalar::lookupOrAddToDict
131 dimensionedScalar::lookupOrAddToDict
167 yStar_(sqrt(k_)*y_/nu() + SMALL),
179 autoCreateLowReNut("nut", mesh_)
182 nut_ = Cmu_*(scalar(1) - exp(-Am_*yStar_))
183 /(scalar(1) - exp(-Aepsilon_*yStar_) + SMALL)*sqr(k_)
184 /(epsilon_ + epsilonSmall_);
185 nut_ = min(nut_, nuRatio()*nu());
186 nut_.correctBoundaryConditions();
192 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
194 tmp<volSymmTensorField> LienLeschzinerLowRe::R() const
196 volTensorField gradU = fvc::grad(U_);
198 return tmp<volSymmTensorField>
200 new volSymmTensorField
210 ((2.0/3.0)*I)*k_ - nut_*twoSymm(gradU),
211 k_.boundaryField().types()
217 tmp<volSymmTensorField> LienLeschzinerLowRe::devReff() const
219 return tmp<volSymmTensorField>
221 new volSymmTensorField
231 -nuEff()*dev(twoSymm(fvc::grad(U_)))
237 tmp<fvVectorMatrix> LienLeschzinerLowRe::divDevReff(volVectorField& U) const
241 - fvm::laplacian(nuEff(), U)
242 //- (fvc::grad(U) & fvc::grad(nuEff()))
243 - fvc::div(nuEff()*fvc::grad(U)().T())
248 bool LienLeschzinerLowRe::read()
250 if (RASModel::read())
252 C1_.readIfPresent(coeffDict());
253 C2_.readIfPresent(coeffDict());
254 sigmak_.readIfPresent(coeffDict());
255 sigmaEps_.readIfPresent(coeffDict());
256 Cmu_.readIfPresent(coeffDict());
257 kappa_.readIfPresent(coeffDict());
258 Am_.readIfPresent(coeffDict());
259 Aepsilon_.readIfPresent(coeffDict());
260 Amu_.readIfPresent(coeffDict());
271 void LienLeschzinerLowRe::correct()
273 // Bound in case of topological change
275 if (mesh_.changing())
278 bound(epsilon_, epsilon0_);
288 if (mesh_.changing())
293 scalar Cmu75 = pow(Cmu_.value(), 0.75);
295 volTensorField gradU = fvc::grad(U_);
298 volScalarField S2 = symm(gradU) && gradU;
300 yStar_ = sqrt(k_)*y_/nu() + SMALL;
301 volScalarField Rt = sqr(k_)/(nu()*epsilon_);
304 (scalar(1) - exp(-Am_*yStar_))
305 /(scalar(1) - exp(-Aepsilon_*yStar_) + SMALL);
307 volScalarField f2 = scalar(1) - 0.3*exp(-sqr(Rt));
309 volScalarField G("RASModel::G", Cmu_*fMu*sqr(k_)/epsilon_*S2);
312 // Dissipation equation
313 tmp<fvScalarMatrix> epsEqn
316 + fvm::div(phi_, epsilon_)
317 + fvm::SuSp(-fvc::div(phi_), epsilon_)
318 - fvm::laplacian(DepsilonEff(), epsilon_)
322 + C2_*f2*Cmu75*sqrt(k_)
323 /(kappa_*y_*(scalar(1) - exp(-Aepsilon_*yStar_)))
324 *exp(-Amu_*sqr(yStar_))*epsilon_
325 - fvm::Sp(C2_*f2*epsilon_/k_, epsilon_)
330 # include "LienLeschzinerLowReSetWallDissipation.H"
331 # include "wallDissipationI.H"
334 bound(epsilon_, epsilon0_);
337 // Turbulent kinetic energy equation
339 tmp<fvScalarMatrix> kEqn
343 + fvm::SuSp(-fvc::div(phi_), k_)
344 - fvm::laplacian(DkEff(), k_)
347 - fvm::Sp(epsilon_/k_, k_)
355 // Re-calculate viscosity
356 nut_ = Cmu_*fMu*sqr(k_)/epsilon_;
357 nut_ = min(nut_, nuRatio()*nu());
358 nut_.correctBoundaryConditions();
362 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
364 } // End namespace RASModels
365 } // End namespace incompressible
366 } // End namespace Foam
368 // ************************************************************************* //