Merge branch 'upstream/OpenFOAM' into master
[freefoam.git] / src / turbulenceModels / incompressible / LES / SpalartAllmarasDDES / SpalartAllmarasDDES.C
blobc7d8ea905a50352adbf837951d72e98be8dc5422
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2008-2009 OpenCFD Ltd.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "SpalartAllmarasDDES.H"
28 #include <OpenFOAM/addToRunTimeSelectionTable.H>
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 namespace Foam
34 namespace incompressible
36 namespace LESModels
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 defineTypeNameAndDebug(SpalartAllmarasDDES, 0);
42 addToRunTimeSelectionTable(LESModel, SpalartAllmarasDDES, dictionary);
44 // * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //
46 tmp<volScalarField> SpalartAllmarasDDES::rd
48     const volScalarField& visc,
49     const volScalarField& S
50 ) const
52     return min
53     (
54         visc
55         /(
56             max
57             (
58                 S,
59                 dimensionedScalar("SMALL", S.dimensions(), SMALL)
60             )*sqr(kappa_*y_)
61           + dimensionedScalar
62             (
63                 "ROOTVSMALL",
64                 dimensionSet(0, 2 , -1, 0, 0),
65                 ROOTVSMALL
66             )
67         ),
68         scalar(10)
69     );
73 tmp<volScalarField> SpalartAllmarasDDES::fd(const volScalarField& S) const
75     return 1 - tanh(pow3(8*rd(nuEff(), S)));
79 tmp<volScalarField> SpalartAllmarasDDES::dTilda(const volScalarField& S) const
81     return max
82     (
83         y_
84       - fd(S)
85        *max(y_ - CDES_*delta(), dimensionedScalar("zero", dimLength, 0)),
86         dimensionedScalar("small", dimLength, SMALL)
87     );
91 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
93 SpalartAllmarasDDES::SpalartAllmarasDDES
95     const volVectorField& U,
96     const surfaceScalarField& phi,
97     transportModel& transport
100     SpalartAllmaras(U, phi, transport, typeName)
104 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
106 } // End namespace LESModels
107 } // End namespace incompressible
108 } // End namespace Foam
110 // ************************ vim: set sw=4 sts=4 et: ************************ //