ShihQuadraticKE.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "ShihQuadraticKE.H"
27 #include "bound.H"
28 #include "wallFvPatch.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace incompressible
37 {
38 namespace RASModels
39 {
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 defineTypeNameAndDebug(ShihQuadraticKE, 0);
44 addToRunTimeSelectionTable(RASModel, ShihQuadraticKE, dictionary);
45 
46 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
47 
49 {
51 }
52 
53 
55 {
56  volSymmTensorField S(symm(gradU));
57  volTensorField W(skew(gradU));
58 
59  volScalarField sBar((k_/epsilon_)*sqrt(2.0)*mag(S));
60  volScalarField wBar((k_/epsilon_)*sqrt(2.0)*mag(W));
61 
62  volScalarField Cmu((2.0/3.0)/(Cmu1_ + sBar + Cmu2_*wBar));
63 
64  nut_ = Cmu*sqr(k_)/epsilon_;
66 
68  k_*sqr(k_/epsilon_)/(Cbeta_ + pow3(sBar))
69  *(
70  Cbeta1_*dev(innerSqr(S))
71  + Cbeta2_*twoSymm(S&W)
72  + Cbeta3_*dev(symm(W&W))
73  );
74 }
75 
76 
77 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
78 
80 (
81  const geometricOneField& alpha,
82  const geometricOneField& rho,
83  const volVectorField& U,
84  const surfaceScalarField& alphaRhoPhi,
85  const surfaceScalarField& phi,
86  const transportModel& transport,
87  const word& type
88 )
89 :
91  (
92  type,
93  alpha,
94  rho,
95  U,
96  alphaRhoPhi,
97  phi,
98  transport
99  ),
100 
101  Ceps1_
102  (
104  (
105  "Ceps1",
106  coeffDict_,
107  1.44
108  )
109  ),
110  Ceps2_
111  (
113  (
114  "Ceps2",
115  coeffDict_,
116  1.92
117  )
118  ),
119  sigmak_
120  (
122  (
123  "sigmak",
124  coeffDict_,
125  1.0
126  )
127  ),
128  sigmaEps_
129  (
131  (
132  "sigmaEps",
133  coeffDict_,
134  1.3
135  )
136  ),
137  Cmu1_
138  (
140  (
141  "Cmu1",
142  coeffDict_,
143  1.25
144  )
145  ),
146  Cmu2_
147  (
149  (
150  "Cmu2",
151  coeffDict_,
152  0.9
153  )
154  ),
155  Cbeta_
156  (
158  (
159  "Cbeta",
160  coeffDict_,
161  1000.0
162  )
163  ),
164  Cbeta1_
165  (
167  (
168  "Cbeta1",
169  coeffDict_,
170  3.0
171  )
172  ),
173  Cbeta2_
174  (
176  (
177  "Cbeta2",
178  coeffDict_,
179  15.0
180  )
181  ),
182  Cbeta3_
183  (
185  (
186  "Cbeta3",
187  coeffDict_,
188  -19.0
189  )
190  ),
191 
192  k_
193  (
194  IOobject
195  (
196  IOobject::groupName("k", alphaRhoPhi.group()),
197  runTime_.timeName(),
198  mesh_,
201  ),
202  mesh_
203  ),
204 
205  epsilon_
206  (
207  IOobject
208  (
209  IOobject::groupName("epsilon", alphaRhoPhi.group()),
210  runTime_.timeName(),
211  mesh_,
214  ),
215  mesh_
216  )
217 {
218  bound(k_, kMin_);
219  bound(epsilon_, epsilonMin_);
220 
221  if (type == typeName)
222  {
223  printCoeffs(type);
224  }
225 }
226 
227 
228 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
229 
231 {
233  {
234  Ceps1_.readIfPresent(coeffDict());
235  Ceps2_.readIfPresent(coeffDict());
236  sigmak_.readIfPresent(coeffDict());
237  sigmaEps_.readIfPresent(coeffDict());
238  Cmu1_.readIfPresent(coeffDict());
239  Cmu2_.readIfPresent(coeffDict());
240  Cbeta_.readIfPresent(coeffDict());
241  Cbeta1_.readIfPresent(coeffDict());
242  Cbeta2_.readIfPresent(coeffDict());
243  Cbeta3_.readIfPresent(coeffDict());
244 
245  return true;
246  }
247  else
248  {
249  return false;
250  }
251 }
252 
253 
255 {
256  if (!turbulence_)
257  {
258  return;
259  }
260 
262 
263  tmp<volTensorField> tgradU = fvc::grad(U_);
264  const volTensorField& gradU = tgradU();
265 
267  (
268  GName(),
269  (nut_*twoSymm(gradU) - nonlinearStress_) && gradU
270  );
271 
272 
273  // Update epsilon and G at the wall
274  epsilon_.boundaryFieldRef().updateCoeffs();
275 
276  // Dissipation equation
277  tmp<fvScalarMatrix> epsEqn
278  (
280  + fvm::div(phi_, epsilon_)
282  ==
283  Ceps1_*G*epsilon_/k_
285  );
286 
287  epsEqn.ref().relax();
288  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
289  solve(epsEqn);
290  bound(epsilon_, epsilonMin_);
291 
292 
293  // Turbulent kinetic energy equation
295  (
296  fvm::ddt(k_)
297  + fvm::div(phi_, k_)
298  - fvm::laplacian(DkEff(), k_)
299  ==
300  G
301  - fvm::Sp(epsilon_/k_, k_)
302  );
303 
304  kEqn.ref().relax();
305  solve(kEqn);
306  bound(k_, kMin_);
307 
308 
309  // Re-calculate viscosity and non-linear stress
310  correctNonlinearStress(gradU);
311 }
312 
313 
314 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
315 
316 } // End namespace RASModels
317 } // End namespace incompressible
318 } // End namespace Foam
319 
320 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
dimensionedTensor skew(const dimensionedTensor &dt)
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
bool readIfPresent(const dictionary &)
Update the value of dimensioned<Type> if found in the dictionary.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Generic dimensioned Type class.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Macros for easy insertion into run-time selection tables.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
phi
Definition: pEqn.H:104
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
virtual bool read()
Read RASProperties dictionary.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
RASModel< momentumTransportModel > RASModel
ShihQuadraticKE(const geometricOneField &alpha, const geometricOneField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
dimensionedSymmTensor innerSqr(const dimensionedSymmTensor &dt)
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Bound the given scalar field if it has gone unbounded.
dimensionedScalar pow3(const dimensionedScalar &ds)
U
Definition: pEqn.H:72
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
addToRunTimeSelectionTable(RASModel, kkLOmega, dictionary)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
Base-class for all transport models used by the incompressible turbulence models. ...
const scalarList W(::W(thermo))
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
void correctBoundaryConditions()
Correct boundary field.
virtual void correctNonlinearStress(const volTensorField &gradU)
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
const dimensionedScalar & G
Newtonian constant of gravitation.
Namespace for OpenFOAM.