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-2016 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 #include "fvOptions.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace LESModels
34 {
35 
36 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
37 
38 template<class BasicTurbulenceModel>
40 {
41  this->nut_ = Ck_*sqrt(k_)*this->delta();
42  this->nut_.correctBoundaryConditions();
43  fv::options::New(this->mesh_).correct(this->nut_);
44 
45  BasicTurbulenceModel::correctNut();
46 }
47 
48 
49 template<class BasicTurbulenceModel>
51 {
52  return tmp<fvScalarMatrix>
53  (
54  new fvScalarMatrix
55  (
56  k_,
57  dimVolume*this->rho_.dimensions()*k_.dimensions()
58  /dimTime
59  )
60  );
61 }
62 
63 
64 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
65 
66 template<class BasicTurbulenceModel>
68 (
69  const alphaField& alpha,
70  const rhoField& rho,
71  const volVectorField& U,
72  const surfaceScalarField& alphaRhoPhi,
73  const surfaceScalarField& phi,
74  const transportModel& transport,
75  const word& propertiesName,
76  const word& type
77 )
78 :
80  (
81  type,
82  alpha,
83  rho,
84  U,
85  alphaRhoPhi,
86  phi,
87  transport,
88  propertiesName
89  ),
90 
91  k_
92  (
93  IOobject
94  (
95  IOobject::groupName("k", this->U_.group()),
96  this->runTime_.timeName(),
97  this->mesh_,
100  ),
101  this->mesh_
102  ),
103 
104  Ck_
105  (
107  (
108  "Ck",
109  this->coeffDict_,
110  0.094
111  )
112  )
113 {
114  bound(k_, this->kMin_);
115 
116  if (type == typeName)
117  {
118  this->printCoeffs(type);
119  }
120 }
121 
122 
123 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
124 
125 template<class BasicTurbulenceModel>
127 {
129  {
130  Ck_.readIfPresent(this->coeffDict());
131 
132  return true;
133  }
134  else
135  {
136  return false;
137  }
138 }
139 
140 
141 template<class BasicTurbulenceModel>
143 {
144  return tmp<volScalarField>
145  (
146  new volScalarField
147  (
148  IOobject
149  (
150  IOobject::groupName("epsilon", this->U_.group()),
151  this->runTime_.timeName(),
152  this->mesh_,
155  ),
156  this->Ce_*k()*sqrt(k())/this->delta()
157  )
158  );
159 }
160 
161 
162 template<class BasicTurbulenceModel>
164 {
165  if (!this->turbulence_)
166  {
167  return;
168  }
169 
170  // Local references
171  const alphaField& alpha = this->alpha_;
172  const rhoField& rho = this->rho_;
173  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
174  const volVectorField& U = this->U_;
175  volScalarField& nut = this->nut_;
176  fv::options& fvOptions(fv::options::New(this->mesh_));
177 
179 
181 
182  tmp<volTensorField> tgradU(fvc::grad(U));
183  volScalarField G(this->GName(), nut*(tgradU() && dev(twoSymm(tgradU()))));
184  tgradU.clear();
185 
187  (
188  fvm::ddt(alpha, rho, k_)
189  + fvm::div(alphaRhoPhi, k_)
190  - fvm::laplacian(alpha*rho*DkEff(), k_)
191  ==
192  alpha*rho*G
193  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
194  - fvm::Sp(this->Ce_*alpha*rho*sqrt(k_)/this->delta(), k_)
195  + kSource()
196  + fvOptions(alpha, rho, k_)
197  );
198 
199  kEqn.ref().relax();
200  fvOptions.constrain(kEqn.ref());
201  solve(kEqn);
202  fvOptions.correct(k_);
203  bound(k_, this->kMin_);
204 
205  correctNut();
206 }
207 
208 
209 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
210 
211 } // End namespace LESModels
212 } // End namespace Foam
213 
214 // ************************************************************************* //
Eddy viscosity LES SGS model base class.
scalar delta
surfaceScalarField & phi
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
U
Definition: pEqn.H:83
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
dimensionedScalar sqrt(const dimensionedScalar &ds)
label k
Boltzmann constant.
Generic dimensioned Type class.
Finite-volume options.
Definition: fvOptions.H:52
One equation eddy-viscosity model.
Definition: kEqn.H:74
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
fv::options & fvOptions
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:230
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
A class for handling words, derived from string.
Definition: word.H:59
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
static word groupName(Name name, const word &group)
BasicTurbulenceModel::transportModel transportModel
Definition: kEqn.H:109
virtual void correct()
Correct eddy-Viscosity and related properties.
Definition: kEqn.C:163
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
virtual bool read()
Read model coefficients if they have changed.
Definition: kEqn.C:126
BasicTurbulenceModel::alphaField alphaField
Definition: kEqn.H:107
virtual void correctNut()
Definition: kEqn.C:39
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:188
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
virtual tmp< fvScalarMatrix > kSource() const
Definition: kEqn.C:50
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
A class for managing temporary objects.
Definition: PtrList.H:54
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:103
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
virtual tmp< volScalarField > epsilon() const
Return sub-grid disipation rate.
Definition: kEqn.C:142
scalar nut
BasicTurbulenceModel::rhoField rhoField
Definition: kEqn.H:108
Namespace for OpenFOAM.