ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / turbulenceModels / compressible / LES / SpalartAllmaras / SpalartAllmaras.C
blob5819741e391e94c62dd936b3517c250d533837f4
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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"
28 #include "wallDist.H"
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 namespace Foam
34 namespace compressible
36 namespace LESModels
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);
75     return
76         (scalar(1) + chi*fv1())
77        *(1/Cv2_)
78        *(3*(scalar(1) + chiByCv2) + sqr(chiByCv2))
79        /pow3(scalar(1) + chiByCv2);
83 tmp<volScalarField> SpalartAllmaras::fw(const volScalarField& Stilda) const
85     volScalarField r
86     (
87         min
88         (
89             nuTilda_
90            /(
91                max
92                (
93                    Stilda,
94                    dimensionedScalar("SMALL", Stilda.dimensions(), SMALL)
95                )
96               *sqr(kappa_*dTilda_)
97             ),
98             scalar(10.0)
99         )
100     );
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),
123     sigmaNut_
124     (
125         dimensioned<scalar>::lookupOrAddToDict
126         (
127             "sigmaNut",
128             coeffDict_,
129             0.66666
130         )
131     ),
132     Prt_
133     (
134         dimensioned<scalar>::lookupOrAddToDict
135         (
136             "Prt",
137             coeffDict_,
138             1.0
139         )
140     ),
142     Cb1_
143     (
144         dimensioned<scalar>::lookupOrAddToDict
145         (
146             "Cb1",
147             coeffDict_,
148             0.1355
149         )
150     ),
151     Cb2_
152     (
153         dimensioned<scalar>::lookupOrAddToDict
154         (
155             "Cb2",
156             coeffDict_,
157             0.622
158         )
159     ),
160     Cv1_
161     (
162         dimensioned<scalar>::lookupOrAddToDict
163         (
164             "Cv1",
165             coeffDict_,
166             7.1
167         )
168     ),
169     Cv2_
170     (
171         dimensioned<scalar>::lookupOrAddToDict
172         (
173             "Cv2",
174             coeffDict_,
175             5.0
176         )
177     ),
178     CDES_
179     (
180         dimensioned<scalar>::lookupOrAddToDict
181         (
182             "CDES",
183             coeffDict_,
184             0.65
185         )
186     ),
187     ck_
188     (
189         dimensioned<scalar>::lookupOrAddToDict
190         (
191             "ck",
192             coeffDict_,
193             0.07
194         )
195     ),
196     kappa_
197     (
198         dimensioned<scalar>::lookupOrAddToDict
199         (
200             "kappa",
201             *this,
202             0.41
203         )
204     ),
205     Cw1_(Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_),
206     Cw2_
207     (
208         dimensioned<scalar>::lookupOrAddToDict
209         (
210             "Cw2",
211             coeffDict_,
212             0.3
213         )
214     ),
215     Cw3_
216     (
217         dimensioned<scalar>::lookupOrAddToDict
218         (
219             "Cw3",
220             coeffDict_,
221             2.0
222         )
223     ),
225     nuTilda_
226     (
227         IOobject
228         (
229             "nuTilda",
230             runTime_.timeName(),
231             mesh_,
232             IOobject::MUST_READ,
233             IOobject::AUTO_WRITE
234         ),
235         mesh_
236     ),
238     dTilda_(min(CDES_*delta(), wallDist(mesh_).y())),
239     muSgs_
240     (
241         IOobject
242         (
243             "muSgs",
244             runTime_.timeName(),
245             mesh_,
246             IOobject::MUST_READ,
247             IOobject::AUTO_WRITE
248         ),
249         mesh_
250     ),
252     alphaSgs_
253     (
254         IOobject
255         (
256             "alphaSgs",
257             runTime_.timeName(),
258             mesh_,
259             IOobject::MUST_READ,
260             IOobject::AUTO_WRITE
261         ),
262         mesh_
263     )
265     updateSubGridScaleFields();
267     printCoeffs();
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
293     return
294     (
295       - fvm::laplacian(muEff(), U) - fvc::div(muEff()*dev2(T(fvc::grad(U))))
296     );
300 void SpalartAllmaras::correct(const tmp<volTensorField>& tgradU)
302     const volTensorField& gradU = tgradU();
303     LESModel::correct(gradU);
305     if (mesh_.changing())
306     {
307         dTilda_ = min(CDES_*delta(), wallDist(mesh_).y());
308     }
310     volScalarField Stilda
311     (
312         fv3()*::sqrt(2.0)*mag(skew(gradU)) + fv2()*nuTilda_/sqr(kappa_*dTilda_)
313     );
315     tmp<fvScalarMatrix> nuTildaEqn
316     (
317         fvm::ddt(rho(), nuTilda_)
318       + fvm::div(phi(), nuTilda_)
319       - fvm::laplacian
320         (
321             (nuTilda_*rho() + mu())/sigmaNut_,
322             nuTilda_,
323             "laplacian(DnuTildaEff,nuTilda)"
324         )
325       - rho()*Cb2_/sigmaNut_*magSqr(fvc::grad(nuTilda_))
326      ==
327         rho()*Cb1_*Stilda*nuTilda_
328       - fvm::Sp(rho()*Cw1_*fw(Stilda)*nuTilda_/sqr(dTilda_), nuTilda_)
329     );
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())
344     {
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);
358         return true;
359     }
360     else
361     {
362         return false;
363     }
367 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
369 } // End namespace LESModels
370 } // End namespace compressible
371 } // End namespace Foam
373 // ************************************************************************* //