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 | Copyright (C) 2011-2016 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& propertiesName,
88  const word& type
89 )
90 :
92  (
93  type,
94  alpha,
95  rho,
96  U,
97  alphaRhoPhi,
98  phi,
99  transport,
100  propertiesName
101  ),
102 
103  Ceps1_
104  (
106  (
107  "Ceps1",
108  coeffDict_,
109  1.44
110  )
111  ),
112  Ceps2_
113  (
115  (
116  "Ceps2",
117  coeffDict_,
118  1.92
119  )
120  ),
121  sigmak_
122  (
124  (
125  "sigmak",
126  coeffDict_,
127  1.0
128  )
129  ),
130  sigmaEps_
131  (
133  (
134  "sigmaEps",
135  coeffDict_,
136  1.3
137  )
138  ),
139  Cmu1_
140  (
142  (
143  "Cmu1",
144  coeffDict_,
145  1.25
146  )
147  ),
148  Cmu2_
149  (
151  (
152  "Cmu2",
153  coeffDict_,
154  0.9
155  )
156  ),
157  Cbeta_
158  (
160  (
161  "Cbeta",
162  coeffDict_,
163  1000.0
164  )
165  ),
166  Cbeta1_
167  (
169  (
170  "Cbeta1",
171  coeffDict_,
172  3.0
173  )
174  ),
175  Cbeta2_
176  (
178  (
179  "Cbeta2",
180  coeffDict_,
181  15.0
182  )
183  ),
184  Cbeta3_
185  (
187  (
188  "Cbeta3",
189  coeffDict_,
190  -19.0
191  )
192  ),
193 
194  k_
195  (
196  IOobject
197  (
198  IOobject::groupName("k", U.group()),
199  runTime_.timeName(),
200  mesh_,
203  ),
204  mesh_
205  ),
206 
207  epsilon_
208  (
209  IOobject
210  (
211  IOobject::groupName("epsilon", U.group()),
212  runTime_.timeName(),
213  mesh_,
216  ),
217  mesh_
218  )
219 {
220  bound(k_, kMin_);
221  bound(epsilon_, epsilonMin_);
222 
223  if (type == typeName)
224  {
225  printCoeffs(type);
226  }
227 }
228 
229 
230 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
231 
233 {
235  {
236  Ceps1_.readIfPresent(coeffDict());
237  Ceps2_.readIfPresent(coeffDict());
238  sigmak_.readIfPresent(coeffDict());
239  sigmaEps_.readIfPresent(coeffDict());
240  Cmu1_.readIfPresent(coeffDict());
241  Cmu2_.readIfPresent(coeffDict());
242  Cbeta_.readIfPresent(coeffDict());
243  Cbeta1_.readIfPresent(coeffDict());
244  Cbeta2_.readIfPresent(coeffDict());
245  Cbeta3_.readIfPresent(coeffDict());
246 
247  return true;
248  }
249  else
250  {
251  return false;
252  }
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  ==
285  Ceps1_*G*epsilon_/k_
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 // ************************************************************************* //
surfaceScalarField & phi
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
U
Definition: pEqn.H:83
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
dimensionedTensor skew(const dimensionedTensor &dt)
const dimensionedScalar G
Newtonian constant of gravitation.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
RASModel< turbulenceModel > RASModel
dimensionedScalar sqrt(const dimensionedScalar &ds)
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.
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)
incompressible::RASModel ::transportModel transportModel
Definition: eddyViscosity.H:75
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
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
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
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.
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:485
dimensionedScalar pow3(const dimensionedScalar &ds)
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
void correctBoundaryConditions()
Correct boundary field.
virtual void correctNonlinearStress(const volTensorField &gradU)
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Namespace for OpenFOAM.