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  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  Cmu_.readIfPresent(coeffDict());
243  kappa_.readIfPresent(coeffDict());
244  Anu_.readIfPresent(coeffDict());
245  Aeps_.readIfPresent(coeffDict());
246  AE_.readIfPresent(coeffDict());
247 
248  return true;
249  }
250  else
251  {
252  return false;
253  }
254 }
255 
256 
258 {
259  if (!turbulence_)
260  {
261  return;
262  }
263 
265 
266  tmp<volTensorField> tgradU = fvc::grad(U_);
268  (
269  GName(),
270  nut_*(tgradU() && twoSymm(tgradU()))
271  );
272  tgradU.clear();
273 
274  // Update epsilon and G at the wall
276 
277  const volScalarField f2(this->f2());
278 
279  // Dissipation equation
280  tmp<fvScalarMatrix> epsEqn
281  (
283  + fvm::div(phi_, epsilon_)
285  ==
288  + E(f2)
289  );
290 
291  epsEqn().relax();
292  epsEqn().boundaryManipulate(epsilon_.boundaryField());
293  solve(epsEqn);
294  bound(epsilon_, epsilonMin_);
295 
296 
297  // Turbulent kinetic energy equation
299  (
300  fvm::ddt(k_)
301  + fvm::div(phi_, k_)
302  - fvm::laplacian(DkEff(), k_)
303  ==
304  G
305  - fvm::Sp(epsilon_/k_, k_)
306  );
307 
308  kEqn().relax();
309  solve(kEqn);
310  bound(k_, kMin_);
311 
312  correctNut();
313 }
314 
315 
316 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
317 
318 } // End namespace RASModels
319 } // End namespace incompressible
320 } // End namespace Foam
321 
322 // ************************************************************************* //
wallDist.H
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::incompressible::RASModels::LienLeschziner::fMu
tmp< volScalarField > fMu() const
Definition: LienLeschziner.C:47
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
Foam::incompressible::RASModels::LienLeschziner::Ceps2_
dimensionedScalar Ceps2_
Definition: LienLeschziner.H:87
Foam::GeometricField::GeometricBoundaryField::updateCoeffs
void updateCoeffs()
Update the boundary condition coefficients.
Definition: GeometricBoundaryField.C:447
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Foam::incompressible::RASModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(RASModel, kkLOmega, dictionary)
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::incompressible::RASModel
RASModel< turbulenceModel > RASModel
Definition: turbulentTransportModel.H:59
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::incompressible::RASModels::LienLeschziner::E
tmp< volScalarField > E(const volScalarField &f2) const
Definition: LienLeschziner.C:65
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::incompressible::RASModels::LienLeschziner::sigmak_
dimensionedScalar sigmak_
Definition: LienLeschziner.H:88
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
Foam::geometricOneField
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Definition: geometricOneField.H:52
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
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:89
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
Foam::incompressible::RASModels::LienLeschziner::epsilon_
volScalarField epsilon_
Definition: LienLeschziner.H:102
Foam::incompressible::RASModels::LienLeschziner::sigmaEps_
dimensionedScalar sigmaEps_
Definition: LienLeschziner.H:89
Foam::MeshObject< fvMesh, UpdateableMeshObject, wallDist >::New
static const wallDist & New(const fvMesh &mesh)
Definition: MeshObject.C:44
LienLeschziner.H
U
U
Definition: pEqn.H:46
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
nu
volScalarField & nu
Definition: readMechanicalProperties.H:179
Foam::incompressible::RASModels::LienLeschziner::Cmu_
dimensionedScalar Cmu_
Definition: LienLeschziner.H:91
Foam::incompressible::RASModels::LienLeschziner::DepsilonEff
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: LienLeschziner.H:159
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
bound.H
Bound the given scalar field if it has gone unbounded.
Foam::incompressible::RASModels::LienLeschziner::kappa_
dimensionedScalar kappa_
Definition: LienLeschziner.H:92
Foam::incompressible::RASModels::defineTypeNameAndDebug
defineTypeNameAndDebug(kkLOmega, 0)
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
Foam::eddyViscosity< incompressible::RASModel >::transportModel
incompressible::RASModel ::transportModel transportModel
Definition: eddyViscosity.H:75
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::dimensioned::readIfPresent
bool readIfPresent(const dictionary &)
Update the value of dimensioned<Type> if found in the dictionary.
Definition: dimensionedType.C:309
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:885
rho
rho
Definition: pEqn.H:3
Foam::incompressible::RASModels::LienLeschziner::Aeps_
dimensionedScalar Aeps_
Definition: LienLeschziner.H:95
Foam::solve
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::incompressible::RASModels::LienLeschziner::f2
tmp< volScalarField > f2() const
Definition: LienLeschziner.C:57
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::incompressible::RASModels::LienLeschziner::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: LienLeschziner.C:234
Foam::incompressible::RASModels::LienLeschziner::Anu_
dimensionedScalar Anu_
Definition: LienLeschziner.H:94
Foam::incompressible::RASModels::LienLeschziner::correctNut
virtual void correctNut()
Definition: LienLeschziner.C:79
Foam::incompressible::RASModels::LienLeschziner::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: LienLeschziner.C:257
Foam::eddyViscosity< incompressible::RASModel >
Foam::eddyViscosity::correct
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
Definition: eddyViscosity.C:132
Foam::eddyViscosity< incompressible::RASModel >::nut_
volScalarField nut_
Definition: eddyViscosity.H:63
Foam::incompressible::RASModels::LienLeschziner::Ceps1_
dimensionedScalar Ceps1_
Definition: LienLeschziner.H:86
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Foam::tmp::clear
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:172
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:93
Foam::incompressible::RASModels::LienLeschziner::y_
const volScalarField & y_
Wall distance.
Definition: LienLeschziner.H:107
Foam::incompressible::RASModels::LienLeschziner::AE_
dimensionedScalar AE_
Definition: LienLeschziner.H:96
Foam::incompressible::RASModels::LienLeschziner::DkEff
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: LienLeschziner.H:150
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::incompressible::RASModels::LienLeschziner::k_
volScalarField k_
Definition: LienLeschziner.H:101