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-2021 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"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace incompressible
39 {
40 namespace RASModels
41 {
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 defineTypeNameAndDebug(ShihQuadraticKE, 0);
46 addToRunTimeSelectionTable(RASModel, ShihQuadraticKE, dictionary);
47 
48 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
49 
51 {
53 }
54 
55 
57 {
58  volSymmTensorField S(symm(gradU));
59  volTensorField W(skew(gradU));
60 
61  volScalarField sBar((k_/epsilon_)*sqrt(2.0)*mag(S));
62  volScalarField wBar((k_/epsilon_)*sqrt(2.0)*mag(W));
63 
64  volScalarField Cmu((2.0/3.0)/(Cmu1_ + sBar + Cmu2_*wBar));
65 
66  nut_ = Cmu*sqr(k_)/epsilon_;
68 
70  k_*sqr(k_/epsilon_)/(Cbeta_ + pow3(sBar))
71  *(
72  Cbeta1_*dev(innerSqr(S))
73  + Cbeta2_*twoSymm(S&W)
74  + Cbeta3_*dev(symm(W&W))
75  );
76 }
77 
78 
79 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
80 
82 (
83  const geometricOneField& alpha,
84  const geometricOneField& rho,
85  const volVectorField& U,
86  const surfaceScalarField& alphaRhoPhi,
87  const surfaceScalarField& phi,
88  const transportModel& transport,
89  const word& propertiesName,
90  const word& type
91 )
92 :
94  (
95  type,
96  alpha,
97  rho,
98  U,
99  alphaRhoPhi,
100  phi,
101  transport,
102  propertiesName
103  ),
104 
105  Ceps1_
106  (
108  (
109  "Ceps1",
110  coeffDict_,
111  1.44
112  )
113  ),
114  Ceps2_
115  (
117  (
118  "Ceps2",
119  coeffDict_,
120  1.92
121  )
122  ),
123  sigmak_
124  (
126  (
127  "sigmak",
128  coeffDict_,
129  1.0
130  )
131  ),
132  sigmaEps_
133  (
135  (
136  "sigmaEps",
137  coeffDict_,
138  1.3
139  )
140  ),
141  Cmu1_
142  (
144  (
145  "Cmu1",
146  coeffDict_,
147  1.25
148  )
149  ),
150  Cmu2_
151  (
153  (
154  "Cmu2",
155  coeffDict_,
156  0.9
157  )
158  ),
159  Cbeta_
160  (
162  (
163  "Cbeta",
164  coeffDict_,
165  1000.0
166  )
167  ),
168  Cbeta1_
169  (
171  (
172  "Cbeta1",
173  coeffDict_,
174  3.0
175  )
176  ),
177  Cbeta2_
178  (
180  (
181  "Cbeta2",
182  coeffDict_,
183  15.0
184  )
185  ),
186  Cbeta3_
187  (
189  (
190  "Cbeta3",
191  coeffDict_,
192  -19.0
193  )
194  ),
195 
196  k_
197  (
198  IOobject
199  (
200  IOobject::groupName("k", alphaRhoPhi.group()),
201  runTime_.timeName(),
202  mesh_,
205  ),
206  mesh_
207  ),
208 
209  epsilon_
210  (
211  IOobject
212  (
213  IOobject::groupName("epsilon", alphaRhoPhi.group()),
214  runTime_.timeName(),
215  mesh_,
218  ),
219  mesh_
220  )
221 {
222  bound(k_, kMin_);
223  bound(epsilon_, epsilonMin_);
224 
225  if (type == typeName)
226  {
227  printCoeffs(type);
228  }
229 }
230 
231 
232 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
233 
235 {
237  {
238  Ceps1_.readIfPresent(coeffDict());
239  Ceps2_.readIfPresent(coeffDict());
240  sigmak_.readIfPresent(coeffDict());
241  sigmaEps_.readIfPresent(coeffDict());
242  Cmu1_.readIfPresent(coeffDict());
243  Cmu2_.readIfPresent(coeffDict());
244  Cbeta_.readIfPresent(coeffDict());
245  Cbeta1_.readIfPresent(coeffDict());
246  Cbeta2_.readIfPresent(coeffDict());
247  Cbeta3_.readIfPresent(coeffDict());
248 
249  return true;
250  }
251 
252  return false;
253 }
254 
255 
257 {
258  if (!turbulence_)
259  {
260  return;
261  }
262 
264 
265  tmp<volTensorField> tgradU = fvc::grad(U_);
266  const volTensorField& gradU = tgradU();
267 
269  (
270  GName(),
271  (nut_*twoSymm(gradU) - nonlinearStress_) && gradU
272  );
273 
274 
275  // Update epsilon and G at the wall
276  epsilon_.boundaryFieldRef().updateCoeffs();
277 
278  // Dissipation equation
279  tmp<fvScalarMatrix> epsEqn
280  (
282  + fvm::div(phi_, epsilon_)
284  ==
287  );
288 
289  epsEqn.ref().relax();
290  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
291  solve(epsEqn);
292  bound(epsilon_, epsilonMin_);
293 
294 
295  // Turbulent kinetic energy equation
297  (
298  fvm::ddt(k_)
299  + fvm::div(phi_, k_)
300  - fvm::laplacian(DkEff(), k_)
301  ==
302  G
303  - fvm::Sp(epsilon_/k_, k_)
304  );
305 
306  kEqn.ref().relax();
307  solve(kEqn);
308  bound(k_, kMin_);
309 
310 
311  // Re-calculate viscosity and non-linear stress
312  correctNonlinearStress(gradU);
313 }
314 
315 
316 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
317 
318 } // End namespace RASModels
319 } // End namespace incompressible
320 } // End namespace Foam
321 
322 // ************************************************************************* //
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:169
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:194
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:65
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:82
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:56
Foam::incompressible::RASModels::ShihQuadraticKE::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: ShihQuadraticKE.C:234
Foam::nonlinearEddyViscosity< incompressible::RASModel >
Foam::incompressible::RASModels::ShihQuadraticKE::Cmu2_
dimensionedScalar Cmu2_
Definition: ShihQuadraticKE.H:88
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
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:42
Foam::incompressible::RASModels::ShihQuadraticKE::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: ShihQuadraticKE.C:256
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
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:50
Foam::IOobject::MUST_READ
Definition: IOobject.H:185
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106