LienLeschziner.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-2023 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 
32 (
33  geometricOneField,
34  geometricOneField,
35  incompressibleMomentumTransportModel
36 )
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 namespace incompressible
43 {
44 namespace RASModels
45 {
46 
47 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
48 
49 defineTypeNameAndDebug(LienLeschziner, 0);
51 (
52  RASincompressibleMomentumTransportModel,
53  LienLeschziner,
54  dictionary
55 );
56 
57 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
58 
59 tmp<volScalarField> LienLeschziner::fMu() const
60 {
61  const volScalarField yStar(sqrt(k_)*y_/nu());
62 
63  return
64  (scalar(1) - exp(-Anu_*yStar))
65  /((scalar(1) + small) - exp(-Aeps_*yStar));
66 }
67 
68 
69 tmp<volScalarField> LienLeschziner::f2() const
70 {
71  tmp<volScalarField> Rt = sqr(k_)/(nu()*epsilon_);
72 
73  return scalar(1) - 0.3*exp(-sqr(Rt));
74 }
75 
76 
77 tmp<volScalarField> LienLeschziner::E(const volScalarField& f2) const
78 {
79  const volScalarField yStar(sqrt(k_)*y_/nu());
80  const volScalarField le
81  (
82  kappa_*y_*((scalar(1) + small) - exp(-Aeps_*yStar))
83  );
84 
85  return
86  (Ceps2_*pow(Cmu_, 0.75))
87  *(f2*sqrt(k_)*epsilon_/le)*exp(-AE_*sqr(yStar));
88 }
89 
90 
92 {
93  nut_ = Cmu_*fMu()*sqr(k_)/epsilon_;
95 }
96 
97 
98 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
99 
101 (
102  const geometricOneField& alpha,
103  const geometricOneField& rho,
104  const volVectorField& U,
105  const surfaceScalarField& alphaRhoPhi,
106  const surfaceScalarField& phi,
107  const viscosity& viscosity,
108  const word& type
109 )
110 :
111  eddyViscosity<incompressible::RASModel>
112  (
113  type,
114  alpha,
115  rho,
116  U,
117  alphaRhoPhi,
118  phi,
119  viscosity
120  ),
121 
122  Ceps1_
123  (
125  (
126  "Ceps1",
127  coeffDict_,
128  1.44
129  )
130  ),
131  Ceps2_
132  (
134  (
135  "Ceps2",
136  coeffDict_,
137  1.92
138  )
139  ),
140  sigmak_
141  (
143  (
144  "sigmak",
145  coeffDict_,
146  1.0
147  )
148  ),
149  sigmaEps_
150  (
152  (
153  "sigmaEps",
154  coeffDict_,
155  1.3
156  )
157  ),
158  Cmu_
159  (
161  (
162  "Cmu",
163  coeffDict_,
164  0.09
165  )
166  ),
167  kappa_
168  (
170  (
171  "kappa",
172  coeffDict_,
173  0.41
174  )
175  ),
176  Anu_
177  (
179  (
180  "Anu",
181  coeffDict_,
182  0.016
183  )
184  ),
185  Aeps_
186  (
188  (
189  "Aeps",
190  coeffDict_,
191  0.263
192  )
193  ),
194  AE_
195  (
197  (
198  "AE",
199  coeffDict_,
200  0.00222
201  )
202  ),
203 
204  k_
205  (
206  IOobject
207  (
208  this->groupName("k"),
209  runTime_.name(),
210  mesh_,
213  ),
214  mesh_
215  ),
216 
217  epsilon_
218  (
219  IOobject
220  (
221  this->groupName("epsilon"),
222  runTime_.name(),
223  mesh_,
226  ),
227  mesh_
228  ),
229 
230  y_(wallDist::New(mesh_).y())
231 {
232  bound(k_, kMin_);
233  bound(epsilon_, epsilonMin_);
234 
235  if (type == typeName)
236  {
237  printCoeffs(type);
238  }
239 }
240 
241 
242 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
243 
245 {
247  {
248  Ceps1_.readIfPresent(coeffDict());
249  Ceps2_.readIfPresent(coeffDict());
250  sigmak_.readIfPresent(coeffDict());
251  sigmaEps_.readIfPresent(coeffDict());
252  Cmu_.readIfPresent(coeffDict());
253  kappa_.readIfPresent(coeffDict());
254  Anu_.readIfPresent(coeffDict());
255  Aeps_.readIfPresent(coeffDict());
256  AE_.readIfPresent(coeffDict());
257 
258  return true;
259  }
260  else
261  {
262  return false;
263  }
264 }
265 
266 
268 {
269  if (!turbulence_)
270  {
271  return;
272  }
273 
275 
276  tmp<volTensorField> tgradU = fvc::grad(U_);
278  (
279  GName(),
280  nut_*(tgradU() && twoSymm(tgradU()))
281  );
282  tgradU.clear();
283 
284  // Update epsilon and G at the wall
286 
287  const volScalarField f2(this->f2());
288 
289  // Dissipation equation
290  tmp<fvScalarMatrix> epsEqn
291  (
293  + fvm::div(phi_, epsilon_)
295  ==
298  + E(f2)
299  );
300 
301  epsEqn.ref().relax();
302  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
303  solve(epsEqn);
304  bound(epsilon_, epsilonMin_);
305 
306 
307  // Turbulent kinetic energy equation
308  tmp<fvScalarMatrix> kEqn
309  (
310  fvm::ddt(k_)
311  + fvm::div(phi_, k_)
312  - fvm::laplacian(DkEff(), k_)
313  ==
314  G
315  - fvm::Sp(epsilon_/k_, k_)
316  );
317 
318  kEqn.ref().relax();
319  solve(kEqn);
320  bound(k_, kMin_);
321 
322  correctNut();
323 }
324 
325 
326 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
327 
328 } // End namespace RASModels
329 } // End namespace incompressible
330 } // End namespace Foam
331 
332 // ************************************************************************* //
scalar y
makeMomentumTransportModelTypes(geometricOneField, geometricOneField, incompressibleMomentumTransportModel)
Bound the given scalar field where it is below the specified minimum.
void updateCoeffs()
Update the boundary condition coefficients.
void correctBoundaryConditions()
Correct boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
static dimensioned< scalar > lookupOrAddToDict(const word &, dictionary &, const dimensionSet &dims=dimless, const scalar &defaultValue=pTraits< scalar >::zero)
Construct from dictionary, with default value.
bool readIfPresent(const dictionary &)
Update the value of dimensioned<Type> if found in the dictionary.
virtual bool read()=0
Re-read model coefficients if they have changed.
Definition: eddyViscosity.C:74
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
tmp< volScalarField > E(const volScalarField &f2) const
const volScalarField & y_
Wall distance.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
LienLeschziner(const geometricOneField &alpha, const geometricOneField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
virtual bool read()
Read RASProperties dictionary.
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
VolField< vector > volVectorField
Definition: volFieldsFwd.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
SurfaceField< scalar > surfaceScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
defineTypeNameAndDebug(combustionModel, 0)
bool bound(volScalarField &, const dimensionedScalar &min)
Bound the given scalar field where it is below the specified min value.
Definition: bound.C:31
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.name(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Creating face flux\n"<< endl;surfaceScalarField phi(IOobject("phi", runTime.name(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), mesh, dimensionedScalar(mesh.Sf().dimensions() *U.dimensions(), 0));autoPtr< viscosityModel > viscosity(viscosityModel::New(mesh))