1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
11 OpenFOAM is free software: you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by
13 the Free Software Foundation, either version 3 of the License, or
14 (at your option) any later version.
16 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "SpalartAllmaras.H"
27 #include "addToRunTimeSelectionTable.H"
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 namespace compressible
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 defineTypeNameAndDebug(SpalartAllmaras, 0);
42 addToRunTimeSelectionTable(LESModel, SpalartAllmaras, dictionary);
45 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
47 void SpalartAllmaras::updateSubGridScaleFields()
49 muSgs_.internalField() = rho()*fv1()*nuTilda_.internalField();
50 muSgs_.correctBoundaryConditions();
52 alphaSgs_ = muSgs_/Prt_;
53 alphaSgs_.correctBoundaryConditions();
57 tmp<volScalarField> SpalartAllmaras::fv1() const
59 volScalarField chi3(pow3(rho()*nuTilda_/mu()));
60 return chi3/(chi3 + pow3(Cv1_));
64 tmp<volScalarField> SpalartAllmaras::fv2() const
66 return 1.0/pow3(scalar(1) + rho()*nuTilda_/(mu()*Cv2_));
70 tmp<volScalarField> SpalartAllmaras::fv3() const
72 volScalarField chi(rho()*nuTilda_/mu());
73 volScalarField chiByCv2((1/Cv2_)*chi);
76 (scalar(1) + chi*fv1())
78 *(3*(scalar(1) + chiByCv2) + sqr(chiByCv2))
79 /pow3(scalar(1) + chiByCv2);
83 tmp<volScalarField> SpalartAllmaras::fw(const volScalarField& Stilda) const
94 dimensionedScalar("SMALL", Stilda.dimensions(), SMALL)
101 r.boundaryField() == 0.0;
103 volScalarField g(r + Cw2_*(pow6(r) - r));
105 return g*pow((1.0 + pow6(Cw3_))/(pow6(g) + pow6(Cw3_)), 1.0/6.0);
109 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
111 SpalartAllmaras::SpalartAllmaras
113 const volScalarField& rho,
114 const volVectorField& U,
115 const surfaceScalarField& phi,
116 const basicThermo& thermoPhysicalModel,
117 const word& turbulenceModelName,
118 const word& modelName
121 LESModel(modelName, rho, U, phi, thermoPhysicalModel, turbulenceModelName),
125 dimensioned<scalar>::lookupOrAddToDict
134 dimensioned<scalar>::lookupOrAddToDict
144 dimensioned<scalar>::lookupOrAddToDict
153 dimensioned<scalar>::lookupOrAddToDict
162 dimensioned<scalar>::lookupOrAddToDict
171 dimensioned<scalar>::lookupOrAddToDict
180 dimensioned<scalar>::lookupOrAddToDict
189 dimensioned<scalar>::lookupOrAddToDict
198 dimensioned<scalar>::lookupOrAddToDict
205 Cw1_(Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_),
208 dimensioned<scalar>::lookupOrAddToDict
217 dimensioned<scalar>::lookupOrAddToDict
238 dTilda_(min(CDES_*delta(), wallDist(mesh_).y())),
265 updateSubGridScaleFields();
271 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
273 tmp<volSymmTensorField> SpalartAllmaras::B() const
275 return ((2.0/3.0)*I)*k() - (muSgs_/rho())*dev(twoSymm(fvc::grad(U())));
279 tmp<volSymmTensorField> SpalartAllmaras::devRhoBeff() const
281 return -muEff()*dev(twoSymm(fvc::grad(U())));
285 tmp<volScalarField> SpalartAllmaras::epsilon() const
287 return 2*muEff()/rho()*magSqr(symm(fvc::grad(U())));
291 tmp<fvVectorMatrix> SpalartAllmaras::divDevRhoBeff(volVectorField& U) const
295 - fvm::laplacian(muEff(), U) - fvc::div(muEff()*dev2(T(fvc::grad(U))))
300 void SpalartAllmaras::correct(const tmp<volTensorField>& tgradU)
302 const volTensorField& gradU = tgradU();
303 LESModel::correct(gradU);
305 if (mesh_.changing())
307 dTilda_ = min(CDES_*delta(), wallDist(mesh_).y());
310 volScalarField Stilda
312 fv3()*::sqrt(2.0)*mag(skew(gradU)) + fv2()*nuTilda_/sqr(kappa_*dTilda_)
315 tmp<fvScalarMatrix> nuTildaEqn
317 fvm::ddt(rho(), nuTilda_)
318 + fvm::div(phi(), nuTilda_)
321 (nuTilda_*rho() + mu())/sigmaNut_,
323 "laplacian(DnuTildaEff,nuTilda)"
325 - rho()*Cb2_/sigmaNut_*magSqr(fvc::grad(nuTilda_))
327 rho()*Cb1_*Stilda*nuTilda_
328 - fvm::Sp(rho()*Cw1_*fw(Stilda)*nuTilda_/sqr(dTilda_), nuTilda_)
331 nuTildaEqn().relax();
332 nuTildaEqn().solve();
334 bound(nuTilda_, dimensionedScalar("zero", nuTilda_.dimensions(), 0.0));
335 nuTilda_.correctBoundaryConditions();
337 updateSubGridScaleFields();
341 bool SpalartAllmaras::read()
343 if (LESModel::read())
345 sigmaNut_.readIfPresent(coeffDict());
346 Prt_.readIfPresent(coeffDict());
347 Cb1_.readIfPresent(coeffDict());
348 Cb2_.readIfPresent(coeffDict());
349 Cw1_ = Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_;
350 Cw2_.readIfPresent(coeffDict());
351 Cw3_.readIfPresent(coeffDict());
352 Cv1_.readIfPresent(coeffDict());
353 Cv2_.readIfPresent(coeffDict());
354 CDES_.readIfPresent(coeffDict());
355 ck_.readIfPresent(coeffDict());
356 kappa_.readIfPresent(*this);
367 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
369 } // End namespace LESModels
370 } // End namespace compressible
371 } // End namespace Foam
373 // ************************************************************************* //