Loading...
Searching...
No Matches
kOmegaSSTDES.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) 2022 Upstream CFD GmbH
10 Copyright (C) 2016-2022 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
30#include "kOmegaSSTDES.H"
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34namespace Foam
36namespace LESModels
37{
38
39// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
40
41template<class BasicTurbulenceModel>
43{
44 // Correct the turbulence viscosity
47 // Correct the turbulence thermal diffusivity
48 BasicTurbulenceModel::correctNut();
49}
50
51
52template<class BasicTurbulenceModel>
54{
55 correctNut(2*magSqr(symm(fvc::grad(this->U_))));
56}
57
58
59template<class BasicTurbulenceModel>
61(
62 const volScalarField& nur,
63 const volScalarField& magGradU
64) const
65{
66 const dimensionedScalar eps(magGradU.dimensions(), SMALL);
67
69 min(nur/(max(magGradU, eps)*sqr(this->kappa_*this->y_)), scalar(10));
70
71 tr.ref().boundaryFieldRef() == 0;
72
73 return tr;
74}
75
76
77template<class BasicTurbulenceModel>
79(
80 const volTensorField& gradU
81) const
82{
85
86 if (this->useSigma_)
87 {
88 volScalarField& S2 = tS2.ref();
89
90 const volScalarField& k = this->k_;
91 const volScalarField& omega = this->omega_;
92
93 const volScalarField CDkOmega
94 (
95 (2*this->alphaOmega2_)*(fvc::grad(k) & fvc::grad(omega))/omega
96 );
97
98 const volScalarField F1(this->F1(CDkOmega));
99
100 const volScalarField CDES(this->CDES(F1));
101 const volScalarField dTilda(this->dTilda(mag(gradU), CDES));
102 const volScalarField lengthScaleRAS(this->lengthScaleRAS());
103 const volScalarField Ssigma(this->Ssigma(gradU));
104
105 S2 =
106 pos(dTilda - lengthScaleRAS)*S2
107 + (scalar(1) - pos(dTilda - lengthScaleRAS))*sqr(Ssigma);
109
110 return tS2;
111}
112
113
114template<class BasicTurbulenceModel>
116(
117 const volScalarField& magGradU,
118 const volScalarField& CDES
119) const
120{
122}
123
124
125template<class BasicTurbulenceModel>
127(
128 const volScalarField& F1,
129 const volTensorField& gradU
130) const
132 volScalarField CDES(this->CDES(F1));
133 return sqrt(this->k_())/dTilda(mag(gradU), CDES)()();
134}
135
136
137template<class BasicTurbulenceModel>
139(
140 const volTensorField& gradU,
141 const volScalarField& S2
142) const
143{
144 if (this->useSigma_)
145 {
146 return S2();
147 }
149 return
151}
152
153
154template<class BasicTurbulenceModel>
156(
157 const volScalarField::Internal& GbyNu0,
160) const
161{
162 return GbyNu0; // Unlimited
163}
164
165
166// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
167
168template<class BasicTurbulenceModel>
169kOmegaSSTDES<BasicTurbulenceModel>::kOmegaSSTDES
170(
171 const alphaField& alpha,
172 const rhoField& rho,
173 const volVectorField& U,
174 const surfaceScalarField& alphaRhoPhi,
175 const surfaceScalarField& phi,
176 const transportModel& transport,
177 const word& propertiesName,
178 const word& type
179)
180:
181 kOmegaSSTBase<DESModel<BasicTurbulenceModel>>
182 (
183 type,
184 alpha,
185 rho,
186 U,
187 alphaRhoPhi,
188 phi,
189 transport,
190 propertiesName
191 ),
192
193 useSigma_
194 (
195 Switch::getOrAddToDict
196 (
197 "useSigma",
198 this->coeffDict_,
199 false
200 )
201 ),
202 kappa_
203 (
204 dimensioned<scalar>::getOrAddToDict
205 (
206 "kappa",
207 this->coeffDict_,
208 0.41
209 )
210 ),
211 CDESkom_
212 (
213 dimensioned<scalar>::getOrAddToDict
214 (
215 "CDESkom",
216 this->coeffDict_,
217 0.82
218 )
219 ),
220 CDESkeps_
221 (
222 dimensioned<scalar>::getOrAddToDict
223 (
224 "CDESkeps",
225 this->coeffDict_,
226 0.60
227 )
228 )
229{
230 // Note: Ctrans coeff is model specific; for k-w = 60
231 this->Ctrans_ =
232 dimensioned<scalar>::getOrAddToDict("Ctrans", this->coeffDict_, 60.0);
233
234 if (type == typeName)
235 {
237 << "This model is not recommended and will be deprecated in future "
238 << "releases. Please consider using the DDES/IDDES versions instead"
239 << endl;
240
241 this->printCoeffs(type);
243}
244
245
246// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
247
248template<class BasicTurbulenceModel>
250{
251 if (kOmegaSSTBase<DESModel<BasicTurbulenceModel>>::read())
252 {
253 useSigma_.readIfPresent("useSigma", this->coeffDict());
254 kappa_.readIfPresent(this->coeffDict());
255 CDESkom_.readIfPresent(this->coeffDict());
256 CDESkeps_.readIfPresent(this->coeffDict());
257
258 return true;
259 }
261 return false;
262}
263
264
265template<class BasicTurbulenceModel>
268{
269 const volScalarField& k = this->k_;
270 const volScalarField& omega = this->omega_;
272 return sqrt(k)/(this->betaStar_*omega);
273}
274
275
276template<class BasicTurbulenceModel>
279(
280 const volScalarField& CDES
281) const
282{
283 return CDES*this->delta();
284}
285
286
287template<class BasicTurbulenceModel>
289{
290 const volScalarField& k = this->k_;
291 const volScalarField& omega = this->omega_;
292 const volVectorField& U = this->U_;
293
294 const volScalarField CDkOmega
295 (
296 (2*this->alphaOmega2_)*(fvc::grad(k) & fvc::grad(omega))/omega
297 );
298
299 const volScalarField F1(this->F1(CDkOmega));
300
302 (
303 IOobject::scopedName("DES", "LESRegion"),
304 neg(dTilda(mag(fvc::grad(U)), CDES(F1)) - lengthScaleRAS())
305 );
306}
307
308
309// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
310
311} // End namespace LESModels
312} // End namespace Foam
313
314// ************************************************************************* //
label k
scalar delta
#define F2(B, C, D)
Definition SHA1.C:150
#define F1(B, C, D)
Definition SHA1.C:149
const dimensionSet & dimensions() const noexcept
Return dimensions.
DimensionedField< scalar, volMesh > Internal
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
Templated abstract base class for DES models.
Definition DESModel.H:55
virtual tmp< volScalarField > dTilda(const volScalarField &magGradU, const volScalarField &CDES) const
Return length scale.
BasicTurbulenceModel::alphaField alphaField
Switch useSigma_
Switch to activate grey-area enhanced sigma-(D)DES.
BasicTurbulenceModel::rhoField rhoField
virtual tmp< volScalarField > lengthScaleLES(const volScalarField &CDES) const
LES length scale.
virtual void correctNut(const volScalarField &S2)
virtual tmp< volScalarField > lengthScaleRAS() const
RAS length scale.
virtual tmp< volScalarField > S2(const volTensorField &gradU) const
Return square of strain rate.
virtual tmp< volScalarField > LESRegion() const
Return the LES field indicator.
virtual tmp< volScalarField::Internal > GbyNu(const volScalarField::Internal &GbyNu0, const volScalarField::Internal &F2, const volScalarField::Internal &S2) const
Return G/nu.
virtual tmp< volScalarField::Internal > GbyNu0(const volTensorField &gradU, const volScalarField &S2) const
Return (G/nu)_0.
virtual tmp< volScalarField > CDES(const volScalarField &F1) const
Blending for CDES parameter.
dimensionedScalar CDESkom_
DES coefficients.
BasicTurbulenceModel::transportModel transportModel
virtual tmp< volScalarField::Internal > epsilonByk(const volScalarField &F1, const volTensorField &gradU) const
Return epsilon/k.
tmp< volScalarField > r(const volScalarField &nur, const volScalarField &magGradU) const
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.
static dimensioned< Type > getOrAddToDict(const word &name, dictionary &dict, const dimensionSet &dims=dimless, const Type &deflt=Type(Zero))
Construct dimensioned from dictionary, with default value.
Base class implementation of the k-omega-SST turbulence model for incompressible and compressible flo...
virtual tmp< volScalarField > F1(const volScalarField &CDkOmega) const
virtual tmp< volScalarField > omega() const
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
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 WarningInFunction
Report a warning using Foam::Warning.
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.
dimensionedScalar pos(const dimensionedScalar &ds)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
GeometricField< tensor, fvPatchField, volMesh > volTensorField
dimensionedScalar sqrt(const dimensionedScalar &ds)
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
dimensionedScalar neg(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
volScalarField & alpha