ShihQuadraticKE.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 "ShihQuadraticKE.H"
30 #include "bound.H"
31 #include "wallFvPatch.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace incompressible
40 {
41 namespace RASModels
42 {
43 
44 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45 
46 defineTypeNameAndDebug(ShihQuadraticKE, 0);
47 addToRunTimeSelectionTable(RASModel, ShihQuadraticKE, dictionary);
48 
49 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
50 
52 {
54 }
55 
56 
58 {
59  volSymmTensorField S(symm(gradU));
60  volTensorField W(skew(gradU));
61 
62  volScalarField sBar((k_/epsilon_)*sqrt(2.0)*mag(S));
63  volScalarField wBar((k_/epsilon_)*sqrt(2.0)*mag(W));
64 
65  volScalarField Cmu((2.0/3.0)/(Cmu1_ + sBar + Cmu2_*wBar));
66 
67  nut_ = Cmu*sqr(k_)/epsilon_;
69 
71  k_*sqr(k_/epsilon_)/(Cbeta_ + pow3(sBar))
72  *(
73  Cbeta1_*dev(innerSqr(S))
74  + Cbeta2_*twoSymm(S&W)
75  + Cbeta3_*dev(symm(W&W))
76  );
77 }
78 
79 
80 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
81 
83 (
84  const geometricOneField& alpha,
85  const geometricOneField& rho,
86  const volVectorField& U,
87  const surfaceScalarField& alphaRhoPhi,
88  const surfaceScalarField& phi,
89  const transportModel& transport,
90  const word& propertiesName,
91  const word& type
92 )
93 :
95  (
96  type,
97  alpha,
98  rho,
99  U,
100  alphaRhoPhi,
101  phi,
102  transport,
103  propertiesName
104  ),
105 
106  Ceps1_
107  (
109  (
110  "Ceps1",
111  coeffDict_,
112  1.44
113  )
114  ),
115  Ceps2_
116  (
118  (
119  "Ceps2",
120  coeffDict_,
121  1.92
122  )
123  ),
124  sigmak_
125  (
127  (
128  "sigmak",
129  coeffDict_,
130  1.0
131  )
132  ),
133  sigmaEps_
134  (
136  (
137  "sigmaEps",
138  coeffDict_,
139  1.3
140  )
141  ),
142  Cmu1_
143  (
145  (
146  "Cmu1",
147  coeffDict_,
148  1.25
149  )
150  ),
151  Cmu2_
152  (
154  (
155  "Cmu2",
156  coeffDict_,
157  0.9
158  )
159  ),
160  Cbeta_
161  (
163  (
164  "Cbeta",
165  coeffDict_,
166  1000.0
167  )
168  ),
169  Cbeta1_
170  (
172  (
173  "Cbeta1",
174  coeffDict_,
175  3.0
176  )
177  ),
178  Cbeta2_
179  (
181  (
182  "Cbeta2",
183  coeffDict_,
184  15.0
185  )
186  ),
187  Cbeta3_
188  (
190  (
191  "Cbeta3",
192  coeffDict_,
193  -19.0
194  )
195  ),
196 
197  k_
198  (
199  IOobject
200  (
201  IOobject::groupName("k", alphaRhoPhi.group()),
202  runTime_.timeName(),
203  mesh_,
206  ),
207  mesh_
208  ),
209 
210  epsilon_
211  (
212  IOobject
213  (
214  IOobject::groupName("epsilon", alphaRhoPhi.group()),
215  runTime_.timeName(),
216  mesh_,
219  ),
220  mesh_
221  )
222 {
223  bound(k_, kMin_);
224  bound(epsilon_, epsilonMin_);
225 
226  if (type == typeName)
227  {
228  printCoeffs(type);
229  }
230 }
231 
232 
233 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
234 
236 {
238  {
239  Ceps1_.readIfPresent(coeffDict());
240  Ceps2_.readIfPresent(coeffDict());
241  sigmak_.readIfPresent(coeffDict());
242  sigmaEps_.readIfPresent(coeffDict());
243  Cmu1_.readIfPresent(coeffDict());
244  Cmu2_.readIfPresent(coeffDict());
245  Cbeta_.readIfPresent(coeffDict());
246  Cbeta1_.readIfPresent(coeffDict());
247  Cbeta2_.readIfPresent(coeffDict());
248  Cbeta3_.readIfPresent(coeffDict());
249 
250  return true;
251  }
252 
253  return false;
254 }
255 
256 
258 {
259  if (!turbulence_)
260  {
261  return;
262  }
263 
265 
266  tmp<volTensorField> tgradU = fvc::grad(U_);
267  const volTensorField& gradU = tgradU();
268 
270  (
271  GName(),
272  (nut_*twoSymm(gradU) - nonlinearStress_) && gradU
273  );
274 
275 
276  // Update epsilon and G at the wall
277  epsilon_.boundaryFieldRef().updateCoeffs();
278 
279  // Dissipation equation
280  tmp<fvScalarMatrix> epsEqn
281  (
283  + fvm::div(phi_, epsilon_)
285  ==
288  );
289 
290  epsEqn.ref().relax();
291  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
292  solve(epsEqn);
293  bound(epsilon_, epsilonMin_);
294 
295 
296  // Turbulent kinetic energy equation
298  (
299  fvm::ddt(k_)
300  + fvm::div(phi_, k_)
301  - fvm::laplacian(DkEff(), k_)
302  ==
303  G
304  - fvm::Sp(epsilon_/k_, k_)
305  );
306 
307  kEqn.ref().relax();
308  solve(kEqn);
309  bound(k_, kMin_);
310 
311 
312  // Re-calculate viscosity and non-linear stress
313  correctNonlinearStress(gradU);
314 }
315 
316 
317 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
318 
319 } // End namespace RASModels
320 } // End namespace incompressible
321 } // End namespace Foam
322 
323 // ************************************************************************* //
Foam::incompressible::RASModels::ShihQuadraticKE::k_
volScalarField k_
Definition: ShihQuadraticKE.H:97
Foam::incompressible::RASModels::ShihQuadraticKE::sigmaEps_
dimensionedScalar sigmaEps_
Definition: ShihQuadraticKE.H:86
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:84
Foam::incompressible::RASModels::ShihQuadraticKE::sigmak_
dimensionedScalar sigmak_
Definition: ShihQuadraticKE.H:85
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Foam::incompressible::RASModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(RASModel, kkLOmega, dictionary)
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::incompressible::RASModel
RASModel< turbulenceModel > RASModel
Definition: turbulentTransportModel.H:63
Foam::skew
dimensionedTensor skew(const dimensionedTensor &dt)
Definition: dimensionedTensor.C:138
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
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::incompressible::RASModels::ShihQuadraticKE::Cmu1_
dimensionedScalar Cmu1_
Definition: ShihQuadraticKE.H:87
Foam::geometricOneField
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Definition: geometricOneField.H:55
wallFvPatch.H
Foam::incompressible::RASModels::ShihQuadraticKE::ShihQuadraticKE
ShihQuadraticKE(const geometricOneField &alpha, const geometricOneField &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: ShihQuadraticKE.C:83
Foam::innerSqr
dimensionedSymmTensor innerSqr(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:62
rho
rho
Definition: readInitialConditions.H:88
Foam::incompressible::RASModels::ShihQuadraticKE::correctNonlinearStress
virtual void correctNonlinearStress(const volTensorField &gradU)
Definition: ShihQuadraticKE.C:57
Foam::incompressible::RASModels::ShihQuadraticKE::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: ShihQuadraticKE.C:235
Foam::nonlinearEddyViscosity< incompressible::RASModel >
Foam::incompressible::RASModels::ShihQuadraticKE::Cmu2_
dimensionedScalar Cmu2_
Definition: ShihQuadraticKE.H:88
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
Foam::incompressible::RASModels::ShihQuadraticKE::Cbeta_
dimensionedScalar Cbeta_
Definition: ShihQuadraticKE.H:89
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::dimensioned::readIfPresent
bool readIfPresent(const dictionary &dict)
Definition: dimensionedType.C:483
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
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.
Foam::solve
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
Foam::incompressible::RASModels::defineTypeNameAndDebug
defineTypeNameAndDebug(kkLOmega, 0)
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:43
Foam::incompressible::RASModels::ShihQuadraticKE::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: ShihQuadraticKE.C:257
Foam::incompressible::RASModels::ShihQuadraticKE::Ceps2_
dimensionedScalar Ceps2_
Definition: ShihQuadraticKE.H:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::incompressible::RASModels::ShihQuadraticKE::epsilon_
volScalarField epsilon_
Definition: ShihQuadraticKE.H:98
Foam::transportModel
Base-class for all transport models used by the incompressible turbulence models.
Definition: transportModel.H:53
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:940
Foam::incompressible::RASModels::ShihQuadraticKE::Cbeta1_
dimensionedScalar Cbeta1_
Definition: ShihQuadraticKE.H:90
Foam::incompressible::RASModels::ShihQuadraticKE::DepsilonEff
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: ShihQuadraticKE.H:148
U
U
Definition: pEqn.H:72
ShihQuadraticKE.H
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::incompressible::RASModels::ShihQuadraticKE::DkEff
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: ShihQuadraticKE.H:139
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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::correct
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
Definition: eddyViscosity.C:134
nutkWallFunctionFvPatchScalarField.H
Foam::incompressible::RASModels::ShihQuadraticKE::Ceps1_
dimensionedScalar Ceps1_
Definition: ShihQuadraticKE.H:83
Foam::incompressible::RASModels::ShihQuadraticKE::Cbeta3_
dimensionedScalar Cbeta3_
Definition: ShihQuadraticKE.H:92
Foam::eddyViscosity< incompressible::RASModel >::nut_
volScalarField nut_
Definition: eddyViscosity.H:66
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Foam::incompressible::RASModels::ShihQuadraticKE::Cbeta2_
dimensionedScalar Cbeta2_
Definition: ShihQuadraticKE.H:91
Foam::GeometricField< tensor, 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::nonlinearEddyViscosity< incompressible::RASModel >::nonlinearStress_
volSymmTensorField nonlinearStress_
Definition: nonlinearEddyViscosity.H:66
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:95
Foam::incompressible::RASModels::ShihQuadraticKE::correctNut
virtual void correctNut()
Definition: ShihQuadraticKE.C:51
Foam::IOobject::MUST_READ
Definition: IOobject.H:120
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106