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 "SpalartAllmarasIDDES.H"
27 #include "addToRunTimeSelectionTable.H"
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 namespace incompressible
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 defineTypeNameAndDebug(SpalartAllmarasIDDES, 0);
41 addToRunTimeSelectionTable(LESModel, SpalartAllmarasIDDES, dictionary);
43 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
45 tmp<volScalarField> SpalartAllmarasIDDES::alpha() const
49 0.25 - y_/static_cast<const volScalarField&>(hmax_()),
55 tmp<volScalarField> SpalartAllmarasIDDES::ft
57 const volScalarField& S
60 return tanh(pow3(sqr(ct_)*rd(nuSgs_, S)));
64 tmp<volScalarField> SpalartAllmarasIDDES::fl
66 const volScalarField& S
69 return tanh(pow(sqr(cl_)*rd(nu(), S), 10));
73 tmp<volScalarField> SpalartAllmarasIDDES::rd
75 const volScalarField& visc,
76 const volScalarField& S
86 dimensionedScalar("SMALL", S.dimensions(), SMALL)
91 dimensionSet(0, 2 , -1, 0, 0),
100 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
102 tmp<volScalarField> SpalartAllmarasIDDES::fd(const volScalarField& S) const
104 return 1 - tanh(pow3(8*rd(nuEff(), S)));
108 tmp<volScalarField> SpalartAllmarasIDDES::dTilda(const volScalarField& S) const
110 const volScalarField alpha(this->alpha());
111 const volScalarField expTerm(exp(sqr(alpha)));
113 tmp<volScalarField> fHill =
114 2*(pos(alpha)*pow(expTerm, -11.09) + neg(alpha)*pow(expTerm, -9.0));
116 tmp<volScalarField> fStep = min(2*pow(expTerm, -9.0), scalar(1));
117 const volScalarField fHyb(max(1 - fd(S), fStep));
118 tmp<volScalarField> fAmp = 1 - max(ft(S), fl(S));
119 tmp<volScalarField> fRestore = max(fHill - 1, scalar(0))*fAmp;
121 // IGNORING ft2 terms
122 const volScalarField Psi
129 (1 - Cb1_/(Cw1_*sqr(kappa_)*fwStar_)*fv2())/max(SMALL, fv1())
136 dimensionedScalar("SMALL", dimLength, SMALL),
137 fHyb*(1 + fRestore*Psi)*y_
138 + (1 - fHyb)*CDES_*Psi*delta()
143 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
145 SpalartAllmarasIDDES::SpalartAllmarasIDDES
147 const volVectorField& U,
148 const surfaceScalarField& phi,
149 transportModel& transport,
150 const word& turbulenceModelName,
151 const word& modelName
154 SpalartAllmaras(U, phi, transport, turbulenceModelName, modelName),
170 this->subDict(typeName + "Coeffs")
175 dimensioned<scalar>::lookupOrAddToDict
184 dimensioned<scalar>::lookupOrAddToDict
193 dimensioned<scalar>::lookupOrAddToDict
203 bool SpalartAllmarasIDDES::read()
205 if (SpalartAllmaras::read())
207 fwStar_.readIfPresent(coeffDict());
208 cl_.readIfPresent(coeffDict());
209 ct_.readIfPresent(coeffDict());
220 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
222 } // End namespace LESModels
223 } // End namespace incompressible
224 } // End namespace Foam
226 // ************************************************************************* //