kEpsilon.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 | Copyright (C) 2011-2015 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 "kEpsilon.H"
27 #include "fvOptions.H"
28 #include "bound.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace RASModels
35 {
36 
37 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38 
39 template<class BasicTurbulenceModel>
41 {
42  this->nut_ = Cmu_*sqr(k_)/epsilon_;
43  this->nut_.correctBoundaryConditions();
44 
45  fv::options& fvOptions(fv::options::New(this->mesh_));
46  fvOptions.correct(this->nut_);
47 
48  BasicTurbulenceModel::correctNut();
49 }
50 
51 
52 template<class BasicTurbulenceModel>
54 {
55  return tmp<fvScalarMatrix>
56  (
57  new fvScalarMatrix
58  (
59  k_,
60  dimVolume*this->rho_.dimensions()*k_.dimensions()
61  /dimTime
62  )
63  );
64 }
65 
66 
67 template<class BasicTurbulenceModel>
69 {
70  return tmp<fvScalarMatrix>
71  (
72  new fvScalarMatrix
73  (
74  epsilon_,
75  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
76  /dimTime
77  )
78  );
79 }
80 
81 
82 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
83 
84 template<class BasicTurbulenceModel>
86 (
87  const alphaField& alpha,
88  const rhoField& rho,
89  const volVectorField& U,
90  const surfaceScalarField& alphaRhoPhi,
91  const surfaceScalarField& phi,
92  const transportModel& transport,
93  const word& propertiesName,
94  const word& type
95 )
96 :
98  (
99  type,
100  alpha,
101  rho,
102  U,
103  alphaRhoPhi,
104  phi,
105  transport,
106  propertiesName
107  ),
108 
109  Cmu_
110  (
112  (
113  "Cmu",
114  this->coeffDict_,
115  0.09
116  )
117  ),
118  C1_
119  (
121  (
122  "C1",
123  this->coeffDict_,
124  1.44
125  )
126  ),
127  C2_
128  (
130  (
131  "C2",
132  this->coeffDict_,
133  1.92
134  )
135  ),
136  C3_
137  (
139  (
140  "C3",
141  this->coeffDict_,
142  -0.33
143  )
144  ),
145  sigmak_
146  (
148  (
149  "sigmak",
150  this->coeffDict_,
151  1.0
152  )
153  ),
154  sigmaEps_
155  (
157  (
158  "sigmaEps",
159  this->coeffDict_,
160  1.3
161  )
162  ),
163 
164  k_
165  (
166  IOobject
167  (
168  IOobject::groupName("k", U.group()),
169  this->runTime_.timeName(),
170  this->mesh_,
173  ),
174  this->mesh_
175  ),
176  epsilon_
177  (
178  IOobject
179  (
180  IOobject::groupName("epsilon", U.group()),
181  this->runTime_.timeName(),
182  this->mesh_,
185  ),
186  this->mesh_
187  )
188 {
189  bound(k_, this->kMin_);
190  bound(epsilon_, this->epsilonMin_);
191 
192  if (type == typeName)
193  {
194  this->printCoeffs(type);
195  }
196 }
197 
198 
199 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
200 
201 template<class BasicTurbulenceModel>
203 {
205  {
206  Cmu_.readIfPresent(this->coeffDict());
207  C1_.readIfPresent(this->coeffDict());
208  C2_.readIfPresent(this->coeffDict());
209  C3_.readIfPresent(this->coeffDict());
210  sigmak_.readIfPresent(this->coeffDict());
211  sigmaEps_.readIfPresent(this->coeffDict());
212 
213  return true;
214  }
215  else
216  {
217  return false;
218  }
219 }
220 
221 
222 template<class BasicTurbulenceModel>
224 {
225  if (!this->turbulence_)
226  {
227  return;
228  }
229 
230  // Local references
231  const alphaField& alpha = this->alpha_;
232  const rhoField& rho = this->rho_;
233  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
234  const volVectorField& U = this->U_;
235  volScalarField& nut = this->nut_;
236  fv::options& fvOptions(fv::options::New(this->mesh_));
237 
239 
241 
242  tmp<volTensorField> tgradU = fvc::grad(U);
243  volScalarField G(this->GName(), nut*(dev(twoSymm(tgradU())) && tgradU()));
244  tgradU.clear();
245 
246  // Update epsilon and G at the wall
247  epsilon_.boundaryField().updateCoeffs();
248 
249  // Dissipation equation
250  tmp<fvScalarMatrix> epsEqn
251  (
252  fvm::ddt(alpha, rho, epsilon_)
253  + fvm::div(alphaRhoPhi, epsilon_)
254  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
255  ==
256  C1_*alpha*rho*G*epsilon_/k_
257  - fvm::SuSp(((2.0/3.0)*C1_ + C3_)*alpha*rho*divU, epsilon_)
258  - fvm::Sp(C2_*alpha*rho*epsilon_/k_, epsilon_)
259  + epsilonSource()
260  + fvOptions(alpha, rho, epsilon_)
261  );
262 
263  epsEqn().relax();
264  fvOptions.constrain(epsEqn());
265  epsEqn().boundaryManipulate(epsilon_.boundaryField());
266  solve(epsEqn);
267  fvOptions.correct(epsilon_);
268  bound(epsilon_, this->epsilonMin_);
269 
270  // Turbulent kinetic energy equation
272  (
273  fvm::ddt(alpha, rho, k_)
274  + fvm::div(alphaRhoPhi, k_)
275  - fvm::laplacian(alpha*rho*DkEff(), k_)
276  ==
277  alpha*rho*G
278  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
279  - fvm::Sp(alpha*rho*epsilon_/k_, k_)
280  + kSource()
281  + fvOptions(alpha, rho, k_)
282  );
283 
284  kEqn().relax();
285  fvOptions.constrain(kEqn());
286  solve(kEqn);
287  fvOptions.correct(k_);
288  bound(k_, this->kMin_);
289 
290  correctNut();
291 }
292 
293 
294 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
295 
296 } // End namespace RASModels
297 } // End namespace Foam
298 
299 // ************************************************************************* //
Foam::RASModels::kEpsilon::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kEpsilon.C:223
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
fvOptions.H
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::RASModel
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:103
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
kEpsilon.H
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::ThermalDiffusivity::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: ThermalDiffusivity.H:55
U
U
Definition: pEqn.H:46
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::RASModels::kEpsilon::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kEpsilon.C:202
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
nut
nut
Definition: createTDFields.H:71
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
bound.H
Bound the given scalar field if it has gone unbounded.
correct
fvOptions correct(rho)
Foam::fvm::SuSp
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::fv::options
Finite-volume options.
Definition: fvOptions.H:52
Foam::RASModels::kEpsilon::correctNut
virtual void correctNut()
Definition: kEpsilon.C:40
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
rho
rho
Definition: pEqn.H:3
Foam::RASModels::kEpsilon::epsilonSource
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: kEpsilon.C:68
Foam::ThermalDiffusivity::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: ThermalDiffusivity.H:57
Foam::RASModels::kEpsilon::kEpsilon
kEpsilon(const kEpsilon &)
Foam::solve
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
divU
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::eddyViscosity
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
Foam::RASModels::kEpsilon::kSource
virtual tmp< fvScalarMatrix > kSource() const
Definition: kEpsilon.C:53
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Foam::tmp::clear
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:172
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:93
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:187
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:16
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:104