BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / LES / SpalartAllmarasDDES / SpalartAllmarasDDES.C
blob06a83c1e8972afdc5afc580ad2def629783ad5a6
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 "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::S(const volTensorField& gradU) const
80     return sqrt(2.0)*mag(symm(gradU));
84 tmp<volScalarField> SpalartAllmarasDDES::dTilda(const volScalarField& S) const
86     return max
87     (
88         y_
89       - fd(S)
90        *max(y_ - CDES_*delta(), dimensionedScalar("zero", dimLength, 0)),
91         dimensionedScalar("small", dimLength, SMALL)
92     );
96 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
98 SpalartAllmarasDDES::SpalartAllmarasDDES
100     const volVectorField& U,
101     const surfaceScalarField& phi,
102     transportModel& transport,
103     const word& turbulenceModelName,
104     const word& modelName
107     SpalartAllmaras(U, phi, transport, turbulenceModelName, modelName)
111 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
113 } // End namespace LESModels
114 } // End namespace incompressible
115 } // End namespace Foam
117 // ************************************************************************* //