LienLeschziner.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-2015 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 "LienLeschziner.H"
27 #include "wallDist.H"
28 #include "bound.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace incompressible
36 {
37 namespace RASModels
38 {
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 defineTypeNameAndDebug(LienLeschziner, 0);
43 addToRunTimeSelectionTable(RASModel, LienLeschziner, dictionary);
44 
45 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
46 
48 {
49  const volScalarField yStar(sqrt(k_)*y_/nu());
50 
51  return
52  (scalar(1) - exp(-Anu_*yStar))
53  /((scalar(1) + SMALL) - exp(-Aeps_*yStar));
54 }
55 
56 
58 {
60 
61  return scalar(1) - 0.3*exp(-sqr(Rt));
62 }
63 
64 
66 {
67  const volScalarField yStar(sqrt(k_)*y_/nu());
68  const volScalarField le
69  (
70  kappa_*y_*((scalar(1) + SMALL) - exp(-Aeps_*yStar))
71  );
72 
73  return
74  (Ceps2_*pow(Cmu_, 0.75))
75  *(f2*sqrt(k_)*epsilon_/le)*exp(-AE_*sqr(yStar));
76 }
77 
78 
80 {
81  nut_ = Cmu_*fMu()*sqr(k_)/epsilon_;
83 }
84 
85 
86 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
87 
89 (
90  const geometricOneField& alpha,
91  const geometricOneField& rho,
92  const volVectorField& U,
93  const surfaceScalarField& alphaRhoPhi,
94  const surfaceScalarField& phi,
95  const transportModel& transport,
96  const word& propertiesName,
97  const word& type
98 )
99 :
101  (
102  type,
103  alpha,
104  rho,
105  U,
106  alphaRhoPhi,
107  phi,
108  transport,
109  propertiesName
110  ),
111 
112  Ceps1_
113  (
115  (
116  "Ceps1",
117  coeffDict_,
118  1.44
119  )
120  ),
121  Ceps2_
122  (
124  (
125  "Ceps2",
126  coeffDict_,
127  1.92
128  )
129  ),
130  sigmak_
131  (
133  (
134  "sigmak",
135  coeffDict_,
136  1.0
137  )
138  ),
139  sigmaEps_
140  (
142  (
143  "sigmaEps",
144  coeffDict_,
145  1.3
146  )
147  ),
148  Cmu_
149  (
151  (
152  "Cmu",
153  coeffDict_,
154  0.09
155  )
156  ),
157  kappa_
158  (
160  (
161  "kappa",
162  coeffDict_,
163  0.41
164  )
165  ),
166  Anu_
167  (
169  (
170  "Anu",
171  coeffDict_,
172  0.016
173  )
174  ),
175  Aeps_
176  (
178  (
179  "Aeps",
180  coeffDict_,
181  0.263
182  )
183  ),
184  AE_
185  (
187  (
188  "AE",
189  coeffDict_,
190  0.00222
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  y_(wallDist::New(mesh_).y())
221 {
222  bound(k_, kMin_);
223  bound(epsilon_, epsilonMin_);
224 
225  if (type == typeName)
226  {
227  correctNut();
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  Cmu_.readIfPresent(coeffDict());
244  kappa_.readIfPresent(coeffDict());
245  Anu_.readIfPresent(coeffDict());
246  Aeps_.readIfPresent(coeffDict());
247  AE_.readIfPresent(coeffDict());
248 
249  return true;
250  }
251  else
252  {
253  return false;
254  }
255 }
256 
257 
259 {
260  if (!turbulence_)
261  {
262  return;
263  }
264 
266 
267  tmp<volTensorField> tgradU = fvc::grad(U_);
269  (
270  GName(),
271  nut_*(tgradU() && twoSymm(tgradU()))
272  );
273  tgradU.clear();
274 
275  // Update epsilon and G at the wall
276  epsilon_.boundaryField().updateCoeffs();
277 
278  const volScalarField f2(this->f2());
279 
280  // Dissipation equation
281  tmp<fvScalarMatrix> epsEqn
282  (
284  + fvm::div(phi_, epsilon_)
286  ==
287  Ceps1_*G*epsilon_/k_
289  + E(f2)
290  );
291 
292  epsEqn().relax();
293  epsEqn().boundaryManipulate(epsilon_.boundaryField());
294  solve(epsEqn);
295  bound(epsilon_, epsilonMin_);
296 
297 
298  // Turbulent kinetic energy equation
300  (
301  fvm::ddt(k_)
302  + fvm::div(phi_, k_)
303  - fvm::laplacian(DkEff(), k_)
304  ==
305  G
306  - fvm::Sp(epsilon_/k_, k_)
307  );
308 
309  kEqn().relax();
310  solve(kEqn);
311  bound(k_, kMin_);
312 
313  correctNut();
314 }
315 
316 
317 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
318 
319 } // End namespace RASModels
320 } // End namespace incompressible
321 } // End namespace Foam
322 
323 // ************************************************************************* //
dimensionedScalar sqrt(const dimensionedScalar &ds)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Bound the given scalar field if it has gone unbounded.
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:174
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
static const wallDist & New(const fvMesh &mesh)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
A class for handling words, derived from string.
Definition: word.H:59
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
dimensionedScalar exp(const dimensionedScalar &ds)
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Namespace for OpenFOAM.
const volScalarField & y_
Wall distance.
static word groupName(Name name, const word &group)
Generic dimensioned Type class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual bool read()
Read RASProperties dictionary.
LienLeschziner(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.
Macros for easy insertion into run-time selection tables.
addToRunTimeSelectionTable(RASModel, kkLOmega, dictionary)
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:589
surfaceScalarField & phi
solverPerformance solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
scalar y
void correctBoundaryConditions()
Correct boundary field.
incompressible::RASModel::transportModel transportModel
Definition: eddyViscosity.H:75
volScalarField & nu
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
tmp< volScalarField > E(const volScalarField &f2) const
dimensionedSymmTensor sqr(const dimensionedVector &dv)
A class for managing temporary objects.
Definition: PtrList.H:118
RASModel< turbulenceModel > RASModel
const dimensionedScalar G
Newtonian constant of gravitation.
U
Definition: pEqn.H:82