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-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
35 namespace Foam
36 {
37 namespace RASModels
38 {
39 
40 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41 
42 template<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 
55 template<class BasicTurbulenceModel>
57 (
58  const alphaField& alpha,
59  const rhoField& rho,
60  const volVectorField& U,
61  const surfaceScalarField& alphaRhoPhi,
62  const surfaceScalarField& phi,
63  const transportModel& transport,
64  const word& propertiesName,
65  const word& type
66 )
67 :
69  (
70  type,
71  alpha,
72  rho,
73  U,
74  alphaRhoPhi,
75  phi,
76  transport,
77  propertiesName
78  ),
79 
80  Cmu_
81  (
83  (
84  "betaStar",
85  this->coeffDict_,
86  0.09
87  )
88  ),
89  beta_
90  (
92  (
93  "beta",
94  this->coeffDict_,
95  0.072
96  )
97  ),
98  gamma_
99  (
101  (
102  "gamma",
103  this->coeffDict_,
104  0.52
105  )
106  ),
107  alphaK_
108  (
110  (
111  "alphaK",
112  this->coeffDict_,
113  0.5
114  )
115  ),
116  alphaOmega_
117  (
119  (
120  "alphaOmega",
121  this->coeffDict_,
122  0.5
123  )
124  ),
125 
126  k_
127  (
128  IOobject
129  (
130  IOobject::groupName("k", alphaRhoPhi.group()),
131  this->runTime_.timeName(),
132  this->mesh_,
135  ),
136  this->mesh_
137  ),
138  omega_
139  (
140  IOobject
141  (
142  IOobject::groupName("omega", alphaRhoPhi.group()),
143  this->runTime_.timeName(),
144  this->mesh_,
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);
157  }
158 }
159 
160 
161 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
162 
163 template<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;
175  }
176 
177  return false;
178 }
179 
180 
181 template<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 
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() && dev(twoSymm(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 
215  // Turbulence specific dissipation rate equation
216  tmp<fvScalarMatrix> omegaEqn
217  (
218  fvm::ddt(alpha, rho, omega_)
219  + fvm::div(alphaRhoPhi, omega_)
220  - fvm::laplacian(alpha*rho*DomegaEff(), omega_)
221  ==
222  gamma_*alpha()*rho()*GbyNu
223  - fvm::SuSp(((2.0/3.0)*gamma_)*alpha()*rho()*divU, omega_)
224  - fvm::Sp(beta_*alpha()*rho()*omega_(), omega_)
225  + fvOptions(alpha, rho, omega_)
226  );
227 
228  omegaEqn.ref().relax();
229  fvOptions.constrain(omegaEqn.ref());
230  omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
231  solve(omegaEqn);
232  fvOptions.correct(omega_);
233  bound(omega_, this->omegaMin_);
234 
235 
236  // Turbulent kinetic energy equation
238  (
239  fvm::ddt(alpha, rho, k_)
240  + fvm::div(alphaRhoPhi, k_)
241  - fvm::laplacian(alpha*rho*DkEff(), k_)
242  ==
243  alpha()*rho()*G
244  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
245  - fvm::Sp(Cmu_*alpha()*rho()*omega_(), k_)
246  + fvOptions(alpha, rho, k_)
247  );
248 
249  kEqn.ref().relax();
250  fvOptions.constrain(kEqn.ref());
251  solve(kEqn);
252  fvOptions.correct(k_);
253  bound(k_, this->kMin_);
254 
255  correctNut();
256 }
257 
258 
259 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
260 
261 } // End namespace RASModels
262 } // End namespace Foam
263 
264 // ************************************************************************* //
Foam::RASModels::kOmega::read
virtual bool read()
Read RASProperties dictionary.
Definition: kOmega.C:164
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::RASModels::kOmega::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: kOmega.H:109
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
fvOptions.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fv::optionList::correct
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
Definition: fvOptionListTemplates.C:355
Foam::tmp::clear
void clear() const noexcept
Definition: tmpI.H:287
Foam::RASModels::kOmega::correctNut
virtual void correctNut()
Definition: kOmega.C:43
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::RASModel
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:52
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:103
Foam::RASModels::kOmega::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: kOmega.H:110
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:35
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
Foam::RASModels::kOmega::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmega.C:182
rho
rho
Definition: readInitialConditions.H:88
Foam::RASModels::kOmega::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: kOmega.H:111
Foam::fvm::SuSp
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
bound.H
Bound the given scalar field if it has gone unbounded.
correct
fvOptions correct(rho)
Foam::solve
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
divU
zeroField divU
Definition: alphaSuSp.H:3
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::fv::options
Finite-volume options.
Definition: fvOptions.H:55
Foam::RASModels::kOmega::kOmega
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:57
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
U
U
Definition: pEqn.H:72
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
kOmega.H
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::eddyViscosity
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:55
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
nut
scalar nut
Definition: evaluateNearWall.H:5
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:95
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:190
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::IOobject::MUST_READ
Definition: IOobject.H:185
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106