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-2021 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 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace incompressible
36 {
37 namespace RASModels
38 {
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 defineTypeNameAndDebug(LamBremhorstKE, 0);
43 addToRunTimeSelectionTable(RASModel, LamBremhorstKE, dictionary);
44 
45 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
46 
47 tmp<volScalarField> LamBremhorstKE::Rt() const
48 {
49  return sqr(k_)/(nu()*epsilon_);
50 }
51 
52 
53 tmp<volScalarField> LamBremhorstKE::fMu(const volScalarField& Rt) const
54 {
55  tmp<volScalarField> Ry(sqrt(k_)*y_/nu());
56  return sqr(scalar(1) - exp(-0.0165*Ry))*(scalar(1) + 20.5/(Rt + small));
57 }
58 
59 
60 tmp<volScalarField> LamBremhorstKE::f1(const volScalarField& fMu) const
61 {
62  return scalar(1) + pow3(0.05/(fMu + small));
63 }
64 
65 
66 tmp<volScalarField> LamBremhorstKE::f2(const volScalarField& Rt) const
67 {
68  return scalar(1) - exp(-sqr(Rt));
69 }
70 
71 
73 {
74  nut_ = Cmu_*fMu*sqr(k_)/epsilon_;
76 }
77 
78 
79 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
80 
82 {
83  correctNut(fMu(Rt()));
84 }
85 
86 
87 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
88 
90 (
91  const geometricOneField& alpha,
92  const geometricOneField& rho,
93  const volVectorField& U,
94  const surfaceScalarField& alphaRhoPhi,
95  const surfaceScalarField& phi,
96  const viscosity& viscosity,
97  const word& type
98 )
99 :
101  (
102  type,
103  alpha,
104  rho,
105  U,
106  alphaRhoPhi,
107  phi,
108  viscosity
109  ),
110 
111  Cmu_
112  (
114  (
115  "Cmu",
116  coeffDict_,
117  0.09
118  )
119  ),
120  Ceps1_
121  (
123  (
124  "Ceps1",
125  coeffDict_,
126  1.44
127  )
128  ),
129  Ceps2_
130  (
132  (
133  "Ceps2",
134  coeffDict_,
135  1.92
136  )
137  ),
138  sigmaEps_
139  (
141  (
142  "alphaEps",
143  coeffDict_,
144  1.3
145  )
146  ),
147 
148  k_
149  (
150  IOobject
151  (
152  IOobject::groupName("k", alphaRhoPhi.group()),
153  runTime_.timeName(),
154  mesh_,
157  ),
158  mesh_
159  ),
160 
161  epsilon_
162  (
163  IOobject
164  (
165  IOobject::groupName("epsilon", alphaRhoPhi.group()),
166  runTime_.timeName(),
167  mesh_,
170  ),
171  mesh_
172  ),
173 
174  y_(wallDist::New(mesh_).y())
175 {
176  bound(k_, kMin_);
177  bound(epsilon_, epsilonMin_);
178 
179  if (type == typeName)
180  {
181  printCoeffs(type);
182  }
183 }
184 
185 
186 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
187 
189 {
191  {
192  Cmu_.readIfPresent(coeffDict());
193  Ceps1_.readIfPresent(coeffDict());
194  Ceps2_.readIfPresent(coeffDict());
195  sigmaEps_.readIfPresent(coeffDict());
196 
197  return true;
198  }
199  else
200  {
201  return false;
202  }
203 }
204 
205 
207 {
208  if (!turbulence_)
209  {
210  return;
211  }
212 
214 
215  tmp<volTensorField> tgradU = fvc::grad(U_);
216  volScalarField G(GName(), nut_*(twoSymm(tgradU()) && tgradU()));
217  tgradU.clear();
218 
219  // Update epsilon and G at the wall
221 
222  const volScalarField Rt(this->Rt());
223  const volScalarField fMu(this->fMu(Rt));
224 
225  // Dissipation equation
226  tmp<fvScalarMatrix> epsEqn
227  (
229  + fvm::div(phi_, epsilon_)
231  ==
232  Ceps1_*f1(fMu)*G*epsilon_/k_
233  - fvm::Sp(Ceps2_*f2(Rt)*epsilon_/k_, epsilon_)
234  );
235 
236  epsEqn.ref().relax();
237  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
238  solve(epsEqn);
239  bound(epsilon_, epsilonMin_);
240 
241  // Turbulent kinetic energy equation
243  (
244  fvm::ddt(k_)
245  + fvm::div(phi_, k_)
246  - fvm::laplacian(DkEff(), k_)
247  ==
248  G - fvm::Sp(epsilon_/k_, k_)
249  );
250 
251  kEqn.ref().relax();
252  solve(kEqn);
253  bound(k_, kMin_);
254 
255  correctNut(fMu);
256 }
257 
258 
259 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
260 
261 } // End namespace RASModels
262 } // End namespace incompressible
263 } // End namespace Foam
264 
265 // ************************************************************************* //
virtual bool read()
Read RASProperties dictionary.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:237
U
Definition: pEqn.H:72
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
dimensionedScalar sqrt(const dimensionedScalar &ds)
bool readIfPresent(const dictionary &)
Update the value of dimensioned<Type> if found in the dictionary.
Generic dimensioned Type class.
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Macros for easy insertion into run-time selection tables.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
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 string.
Definition: word.H:59
static word groupName(Name name, const word &group)
Info<< "Reading field p\"<< endl;volScalarField p(IOobject("p", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);pressureReference pressureReference(p, simple.dict());mesh.schemes().setFluxRequired(p.name());Info<< "Reading field pa\"<< endl;volScalarField pa(IOobject("pa", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field Ua\"<< endl;volVectorField Ua(IOobject("Ua", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);# 65 "/home/ubuntu/OpenFOAM-10/applications/solvers/incompressible/adjointShapeOptimisationFoam/createFields.H" 2label paRefCell=0;scalar paRefValue=0.0;setRefCell(pa, simple.dict(), paRefCell, paRefValue);mesh.schemes().setFluxRequired(pa.name());autoPtr< viscosityModel > viscosity(viscosityModel::New(mesh))
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
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.
phi
Definition: correctPhi.H:3
RASModel< momentumTransportModel > RASModel
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:46
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:106
dimensionedScalar pow3(const dimensionedScalar &ds)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
addToRunTimeSelectionTable(RASModel, kkLOmega, dictionary)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
void correctBoundaryConditions()
Correct boundary field.
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
Namespace for OpenFOAM.