LamBremhorstKE.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2020,2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "LamBremhorstKE.H"
30 #include "wallDist.H"
31 #include "bound.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace incompressible
39 {
40 namespace RASModels
41 {
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 defineTypeNameAndDebug(LamBremhorstKE, 0);
46 addToRunTimeSelectionTable(RASModel, LamBremhorstKE, dictionary);
47 
48 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
49 
50 tmp<volScalarField> LamBremhorstKE::Rt() const
51 {
52  return sqr(k_)/(nu()*epsilon_);
53 }
54 
55 
56 tmp<volScalarField> LamBremhorstKE::fMu(const volScalarField& Rt) const
57 {
58  tmp<volScalarField> Ry(sqrt(k_)*y_/nu());
59  return sqr(scalar(1) - exp(-0.0165*Ry))*(scalar(1) + 20.5/(Rt + SMALL));
60 }
61 
62 
63 tmp<volScalarField> LamBremhorstKE::f1(const volScalarField& fMu) const
64 {
65  return scalar(1) + pow3(0.05/(fMu + SMALL));
66 }
67 
68 
69 tmp<volScalarField> LamBremhorstKE::f2(const volScalarField& Rt) const
70 {
71  return scalar(1) - exp(-sqr(Rt));
72 }
73 
74 
76 {
79 }
80 
81 
82 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
83 
85 {
86  correctNut(fMu(Rt()));
87 }
88 
89 
90 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
91 
92 LamBremhorstKE::LamBremhorstKE
93 (
94  const geometricOneField& alpha,
95  const geometricOneField& rho,
96  const volVectorField& U,
97  const surfaceScalarField& alphaRhoPhi,
98  const surfaceScalarField& phi,
99  const transportModel& transport,
100  const word& propertiesName,
101  const word& type
102 )
103 :
104  eddyViscosity<incompressible::RASModel>
105  (
106  type,
107  alpha,
108  rho,
109  U,
110  alphaRhoPhi,
111  phi,
112  transport,
113  propertiesName
114  ),
115 
116  Cmu_
117  (
118  dimensioned<scalar>::getOrAddToDict
119  (
120  "Cmu",
121  coeffDict_,
122  0.09
123  )
124  ),
125  Ceps1_
126  (
127  dimensioned<scalar>::getOrAddToDict
128  (
129  "Ceps1",
130  coeffDict_,
131  1.44
132  )
133  ),
134  Ceps2_
135  (
136  dimensioned<scalar>::getOrAddToDict
137  (
138  "Ceps2",
139  coeffDict_,
140  1.92
141  )
142  ),
143  sigmaEps_
144  (
145  dimensioned<scalar>::getOrAddToDict
146  (
147  "alphaEps",
148  coeffDict_,
149  1.3
150  )
151  ),
152 
153  k_
154  (
155  IOobject
156  (
157  IOobject::groupName("k", alphaRhoPhi.group()),
158  runTime_.timeName(),
159  mesh_,
160  IOobject::MUST_READ,
161  IOobject::AUTO_WRITE
162  ),
163  mesh_
164  ),
165 
166  epsilon_
167  (
168  IOobject
169  (
170  IOobject::groupName("epsilon", alphaRhoPhi.group()),
171  runTime_.timeName(),
172  mesh_,
173  IOobject::MUST_READ,
174  IOobject::AUTO_WRITE
175  ),
176  mesh_
177  ),
178 
179  y_(wallDist::New(mesh_).y())
180 {
181  bound(k_, kMin_);
182  bound(epsilon_, epsilonMin_);
183 
184  if (type == typeName)
185  {
186  printCoeffs(type);
187  }
188 }
189 
190 
191 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
192 
194 {
196  {
197  Cmu_.readIfPresent(coeffDict());
198  Ceps1_.readIfPresent(coeffDict());
199  Ceps2_.readIfPresent(coeffDict());
200  sigmaEps_.readIfPresent(coeffDict());
201 
202  return true;
203  }
204 
205  return false;
206 }
207 
208 
210 {
211  if (!turbulence_)
212  {
213  return;
214  }
215 
217 
218  tmp<volTensorField> tgradU = fvc::grad(U_);
219  volScalarField G(GName(), nut_*(twoSymm(tgradU()) && tgradU()));
220  tgradU.clear();
221 
222  // Update epsilon and G at the wall
224  // Push any changed cell values to coupled neighbours
226 
227  const volScalarField Rt(this->Rt());
228  const volScalarField fMu(this->fMu(Rt));
229 
230  // Dissipation equation
231  tmp<fvScalarMatrix> epsEqn
232  (
234  + fvm::div(phi_, epsilon_)
236  ==
237  Ceps1_*f1(fMu)*G*epsilon_/k_
238  - fvm::Sp(Ceps2_*f2(Rt)*epsilon_/k_, epsilon_)
239  );
240 
241  epsEqn.ref().relax();
242  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
243  solve(epsEqn);
244  bound(epsilon_, epsilonMin_);
245 
246  // Turbulent kinetic energy equation
248  (
249  fvm::ddt(k_)
250  + fvm::div(phi_, k_)
251  - fvm::laplacian(DkEff(), k_)
252  ==
253  G - fvm::Sp(epsilon_/k_, k_)
254  );
255 
256  kEqn.ref().relax();
257  solve(kEqn);
258  bound(k_, kMin_);
259 
260  // Update nut with latest available k,epsilon
261  correctNut();
262 }
263 
264 
265 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
266 
267 } // End namespace RASModels
268 } // End namespace incompressible
269 } // End namespace Foam
270 
271 // ************************************************************************* //
virtual bool read()
Re-read model coefficients if they have changed.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
const dimensionedScalar G
Newtonian constant of gravitation.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
bool readIfPresent(const dictionary &dict)
Update the value of dimensioned<Type> if found in the dictionary, lookup in dictionary with the name(...
An abstract base class for patches that couple regions of the computational domain e...
RASModel< turbulenceModel > RASModel
dimensionedScalar sqrt(const dimensionedScalar &ds)
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:51
Macros for easy insertion into run-time selection tables.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
constexpr const char *const group
Group name for atomic constants.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
word timeName
Definition: getTimeIndex.H:3
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
scalar y
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
dimensionedScalar exp(const dimensionedScalar &ds)
A class for handling words, derived from Foam::string.
Definition: word.H:63
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
void evaluateCoupled()
Evaluate boundary conditions on a subset of coupled patches.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
void updateCoeffs()
Update the boundary condition coefficients.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
Bound the given scalar field if it has gone unbounded.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
const volScalarField & y_
Wall distance.
tensor Ry(const scalar omega)
Rotational transformation tensor about the y-axis by omega radians.
Definition: transform.H:99
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
dimensionedScalar pow3(const dimensionedScalar &ds)
U
Definition: pEqn.H:72
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:29
addToRunTimeSelectionTable(RASModel, kkLOmega, dictionary)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
Base-class for all transport models used by the incompressible turbulence models. ...
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition: tmpI.H:289
Interface to run-time selectable methods to calculate the distance-to-wall and normal-to-wall fields...
Definition: wallDist.H:71
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
volScalarField & nu
Namespace for OpenFOAM.