Loading...
Searching...
No Matches
kOmega.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2019-2023 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 "kOmega.H"
30#include "fvOptions.H"
31#include "bound.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
37namespace RASModels
38{
39
40// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41
42template<class BasicTurbulenceModel>
44{
45 this->nut_ = k_/omega_;
46 this->nut_.correctBoundaryConditions();
47 fv::options::New(this->mesh_).correct(this->nut_);
48
49 BasicTurbulenceModel::correctNut();
50}
51
52
53// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54
55template<class BasicTurbulenceModel>
57(
58 const alphaField& alpha,
59 const rhoField& rho,
60 const volVectorField& U,
61 const surfaceScalarField& alphaRhoPhi,
63 const transportModel& transport,
64 const word& propertiesName,
65 const word& type
66)
67:
68 eddyViscosity<RASModel<BasicTurbulenceModel>>
69 (
70 type,
71 alpha,
72 rho,
73 U,
74 alphaRhoPhi,
75 phi,
76 transport,
77 propertiesName
78 ),
79
80 Cmu_
81 (
82 dimensioned<scalar>::getOrAddToDict
83 (
84 "betaStar",
85 this->coeffDict_,
86 0.09
87 )
88 ),
89 beta_
90 (
91 dimensioned<scalar>::getOrAddToDict
92 (
93 "beta",
94 this->coeffDict_,
95 0.072
96 )
97 ),
98 gamma_
99 (
100 dimensioned<scalar>::getOrAddToDict
101 (
102 "gamma",
103 this->coeffDict_,
104 0.52
105 )
106 ),
107 alphaK_
108 (
109 dimensioned<scalar>::getOrAddToDict
110 (
111 "alphaK",
112 this->coeffDict_,
113 0.5
114 )
115 ),
116 alphaOmega_
117 (
118 dimensioned<scalar>::getOrAddToDict
119 (
120 "alphaOmega",
121 this->coeffDict_,
122 0.5
123 )
124 ),
125
126 k_
127 (
129 (
130 IOobject::groupName("k", alphaRhoPhi.group()),
131 this->runTime_.timeName(),
132 this->mesh_,
133 IOobject::MUST_READ,
134 IOobject::AUTO_WRITE
135 ),
136 this->mesh_
137 ),
138 omega_
139 (
141 (
142 IOobject::groupName("omega", alphaRhoPhi.group()),
143 this->runTime_.timeName(),
144 this->mesh_,
145 IOobject::MUST_READ,
146 IOobject::AUTO_WRITE
147 ),
148 this->mesh_
149 )
150{
151 bound(k_, this->kMin_);
152 bound(omega_, this->omegaMin_);
153
154 if (type == typeName)
155 {
156 this->printCoeffs(type);
158}
159
160
161// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
162
163template<class BasicTurbulenceModel>
165{
167 {
168 Cmu_.readIfPresent(this->coeffDict());
169 beta_.readIfPresent(this->coeffDict());
170 gamma_.readIfPresent(this->coeffDict());
171 alphaK_.readIfPresent(this->coeffDict());
172 alphaOmega_.readIfPresent(this->coeffDict());
173
174 return true;
176
177 return false;
178}
179
180
181template<class BasicTurbulenceModel>
183{
184 if (!this->turbulence_)
185 {
186 return;
187 }
188
189 // Local references
190 const alphaField& alpha = this->alpha_;
191 const rhoField& rho = this->rho_;
192 const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
193 const volVectorField& U = this->U_;
194 const volScalarField& nut = this->nut_;
195 fv::options& fvOptions(fv::options::New(this->mesh_));
196
197 eddyViscosity<RASModel<BasicTurbulenceModel>>::correct();
198
200 (
201 fvc::div(fvc::absolute(this->phi(), U))().v()
202 );
203
204 tmp<volTensorField> tgradU = fvc::grad(U);
205 const volScalarField::Internal GbyNu
206 (
207 tgradU().v() && devTwoSymm(tgradU().v())
208 );
209 const volScalarField::Internal G(this->GName(), nut()*GbyNu);
210 tgradU.clear();
211
212 // Update omega and G at the wall
213 omega_.boundaryFieldRef().updateCoeffs();
214 // Push any changed cell values to coupled neighbours
215 omega_.boundaryFieldRef().template evaluateCoupled<coupledFvPatch>();
216
217 // Turbulence specific dissipation rate equation
218 tmp<fvScalarMatrix> omegaEqn
219 (
220 fvm::ddt(alpha, rho, omega_)
221 + fvm::div(alphaRhoPhi, omega_)
222 - fvm::laplacian(alpha*rho*DomegaEff(), omega_)
223 ==
224 gamma_*alpha()*rho()*GbyNu
225 - fvm::SuSp(((2.0/3.0)*gamma_)*alpha()*rho()*divU, omega_)
226 - fvm::Sp(beta_*alpha()*rho()*omega_(), omega_)
227 + fvOptions(alpha, rho, omega_)
228 );
229
230 omegaEqn.ref().relax();
231 fvOptions.constrain(omegaEqn.ref());
232 omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
233 solve(omegaEqn);
234 fvOptions.correct(omega_);
235 bound(omega_, this->omegaMin_);
236
237
238 // Turbulent kinetic energy equation
239 tmp<fvScalarMatrix> kEqn
240 (
241 fvm::ddt(alpha, rho, k_)
242 + fvm::div(alphaRhoPhi, k_)
243 - fvm::laplacian(alpha*rho*DkEff(), k_)
244 ==
245 alpha()*rho()*G
246 - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
247 - fvm::Sp(Cmu_*alpha()*rho()*omega_(), k_)
248 + fvOptions(alpha, rho, k_)
249 );
250
251 kEqn.ref().relax();
252 fvOptions.constrain(kEqn.ref());
253 solve(kEqn);
254 fvOptions.correct(k_);
255 bound(k_, this->kMin_);
256
257 correctNut();
258}
259
260
261// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
262
263} // End namespace RASModels
264} // End namespace Foam
265
266// ************************************************************************* //
Bound the given scalar field if it has gone unbounded.
fv::options & fvOptions
DimensionedField< scalar, volMesh > Internal
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Templated abstract base class for RAS turbulence models.
Definition RASModel.H:81
BasicTurbulenceModel::alphaField alphaField
Definition kOmega.H:107
BasicTurbulenceModel::rhoField rhoField
Definition kOmega.H:108
volScalarField k_
Definition kOmega.H:96
tmp< volScalarField > DomegaEff() const
Return the effective diffusivity for omega.
Definition kOmega.H:167
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition kOmega.C:175
dimensionedScalar alphaOmega_
Definition kOmega.H:91
kOmega(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
Definition kOmega.C:50
dimensionedScalar gamma_
Definition kOmega.H:89
volScalarField omega_
Definition kOmega.H:97
virtual void correctNut()
Definition kOmega.C:36
BasicTurbulenceModel::transportModel transportModel
Definition kOmega.H:109
dimensionedScalar Cmu_
Definition kOmega.H:87
dimensionedScalar alphaK_
Definition kOmega.H:90
virtual bool read()
Read RASProperties dictionary.
Definition kOmega.C:157
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition kOmega.H:152
dimensionedScalar beta_
Definition kOmega.H:88
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
eddyViscosity(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
Finite-volume options, which is an IOdictionary of values and a fv::optionList.
Definition fvOptions.H:69
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
Definition fvOptions.C:116
A class for managing temporary objects.
Definition tmp.H:75
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition tmpI.H:289
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
thermo correct()
scalar nut
zeroField divU
Definition alphaSuSp.H:3
word timeName
Definition getTimeIndex.H:3
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition fvcMeshPhi.C:183
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvmDiv.C:41
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvmDdt.C:41
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition bound.C:29
GeometricField< scalar, fvPatchField, volMesh > volScalarField
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
volScalarField & alpha
CEqn solve()