Loading...
Searching...
No Matches
kOmegaSSTIDDES.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2015 OpenFOAM Foundation
9 Copyright (C) 2015-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "kOmegaSSTIDDES.H"
30
31// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32
33namespace Foam
34{
35namespace LESModels
36{
37
38// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
39
40template<class BasicTurbulenceModel>
41const IDDESDelta& kOmegaSSTIDDES<BasicTurbulenceModel>::setDelta() const
42{
43 if (!isA<IDDESDelta>(this->delta_()))
44 {
46 << "The delta function must be set to a " << IDDESDelta::typeName
47 << " -based model" << exit(FatalError);
48 }
49
50 return refCast<const IDDESDelta>(this->delta_());
51}
52
53
54template<class BasicTurbulenceModel>
55tmp<volScalarField> kOmegaSSTIDDES<BasicTurbulenceModel>::alpha() const
56{
57 return max(0.25 - this->y_/IDDESDelta_.hmax(), scalar(-5));
58}
59
60
61template<class BasicTurbulenceModel>
62tmp<volScalarField> kOmegaSSTIDDES<BasicTurbulenceModel>::ft
63(
64 const volScalarField& magGradU
65) const
66{
67 return tanh(pow3(sqr(Ct_)*this->r(this->nut_, magGradU)));
68}
69
70
71template<class BasicTurbulenceModel>
72tmp<volScalarField> kOmegaSSTIDDES<BasicTurbulenceModel>::fl
73(
74 const volScalarField& magGradU
75) const
76{
77 return tanh(pow(sqr(Cl_)*this->r(this->nu(), magGradU), 10));
78}
79
80
81template<class BasicTurbulenceModel>
82tmp<volScalarField> kOmegaSSTIDDES<BasicTurbulenceModel>::fdt
83(
84 const volScalarField& magGradU
85) const
86{
87 return 1 - tanh(pow(Cdt1_*this->r(this->nut_, magGradU), Cdt2_));
88}
89
90
91// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
92
93template<class BasicTurbulenceModel>
95(
96 const volScalarField& magGradU,
97 const volScalarField& CDES
98) const
99{
100 const volScalarField lRAS(this->lengthScaleRAS());
101 const volScalarField lLES(this->lengthScaleLES(CDES));
102
103 const volScalarField alpha(this->alpha());
104 const volScalarField expTerm(exp(sqr(alpha)));
105
106 tmp<volScalarField> fB = min(2*pow(expTerm, -9.0), scalar(1));
107 const volScalarField fdTilda(max(1 - fdt(magGradU), fB));
108
109 if (fe_)
110 {
112 2*lerp(pow(expTerm, -9.), pow(expTerm, -11.09), pos0(alpha));
113 tmp<volScalarField> fe2 = 1 - max(ft(magGradU), fl(magGradU));
114 tmp<volScalarField> fe = max(fe1 - 1, scalar(0))*fe2;
115
116 // Original formulation from Shur et al. paper (2008)
117 return max
118 (
119 fdTilda*(1 + fe)*lRAS + (1 - fdTilda)*lLES,
121 );
122 }
123
124 // Simplified formulation from Gritskevich et al. paper (2011) where fe = 0
125 return max
126 (
127 lerp(lLES, lRAS, fdTilda),
129 );
130}
131
132
133// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
134
135template<class BasicTurbulenceModel>
136kOmegaSSTIDDES<BasicTurbulenceModel>::kOmegaSSTIDDES
137(
138 const alphaField& alpha,
139 const rhoField& rho,
140 const volVectorField& U,
141 const surfaceScalarField& alphaRhoPhi,
142 const surfaceScalarField& phi,
143 const transportModel& transport,
144 const word& propertiesName,
145 const word& type
146)
147:
148 kOmegaSSTDES<BasicTurbulenceModel>
149 (
150 alpha,
151 rho,
152 U,
153 alphaRhoPhi,
154 phi,
155 transport,
156 propertiesName,
157 type
158 ),
159
160 Cdt1_
161 (
162 dimensioned<scalar>::getOrAddToDict
163 (
164 "Cdt1",
165 this->coeffDict_,
166 20
167 )
168 ),
169 Cdt2_
170 (
171 dimensioned<scalar>::getOrAddToDict
172 (
173 "Cdt2",
174 this->coeffDict_,
175 3
176 )
177 ),
178 Cl_
179 (
180 dimensioned<scalar>::getOrAddToDict
181 (
182 "Cl",
183 this->coeffDict_,
184 5
185 )
186 ),
187 Ct_
188 (
189 dimensioned<scalar>::getOrAddToDict
190 (
191 "Ct",
192 this->coeffDict_,
193 1.87
194 )
195 ),
196 fe_
197 (
198 Switch::getOrAddToDict
199 (
200 "fe",
201 this->coeffDict_,
202 true
203 )
204 ),
205
206 IDDESDelta_(setDelta())
207{
208 if (type == typeName)
209 {
210 this->printCoeffs(type);
212}
213
214
215// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
216
217template<class BasicTurbulenceModel>
219{
221 {
222 Cdt1_.readIfPresent(this->coeffDict());
223 Cdt2_.readIfPresent(this->coeffDict());
224 Cl_.readIfPresent(this->coeffDict());
225 Ct_.readIfPresent(this->coeffDict());
226
227 return true;
229
230 return false;
231}
232
233
234template<class BasicTurbulenceModel>
236{
237 const volScalarField alpha(this->alpha());
238 const volScalarField expTerm(exp(sqr(alpha)));
239
240 tmp<volScalarField> fB = min(2*pow(expTerm, -9.0), scalar(1));
241 return max(1 - fdt(mag(fvc::grad(this->U_))), fB);
242}
243
244
245// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
246
247} // End namespace LESModels
248} // End namespace Foam
249
250// ************************************************************************* //
k-omega-SST DES turbulence model for incompressible and compressible flows.
virtual tmp< volScalarField > lengthScaleLES(const volScalarField &CDES) const
LES length scale.
virtual tmp< volScalarField > lengthScaleRAS() const
RAS length scale.
virtual tmp< volScalarField > CDES(const volScalarField &F1) const
Blending for CDES parameter.
tmp< volScalarField > r(const volScalarField &nur, const volScalarField &magGradU) const
virtual bool read()
Re-read model coefficients if they have changed.
virtual tmp< volScalarField > dTilda(const volScalarField &magGradU, const volScalarField &CDES) const
Return the length scale.
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::rhoField rhoField
virtual tmp< volScalarField > fd() const
Return the shielding function.
BasicTurbulenceModel::transportModel transportModel
const IDDESDelta & IDDESDelta_
IDDES delta.
virtual bool read()
Re-read model coefficients if they have changed.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
Generic dimensioned Type class.
A class for managing temporary objects.
Definition tmp.H:75
Base-class for all transport models used by the incompressible turbulence models.
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Namespace for LES SGS models.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
dimensionedScalar exp(const dimensionedScalar &ds)
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar tanh(const dimensionedScalar &ds)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
volScalarField & nu
volScalarField & alpha