kEqn.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 "kEqn.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 namespace LESModels
33 {
34 
35 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
36 
37 template<class BasicTurbulenceModel>
39 {
40  this->nut_ = Ck_*sqrt(k_)*this->delta();
41  this->nut_.correctBoundaryConditions();
42 
43  BasicTurbulenceModel::correctNut();
44 }
45 
46 
47 template<class BasicTurbulenceModel>
49 {
50  return tmp<fvScalarMatrix>
51  (
52  new fvScalarMatrix
53  (
54  k_,
55  dimVolume*this->rho_.dimensions()*k_.dimensions()
56  /dimTime
57  )
58  );
59 }
60 
61 
62 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
63 
64 template<class BasicTurbulenceModel>
66 (
67  const alphaField& alpha,
68  const rhoField& rho,
69  const volVectorField& U,
70  const surfaceScalarField& alphaRhoPhi,
71  const surfaceScalarField& phi,
72  const transportModel& transport,
73  const word& propertiesName,
74  const word& type
75 )
76 :
78  (
79  type,
80  alpha,
81  rho,
82  U,
83  alphaRhoPhi,
84  phi,
85  transport,
86  propertiesName
87  ),
88 
89  k_
90  (
91  IOobject
92  (
93  IOobject::groupName("k", this->U_.group()),
94  this->runTime_.timeName(),
95  this->mesh_,
98  ),
99  this->mesh_
100  ),
101 
102  Ck_
103  (
105  (
106  "Ck",
107  this->coeffDict_,
108  0.094
109  )
110  )
111 {
112  bound(k_, this->kMin_);
113 
114  if (type == typeName)
115  {
116  correctNut();
117  this->printCoeffs(type);
118  }
119 }
120 
121 
122 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
123 
124 template<class BasicTurbulenceModel>
126 {
128  {
129  Ck_.readIfPresent(this->coeffDict());
130 
131  return true;
132  }
133  else
134  {
135  return false;
136  }
137 }
138 
139 
140 template<class BasicTurbulenceModel>
142 {
143  return tmp<volScalarField>
144  (
145  new volScalarField
146  (
147  IOobject
148  (
149  IOobject::groupName("epsilon", this->U_.group()),
150  this->runTime_.timeName(),
151  this->mesh_,
154  ),
155  this->Ce_*k()*sqrt(k())/this->delta()
156  )
157  );
158 }
159 
160 
161 template<class BasicTurbulenceModel>
163 {
164  if (!this->turbulence_)
165  {
166  return;
167  }
168 
169  // Local references
170  const alphaField& alpha = this->alpha_;
171  const rhoField& rho = this->rho_;
172  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
173  const volVectorField& U = this->U_;
174  volScalarField& nut = this->nut_;
175 
177 
179 
180  tmp<volTensorField> tgradU(fvc::grad(U));
181  volScalarField G(this->GName(), nut*(tgradU() && dev(twoSymm(tgradU()))));
182  tgradU.clear();
183 
185  (
186  fvm::ddt(alpha, rho, k_)
187  + fvm::div(alphaRhoPhi, k_)
188  - fvm::laplacian(alpha*rho*DkEff(), k_)
189  ==
190  alpha*rho*G
191  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
192  - fvm::Sp(this->Ce_*alpha*rho*sqrt(k_)/this->delta(), k_)
193  + kSource()
194  );
195 
196  kEqn().relax();
197  solve(kEqn);
198  bound(k_, this->kMin_);
199 
200  correctNut();
201 }
202 
203 
204 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
205 
206 } // End namespace LESModels
207 } // End namespace Foam
208 
209 // ************************************************************************* //
dimensionedScalar sqrt(const dimensionedScalar &ds)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
One equation eddy-viscosity model.
Definition: kEqn.H:74
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
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
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
virtual tmp< fvScalarMatrix > kSource() const
Definition: kEqn.C:48
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
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
virtual void correct()
Correct eddy-Viscosity and related properties.
Definition: kEqn.C:162
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
static word groupName(Name name, const word &group)
BasicTurbulenceModel::transportModel transportModel
Definition: kEqn.H:109
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
Generic dimensioned Type class.
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 > &)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
scalar nut
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual bool read()
Read model coefficients if they have changed.
Definition: kEqn.C:125
label k
Boltzmann constant.
scalar delta
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
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
virtual void correctNut()
Definition: kEqn.C:38
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Eddy viscosity LES SGS model base class.
virtual tmp< volScalarField > epsilon() const
Return sub-grid disipation rate.
Definition: kEqn.C:141
BasicTurbulenceModel::rhoField rhoField
Definition: kEqn.H:108
BasicTurbulenceModel::alphaField alphaField
Definition: kEqn.H:107
A class for managing temporary objects.
Definition: PtrList.H:118
const dimensionedScalar G
Newtonian constant of gravitation.
U
Definition: pEqn.H:82