LienCubicKE.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 "LienCubicKE.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(LienCubicKE, 0);
51 (
52  RASincompressibleMomentumTransportModel,
53  LienCubicKE,
54  dictionary
55 );
56 
57 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
58 
59 tmp<volScalarField> LienCubicKE::fMu() const
60 {
61  const volScalarField yStar(sqrt(k_)*y_/nu());
62 
63  return
64  (scalar(1) - exp(-Anu_*yStar))
65  *(scalar(1) + (2*kappa_/(pow(Cmu_, 0.75))/(yStar + small)));
66 }
67 
68 
69 tmp<volScalarField> LienCubicKE::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> LienCubicKE::E(const volScalarField& f2) const
78 {
79  const volScalarField yStar(sqrt(k_)*y_/nu());
80  const volScalarField le
81  (
82  kappa_*y_/(scalar(1) + (2*kappa_/(pow(Cmu_, 0.75))/(yStar + small)))
83  );
84 
85  return
86  (Ceps2_*pow(Cmu_, 0.75))
87  *(f2*sqrt(k_)*epsilon_/le)*exp(-AE_*sqr(yStar));
88 }
89 
90 
92 {
94 }
95 
96 
98 {
99  volSymmTensorField S(symm(gradU));
100  volTensorField W(skew(gradU));
101 
102  volScalarField sBar((k_/epsilon_)*sqrt(2.0)*mag(S));
103  volScalarField wBar((k_/epsilon_)*sqrt(2.0)*mag(W));
104 
105  volScalarField Cmu((2.0/3.0)/(Cmu1_ + sBar + Cmu2_*wBar));
106  volScalarField fMu(this->fMu());
107 
108  nut_ = Cmu*fMu*sqr(k_)/epsilon_;
110 
112  fMu*k_
113  *(
114  // Quadratic terms
115  sqr(k_/epsilon_)/(Cbeta_ + pow3(sBar))
116  *(
118  + Cbeta2_*twoSymm(S&W)
119  + Cbeta3_*dev(symm(W&W))
120  )
121 
122  // Cubic terms
123  - pow3(Cmu*k_/epsilon_)
124  *(
126  + Cgamma4_*twoSymm((innerSqr(S)&W))
127  )
128  );
129 }
130 
131 
132 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
133 
135 (
136  const geometricOneField& alpha,
137  const geometricOneField& rho,
138  const volVectorField& U,
139  const surfaceScalarField& alphaRhoPhi,
140  const surfaceScalarField& phi,
141  const viscosity& viscosity,
142  const word& type
143 )
144 :
145  nonlinearEddyViscosity<incompressible::RASModel>
146  (
147  type,
148  alpha,
149  rho,
150  U,
151  alphaRhoPhi,
152  phi,
153  viscosity
154  ),
155 
156  Ceps1_
157  (
159  (
160  "Ceps1",
161  coeffDict_,
162  1.44
163  )
164  ),
165  Ceps2_
166  (
168  (
169  "Ceps2",
170  coeffDict_,
171  1.92
172  )
173  ),
174  sigmak_
175  (
177  (
178  "sigmak",
179  coeffDict_,
180  1.0
181  )
182  ),
183  sigmaEps_
184  (
186  (
187  "sigmaEps",
188  coeffDict_,
189  1.3
190  )
191  ),
192  Cmu1_
193  (
195  (
196  "Cmu1",
197  coeffDict_,
198  1.25
199  )
200  ),
201  Cmu2_
202  (
204  (
205  "Cmu2",
206  coeffDict_,
207  0.9
208  )
209  ),
210  Cbeta_
211  (
213  (
214  "Cbeta",
215  coeffDict_,
216  1000.0
217  )
218  ),
219  Cbeta1_
220  (
222  (
223  "Cbeta1",
224  coeffDict_,
225  3.0
226  )
227  ),
228  Cbeta2_
229  (
231  (
232  "Cbeta2",
233  coeffDict_,
234  15.0
235  )
236  ),
237  Cbeta3_
238  (
240  (
241  "Cbeta3",
242  coeffDict_,
243  -19.0
244  )
245  ),
246  Cgamma1_
247  (
249  (
250  "Cgamma1",
251  coeffDict_,
252  16.0
253  )
254  ),
255  Cgamma2_
256  (
258  (
259  "Cgamma2",
260  coeffDict_,
261  16.0
262  )
263  ),
264  Cgamma4_
265  (
267  (
268  "Cgamma4",
269  coeffDict_,
270  -80.0
271  )
272  ),
273  Cmu_
274  (
276  (
277  "Cmu",
278  coeffDict_,
279  0.09
280  )
281  ),
282  kappa_
283  (
285  (
286  "kappa",
287  coeffDict_,
288  0.41
289  )
290  ),
291  Anu_
292  (
294  (
295  "Anu",
296  coeffDict_,
297  0.0198
298  )
299  ),
300  AE_
301  (
303  (
304  "AE",
305  coeffDict_,
306  0.00375
307  )
308  ),
309 
310  k_
311  (
312  IOobject
313  (
314  this->groupName("k"),
315  runTime_.name(),
316  mesh_,
319  ),
320  mesh_
321  ),
322 
323  epsilon_
324  (
325  IOobject
326  (
327  this->groupName("epsilon"),
328  runTime_.name(),
329  mesh_,
332  ),
333  mesh_
334  ),
335 
336  y_(wallDist::New(mesh_).y())
337 {
338  bound(k_, kMin_);
339  bound(epsilon_, epsilonMin_);
340 
341  if (type == typeName)
342  {
343  printCoeffs(type);
344  }
345 }
346 
347 
348 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
349 
350 bool LienCubicKE::read()
351 {
353  {
354  Ceps1_.readIfPresent(coeffDict());
355  Ceps2_.readIfPresent(coeffDict());
356  sigmak_.readIfPresent(coeffDict());
357  sigmaEps_.readIfPresent(coeffDict());
358  Cmu1_.readIfPresent(coeffDict());
359  Cmu2_.readIfPresent(coeffDict());
360  Cbeta_.readIfPresent(coeffDict());
361  Cbeta1_.readIfPresent(coeffDict());
362  Cbeta2_.readIfPresent(coeffDict());
363  Cbeta3_.readIfPresent(coeffDict());
364  Cgamma1_.readIfPresent(coeffDict());
365  Cgamma2_.readIfPresent(coeffDict());
366  Cgamma4_.readIfPresent(coeffDict());
367  Cmu_.readIfPresent(coeffDict());
368  kappa_.readIfPresent(coeffDict());
369  Anu_.readIfPresent(coeffDict());
370  AE_.readIfPresent(coeffDict());
371 
372  return true;
373  }
374  else
375  {
376  return false;
377  }
378 }
379 
380 
382 {
383  if (!turbulence_)
384  {
385  return;
386  }
387 
389 
390  tmp<volTensorField> tgradU = fvc::grad(U_);
391  const volTensorField& gradU = tgradU();
392 
394  (
395  GName(),
396  (nut_*twoSymm(gradU) - nonlinearStress_) && gradU
397  );
398 
399 
400  // Update epsilon and G at the wall
402 
403  const volScalarField f2(this->f2());
404 
405  // Dissipation equation
406  tmp<fvScalarMatrix> epsEqn
407  (
409  + fvm::div(phi_, epsilon_)
411  ==
414  + E(f2)
415  );
416 
417  epsEqn.ref().relax();
418  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
419  solve(epsEqn);
420  bound(epsilon_, epsilonMin_);
421 
422 
423  // Turbulent kinetic energy equation
424  tmp<fvScalarMatrix> kEqn
425  (
426  fvm::ddt(k_)
427  + fvm::div(phi_, k_)
428  - fvm::laplacian(DkEff(), k_)
429  ==
430  G
431  - fvm::Sp(epsilon_/k_, k_)
432  );
433 
434  kEqn.ref().relax();
435  solve(kEqn);
436  bound(k_, kMin_);
437 
438 
439  // Re-calculate viscosity and non-linear stress
440  correctNonlinearStress(gradU);
441 }
442 
443 
444 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
445 
446 } // End namespace RASModels
447 } // End namespace incompressible
448 } // End namespace Foam
449 
450 // ************************************************************************* //
scalar y
makeMomentumTransportModelTypes(geometricOneField, geometricOneField, incompressibleMomentumTransportModel)
Definition: LienCubicKE.C:32
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.
LienCubicKE(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.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: LienCubicKE.H:158
tmp< volScalarField > f2() const
tmp< volScalarField > E(const volScalarField &f2) const
const volScalarField & y_
Wall distance.
Definition: LienCubicKE.H:114
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: LienCubicKE.H:168
virtual void correctNonlinearStress(const volTensorField &gradU)
virtual bool read()
Read RASProperties dictionary.
tmp< volScalarField > fMu() const
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 > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
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.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensionedScalar exp(const dimensionedScalar &ds)
VolField< vector > volVectorField
Definition: volFieldsFwd.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
SurfaceField< scalar > surfaceScalarField
VolField< tensor > volTensorField
Definition: volFieldsFwd.H:65
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
scalarList W(const fluidMulticomponentThermo &thermo)
dimensionedScalar sqrt(const dimensionedScalar &ds)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
dimensioned< scalar > mag(const dimensioned< Type > &)
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
dimensionedSymmTensor innerSqr(const dimensionedSymmTensor &dt)
VolField< symmTensor > volSymmTensorField
Definition: volFieldsFwd.H:64
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dimensionedTensor skew(const dimensionedTensor &dt)
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))