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::boundEpsilon()
60 {
61  tmp<volScalarField> tCmuk2(Cmu_*sqr(k_));
62  epsilon_ = max(epsilon_, tCmuk2()/(nutMaxCoeff_*nu()));
63  return tCmuk2;
64 }
65 
66 
67 tmp<volScalarField> LienLeschziner::fMu() const
68 {
69  const volScalarField yStar(sqrt(k_)*y()/nu());
70 
71  return
72  (scalar(1) - exp(-Anu_*yStar))
73  /((scalar(1) + small) - exp(-Aeps_*yStar));
74 }
75 
76 
77 tmp<volScalarField> LienLeschziner::f2() const
78 {
79  tmp<volScalarField> Rt = sqr(k_)/(nu()*epsilon_);
80 
81  return scalar(1) - 0.3*exp(-sqr(Rt));
82 }
83 
84 
85 tmp<volScalarField> LienLeschziner::E(const volScalarField& f2) const
86 {
87  const volScalarField yStar(sqrt(k_)*y()/nu());
88  const volScalarField le
89  (
90  kappa_*y()*((scalar(1) + small) - exp(-Aeps_*yStar))
91  );
92 
93  return
94  (Ceps2_*pow(Cmu_, 0.75))
95  *(f2*sqrt(k_)*epsilon_/le)*exp(-AE_*sqr(yStar));
96 }
97 
98 
100 {
103 }
104 
105 
106 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
107 
109 (
110  const geometricOneField& alpha,
111  const geometricOneField& rho,
112  const volVectorField& U,
113  const surfaceScalarField& alphaRhoPhi,
114  const surfaceScalarField& phi,
115  const viscosity& viscosity,
116  const word& type
117 )
118 :
119  eddyViscosity<incompressible::RASModel>
120  (
121  type,
122  alpha,
123  rho,
124  U,
125  alphaRhoPhi,
126  phi,
127  viscosity
128  ),
129 
130  Ceps1_
131  (
133  (
134  "Ceps1",
135  coeffDict_,
136  1.44
137  )
138  ),
139  Ceps2_
140  (
142  (
143  "Ceps2",
144  coeffDict_,
145  1.92
146  )
147  ),
148  sigmak_
149  (
151  (
152  "sigmak",
153  coeffDict_,
154  1.0
155  )
156  ),
157  sigmaEps_
158  (
160  (
161  "sigmaEps",
162  coeffDict_,
163  1.3
164  )
165  ),
166  Cmu_
167  (
169  (
170  "Cmu",
171  coeffDict_,
172  0.09
173  )
174  ),
175  kappa_
176  (
178  (
179  "kappa",
180  coeffDict_,
181  0.41
182  )
183  ),
184  Anu_
185  (
187  (
188  "Anu",
189  coeffDict_,
190  0.016
191  )
192  ),
193  Aeps_
194  (
196  (
197  "Aeps",
198  coeffDict_,
199  0.263
200  )
201  ),
202  AE_
203  (
205  (
206  "AE",
207  coeffDict_,
208  0.00222
209  )
210  ),
211 
212  k_
213  (
214  IOobject
215  (
216  this->groupName("k"),
217  runTime_.name(),
218  mesh_,
221  ),
222  mesh_
223  ),
224 
225  epsilon_
226  (
227  IOobject
228  (
229  this->groupName("epsilon"),
230  runTime_.name(),
231  mesh_,
234  ),
235  mesh_
236  )
237 {
238  bound(k_, kMin_);
239  boundEpsilon();
240 
241  if (type == typeName)
242  {
243  printCoeffs(type);
244  }
245 }
246 
247 
248 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
249 
251 {
253  {
254  Ceps1_.readIfPresent(coeffDict());
255  Ceps2_.readIfPresent(coeffDict());
256  sigmak_.readIfPresent(coeffDict());
257  sigmaEps_.readIfPresent(coeffDict());
258  Cmu_.readIfPresent(coeffDict());
259  kappa_.readIfPresent(coeffDict());
260  Anu_.readIfPresent(coeffDict());
261  Aeps_.readIfPresent(coeffDict());
262  AE_.readIfPresent(coeffDict());
263 
264  return true;
265  }
266  else
267  {
268  return false;
269  }
270 }
271 
272 
274 {
275  if (!turbulence_)
276  {
277  return;
278  }
279 
281 
282  tmp<volTensorField> tgradU = fvc::grad(U_);
284  (
285  GName(),
286  nut_*(tgradU() && twoSymm(tgradU()))
287  );
288  tgradU.clear();
289 
290  // Update epsilon and G at the wall
292 
293  const volScalarField f2(this->f2());
294 
295  // Dissipation equation
296  tmp<fvScalarMatrix> epsEqn
297  (
299  + fvm::div(phi_, epsilon_)
301  ==
304  + E(f2)
305  );
306 
307  epsEqn.ref().relax();
308  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
309  solve(epsEqn);
310  boundEpsilon();
311 
312 
313  // Turbulent kinetic energy equation
314  tmp<fvScalarMatrix> kEqn
315  (
316  fvm::ddt(k_)
317  + fvm::div(phi_, k_)
318  - fvm::laplacian(DkEff(), k_)
319  ==
320  G
321  - fvm::Sp(epsilon_/k_, k_)
322  );
323 
324  kEqn.ref().relax();
325  solve(kEqn);
326  bound(k_, kMin_);
327 
328  correctNut();
329 }
330 
331 
332 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
333 
334 } // End namespace RASModels
335 } // End namespace incompressible
336 } // End namespace Foam
337 
338 // ************************************************************************* //
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.
virtual void correctNut()
Correct the eddy-viscosity nut.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
tmp< volScalarField > boundEpsilon()
Bound epsilon and return Cmu*sqr(k) for nut.
tmp< volScalarField > E(const volScalarField &f2) const
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 > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
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:65
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
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:64
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
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.
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))