Forward compatibility: flex
[foam-extend-3.2.git] / src / turbulenceModels / incompressible / LES / SpalartAllmarasDDES / SpalartAllmarasDDES.C
blobc9c60688b16150ce3b6946d76bae729d1c7d716e
1 /*---------------------------------------------------------------------------*\
2   =========                 |
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 -------------------------------------------------------------------------------
8 License
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 "SpalartAllmarasDDES.H"
27 #include "addToRunTimeSelectionTable.H"
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 namespace Foam
33 namespace incompressible
35 namespace LESModels
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 defineTypeNameAndDebug(SpalartAllmarasDDES, 0);
41 addToRunTimeSelectionTable(LESModel, SpalartAllmarasDDES, dictionary);
43 // * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //
45 tmp<volScalarField> SpalartAllmarasDDES::rd
47     const volScalarField& visc,
48     const volScalarField& S
49 ) const
51     return min
52     (
53         visc
54         /(
55             max
56             (
57                 S,
58                 dimensionedScalar("SMALL", S.dimensions(), SMALL)
59             )*sqr(kappa_*y_)
60           + dimensionedScalar
61             (
62                 "ROOTVSMALL",
63                 dimensionSet(0, 2 , -1, 0, 0),
64                 ROOTVSMALL
65             )
66         ),
67         scalar(10)
68     );
72 tmp<volScalarField> SpalartAllmarasDDES::fd(const volScalarField& S) const
74     return 1 - tanh(pow3(8*rd(nuEff(), S)));
78 tmp<volScalarField> SpalartAllmarasDDES::dTilda(const volScalarField& S) const
80     return max
81     (
82         y_
83       - fd(S)
84        *max(y_ - CDES_*delta(), dimensionedScalar("zero", dimLength, 0)),
85         dimensionedScalar("small", dimLength, SMALL)
86     );
90 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
92 SpalartAllmarasDDES::SpalartAllmarasDDES
94     const volVectorField& U,
95     const surfaceScalarField& phi,
96     transportModel& transport
99     SpalartAllmaras(U, phi, transport, typeName)
103 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
105 } // End namespace LESModels
106 } // End namespace incompressible
107 } // End namespace Foam
109 // ************************************************************************* //