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