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 | 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 "LienLeschziner.H"
30 #include "wallDist.H"
31 #include "bound.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace incompressible
39 {
40 namespace RASModels
41 {
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 defineTypeNameAndDebug(LienLeschziner, 0);
46 addToRunTimeSelectionTable(RASModel, LienLeschziner, dictionary);
47 
48 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
49 
51 {
52  const volScalarField yStar(sqrt(k_)*y_/nu());
53 
54  return
55  (scalar(1) - exp(-Anu_*yStar))
56  /((scalar(1) + SMALL) - exp(-Aeps_*yStar));
57 }
58 
59 
61 {
63 
64  return scalar(1) - 0.3*exp(-sqr(Rt));
65 }
66 
67 
69 {
70  const volScalarField yStar(sqrt(k_)*y_/nu());
71  const volScalarField le
72  (
73  kappa_*y_*((scalar(1) + SMALL) - exp(-Aeps_*yStar))
74  );
75 
76  return
77  (Ceps2_*pow(Cmu_, 0.75))
78  *(f2*sqrt(k_)*epsilon_/le)*exp(-AE_*sqr(yStar));
79 }
80 
81 
83 {
84  nut_ = Cmu_*fMu()*sqr(k_)/epsilon_;
86 }
87 
88 
89 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
90 
92 (
93  const geometricOneField& alpha,
94  const geometricOneField& rho,
95  const volVectorField& U,
96  const surfaceScalarField& alphaRhoPhi,
97  const surfaceScalarField& phi,
98  const transportModel& transport,
99  const word& propertiesName,
100  const word& type
101 )
102 :
104  (
105  type,
106  alpha,
107  rho,
108  U,
109  alphaRhoPhi,
110  phi,
111  transport,
112  propertiesName
113  ),
114 
115  Ceps1_
116  (
118  (
119  "Ceps1",
120  coeffDict_,
121  1.44
122  )
123  ),
124  Ceps2_
125  (
127  (
128  "Ceps2",
129  coeffDict_,
130  1.92
131  )
132  ),
133  sigmak_
134  (
136  (
137  "sigmak",
138  coeffDict_,
139  1.0
140  )
141  ),
142  sigmaEps_
143  (
145  (
146  "sigmaEps",
147  coeffDict_,
148  1.3
149  )
150  ),
151  Cmu_
152  (
154  (
155  "Cmu",
156  coeffDict_,
157  0.09
158  )
159  ),
160  kappa_
161  (
163  (
164  "kappa",
165  coeffDict_,
166  0.41
167  )
168  ),
169  Anu_
170  (
172  (
173  "Anu",
174  coeffDict_,
175  0.016
176  )
177  ),
178  Aeps_
179  (
181  (
182  "Aeps",
183  coeffDict_,
184  0.263
185  )
186  ),
187  AE_
188  (
190  (
191  "AE",
192  coeffDict_,
193  0.00222
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  y_(wallDist::New(mesh_).y())
224 {
225  bound(k_, kMin_);
226  bound(epsilon_, epsilonMin_);
227 
228  if (type == typeName)
229  {
230  printCoeffs(type);
231  }
232 }
233 
234 
235 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
236 
238 {
240  {
241  Ceps1_.readIfPresent(coeffDict());
242  Ceps2_.readIfPresent(coeffDict());
243  sigmak_.readIfPresent(coeffDict());
244  sigmaEps_.readIfPresent(coeffDict());
245  Cmu_.readIfPresent(coeffDict());
246  kappa_.readIfPresent(coeffDict());
247  Anu_.readIfPresent(coeffDict());
248  Aeps_.readIfPresent(coeffDict());
249  AE_.readIfPresent(coeffDict());
250 
251  return true;
252  }
253 
254  return false;
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_.boundaryFieldRef().updateCoeffs();
277 
278  const volScalarField f2(this->f2());
279 
280  // Dissipation equation
281  tmp<fvScalarMatrix> epsEqn
282  (
284  + fvm::div(phi_, epsilon_)
286  ==
289  + E(f2)
290  );
291 
292  epsEqn.ref().relax();
293  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
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.ref().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 // ************************************************************************* //
wallDist.H
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::incompressible::RASModels::LienLeschziner::fMu
tmp< volScalarField > fMu() const
Definition: LienLeschziner.C:50
Foam::incompressible::RASModels::LienLeschziner::Ceps2_
dimensionedScalar Ceps2_
Definition: LienLeschziner.H:90
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::tmp::clear
void clear() const noexcept
Definition: tmpI.H:287
Foam::incompressible::RASModels::LienLeschziner::E
tmp< volScalarField > E(const volScalarField &f2) const
Definition: LienLeschziner.C:68
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::incompressible::RASModels::LienLeschziner::sigmak_
dimensionedScalar sigmak_
Definition: LienLeschziner.H:91
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:35
Foam::MeshObject< fvMesh, UpdateableMeshObject, wallDist >::New
static const wallDist & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::geometricOneField
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Definition: geometricOneField.H:55
Foam::incompressible::RASModels::LienLeschziner::LienLeschziner
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.
Definition: LienLeschziner.C:92
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::incompressible::RASModels::LienLeschziner::epsilon_
volScalarField epsilon_
Definition: LienLeschziner.H:105
rho
rho
Definition: readInitialConditions.H:88
Foam::incompressible::RASModels::LienLeschziner::sigmaEps_
dimensionedScalar sigmaEps_
Definition: LienLeschziner.H:92
LienLeschziner.H
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::incompressible::RASModels::LienLeschziner::Cmu_
dimensionedScalar Cmu_
Definition: LienLeschziner.H:94
Foam::incompressible::RASModels::LienLeschziner::DepsilonEff
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: LienLeschziner.H:161
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
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::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::incompressible::RASModels::LienLeschziner::kappa_
dimensionedScalar kappa_
Definition: LienLeschziner.H:95
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
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
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
U
U
Definition: pEqn.H:72
Foam::incompressible::RASModels::LienLeschziner::Aeps_
dimensionedScalar Aeps_
Definition: LienLeschziner.H:98
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::LienLeschziner::f2
tmp< volScalarField > f2() const
Definition: LienLeschziner.C:60
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::incompressible::RASModels::LienLeschziner::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: LienLeschziner.C:237
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::incompressible::RASModels::LienLeschziner::Anu_
dimensionedScalar Anu_
Definition: LienLeschziner.H:97
Foam::incompressible::RASModels::LienLeschziner::correctNut
virtual void correctNut()
Definition: LienLeschziner.C:82
Foam::incompressible::RASModels::LienLeschziner::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: LienLeschziner.C:258
Foam::eddyViscosity< incompressible::RASModel >
Foam::eddyViscosity::correct
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
Definition: eddyViscosity.C:134
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::LienLeschziner::Ceps1_
dimensionedScalar Ceps1_
Definition: LienLeschziner.H:89
Foam::GeometricField< scalar, 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::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:95
Foam::incompressible::RASModels::LienLeschziner::y_
const volScalarField & y_
Wall distance.
Definition: LienLeschziner.H:110
Foam::incompressible::RASModels::LienLeschziner::AE_
dimensionedScalar AE_
Definition: LienLeschziner.H:99
Foam::incompressible::RASModels::LienLeschziner::DkEff
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: LienLeschziner.H:152
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::IOobject::MUST_READ
Definition: IOobject.H:185
Foam::incompressible::RASModels::LienLeschziner::k_
volScalarField k_
Definition: LienLeschziner.H:104