kOmegaSSTSAS.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-2017 OpenFOAM Foundation
9 Copyright (C) 2019-2020 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 "kOmegaSSTSAS.H"
30
31// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32
33namespace Foam
34{
35namespace RASModels
36{
37
38// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
39
40template<class BasicTurbulenceModel>
42(
46) const
47{
49 (
50 sqrt(this->k_())/(pow025(this->betaStar_)*this->omega_())
51 );
52
54 (
55 max
56 (
57 kappa_*sqrt(S2)
58 /(
59 mag(fvc::laplacian(this->U_))()()
61 (
62 "ROOTVSMALL",
63 dimensionSet(0, -1, -1, 0, 0),
64 ROOTVSMALL
65 )
66 ),
67 Cs_*sqrt(kappa_*zeta2_/(beta/this->betaStar_ - gamma))*delta()()
68 )
69 );
70
71 return fvm::Su
72 (
73 this->alpha_()*this->rho_()
74 *min
75 (
76 max
77 (
78 zeta2_*kappa_*S2*sqr(L/Lvk)
79 - (2*C_/sigmaPhi_)*this->k_()
80 *max
81 (
82 magSqr(fvc::grad(this->omega_)()())/sqr(this->omega_()),
83 magSqr(fvc::grad(this->k_)()())/sqr(this->k_())
84 ),
85 dimensionedScalar(dimensionSet(0, 0, -2, 0, 0), Zero)
86 ),
87 // Limit SAS production of omega for numerical stability,
88 // particularly during start-up
89 this->omega_()/(0.1*this->omega_.time().deltaT())
90 ),
91 this->omega_
92 );
93}
94
95
96// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
97
98template<class BasicTurbulenceModel>
100(
101 const alphaField& alpha,
102 const rhoField& rho,
103 const volVectorField& U,
104 const surfaceScalarField& alphaRhoPhi,
105 const surfaceScalarField& phi,
106 const transportModel& transport,
107 const word& propertiesName,
108 const word& type
109)
110:
111 kOmegaSST<BasicTurbulenceModel>
112 (
113 alpha,
114 rho,
115 U,
116 alphaRhoPhi,
117 phi,
118 transport,
119 propertiesName,
120 type
121 ),
122
123 Cs_
124 (
125 dimensioned<scalar>::getOrAddToDict
126 (
127 "Cs",
128 this->coeffDict_,
129 0.11
130 )
131 ),
132 kappa_
133 (
134 dimensioned<scalar>::getOrAddToDict
135 (
136 "kappa",
137 this->coeffDict_,
138 0.41
139 )
140 ),
141 zeta2_
142 (
143 dimensioned<scalar>::getOrAddToDict
144 (
145 "zeta2",
146 this->coeffDict_,
147 3.51
148 )
149 ),
150 sigmaPhi_
151 (
152 dimensioned<scalar>::getOrAddToDict
153 (
154 "sigmaPhi",
155 this->coeffDict_,
156 2.0/3.0
157 )
158 ),
159 C_
160 (
161 dimensioned<scalar>::getOrAddToDict
162 (
163 "C",
164 this->coeffDict_,
165 2
166 )
167 ),
168
169 delta_
170 (
172 (
173 IOobject::groupName("delta", alphaRhoPhi.group()),
174 *this,
175 this->coeffDict_
176 )
177 )
178{
179 if (type == typeName)
180 {
181 this->correctNut();
182 this->printCoeffs(type);
183 }
184}
185
186
187// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
188
189template<class BasicTurbulenceModel>
191{
193 {
194 Cs_.readIfPresent(this->coeffDict());
195 kappa_.readIfPresent(this->coeffDict());
196 sigmaPhi_.readIfPresent(this->coeffDict());
197 zeta2_.readIfPresent(this->coeffDict());
198 C_.readIfPresent(this->coeffDict());
199
200 return true;
201 }
202
203 return false;
204}
205
206
207// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
208
209} // End namespace RASModels
210} // End namespace Foam
211
212// ************************************************************************* //
scalar delta
surfaceScalarField & phi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Abstract base class for LES deltas.
Definition: LESdelta.H:54
Scale-adaptive URAS model based on the k-omega-SST RAS model.
Definition: kOmegaSSTSAS.H:104
BasicTurbulenceModel::alphaField alphaField
Definition: kOmegaSSTSAS.H:146
BasicTurbulenceModel::rhoField rhoField
Definition: kOmegaSSTSAS.H:147
virtual tmp< fvScalarMatrix > Qsas(const volScalarField::Internal &S2, const volScalarField::Internal &gamma, const volScalarField::Internal &beta) const
SAS omega source.
Definition: kOmegaSSTSAS.C:42
BasicTurbulenceModel::transportModel transportModel
Definition: kOmegaSSTSAS.H:148
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kOmegaSSTSAS.C:190
Implementation of the k-omega-SST turbulence model for incompressible and compressible flows.
Definition: kOmegaSST.H:132
virtual void correctNut()
Definition: kOmegaSST.C:54
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
const scalar gamma
Definition: EEqn.H:9
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:47
zeroField Su(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedSymmTensor sqr(const dimensionedVector &dv)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
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:33
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensionedScalar pow025(const dimensionedScalar &ds)
volScalarField & alpha
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
const vector L(dict.get< vector >("L"))