LamBremhorstKE.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 "LamBremhorstKE.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(LamBremhorstKE, 0);
51 (
52  RASincompressibleMomentumTransportModel,
53  LamBremhorstKE,
54  dictionary
55 );
56 
57 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
58 
59 tmp<volScalarField> LamBremhorstKE::Rt() const
60 {
61  return sqr(k_)/(nu()*epsilon_);
62 }
63 
64 
65 tmp<volScalarField> LamBremhorstKE::fMu(const volScalarField& Rt) const
66 {
67  tmp<volScalarField> Ry(sqrt(k_)*y_/nu());
68  return sqr(scalar(1) - exp(-0.0165*Ry))*(scalar(1) + 20.5/(Rt + small));
69 }
70 
71 
72 tmp<volScalarField> LamBremhorstKE::f1(const volScalarField& fMu) const
73 {
74  return scalar(1) + pow3(0.05/(fMu + small));
75 }
76 
77 
78 tmp<volScalarField> LamBremhorstKE::f2(const volScalarField& Rt) const
79 {
80  return scalar(1) - exp(-sqr(Rt));
81 }
82 
83 
85 {
86  nut_ = Cmu_*fMu*sqr(k_)/epsilon_;
88 }
89 
90 
91 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
92 
94 {
95  correctNut(fMu(Rt()));
96 }
97 
98 
99 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
100 
102 (
103  const geometricOneField& alpha,
104  const geometricOneField& rho,
105  const volVectorField& U,
106  const surfaceScalarField& alphaRhoPhi,
107  const surfaceScalarField& phi,
108  const viscosity& viscosity,
109  const word& type
110 )
111 :
112  eddyViscosity<incompressible::RASModel>
113  (
114  type,
115  alpha,
116  rho,
117  U,
118  alphaRhoPhi,
119  phi,
120  viscosity
121  ),
122 
123  Cmu_
124  (
126  (
127  "Cmu",
128  coeffDict_,
129  0.09
130  )
131  ),
132  Ceps1_
133  (
135  (
136  "Ceps1",
137  coeffDict_,
138  1.44
139  )
140  ),
141  Ceps2_
142  (
144  (
145  "Ceps2",
146  coeffDict_,
147  1.92
148  )
149  ),
150  sigmaEps_
151  (
153  (
154  "alphaEps",
155  coeffDict_,
156  1.3
157  )
158  ),
159 
160  k_
161  (
162  IOobject
163  (
164  this->groupName("k"),
165  runTime_.name(),
166  mesh_,
169  ),
170  mesh_
171  ),
172 
173  epsilon_
174  (
175  IOobject
176  (
177  this->groupName("epsilon"),
178  runTime_.name(),
179  mesh_,
182  ),
183  mesh_
184  ),
185 
186  y_(wallDist::New(mesh_).y())
187 {
188  bound(k_, kMin_);
189  bound(epsilon_, epsilonMin_);
190 
191  if (type == typeName)
192  {
193  printCoeffs(type);
194  }
195 }
196 
197 
198 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
199 
201 {
203  {
204  Cmu_.readIfPresent(coeffDict());
205  Ceps1_.readIfPresent(coeffDict());
206  Ceps2_.readIfPresent(coeffDict());
207  sigmaEps_.readIfPresent(coeffDict());
208 
209  return true;
210  }
211  else
212  {
213  return false;
214  }
215 }
216 
217 
219 {
220  if (!turbulence_)
221  {
222  return;
223  }
224 
226 
227  tmp<volTensorField> tgradU = fvc::grad(U_);
228  volScalarField G(GName(), nut_*(twoSymm(tgradU()) && tgradU()));
229  tgradU.clear();
230 
231  // Update epsilon and G at the wall
233 
234  const volScalarField Rt(this->Rt());
235  const volScalarField fMu(this->fMu(Rt));
236 
237  // Dissipation equation
238  tmp<fvScalarMatrix> epsEqn
239  (
241  + fvm::div(phi_, epsilon_)
243  ==
244  Ceps1_*f1(fMu)*G*epsilon_/k_
245  - fvm::Sp(Ceps2_*f2(Rt)*epsilon_/k_, epsilon_)
246  );
247 
248  epsEqn.ref().relax();
249  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
250  solve(epsEqn);
251  bound(epsilon_, epsilonMin_);
252 
253  // Turbulent kinetic energy equation
254  tmp<fvScalarMatrix> kEqn
255  (
256  fvm::ddt(k_)
257  + fvm::div(phi_, k_)
258  - fvm::laplacian(DkEff(), k_)
259  ==
260  G - fvm::Sp(epsilon_/k_, k_)
261  );
262 
263  kEqn.ref().relax();
264  solve(kEqn);
265  bound(k_, kMin_);
266 
267  correctNut(fMu);
268 }
269 
270 
271 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
272 
273 } // End namespace RASModels
274 } // End namespace incompressible
275 } // End namespace Foam
276 
277 // ************************************************************************* //
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.
LamBremhorstKE(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.
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.
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)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
SurfaceField< scalar > surfaceScalarField
tensor Ry(const scalar &omega)
Rotational transformation tensor about the y-axis by omega radians.
Definition: transform.H:106
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))