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  this->printCoeffs(type);
117  }
118 }
119 
120 
121 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
122 
123 template<class BasicTurbulenceModel>
125 {
127  {
128  Ck_.readIfPresent(this->coeffDict());
129 
130  return true;
131  }
132  else
133  {
134  return false;
135  }
136 }
137 
138 
139 template<class BasicTurbulenceModel>
141 {
142  return tmp<volScalarField>
143  (
144  new volScalarField
145  (
146  IOobject
147  (
148  IOobject::groupName("epsilon", this->U_.group()),
149  this->runTime_.timeName(),
150  this->mesh_,
153  ),
154  this->Ce_*k()*sqrt(k())/this->delta()
155  )
156  );
157 }
158 
159 
160 template<class BasicTurbulenceModel>
162 {
163  if (!this->turbulence_)
164  {
165  return;
166  }
167 
168  // Local references
169  const alphaField& alpha = this->alpha_;
170  const rhoField& rho = this->rho_;
171  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
172  const volVectorField& U = this->U_;
173  volScalarField& nut = this->nut_;
174 
176 
178 
180  volScalarField G(this->GName(), nut*(tgradU() && dev(twoSymm(tgradU()))));
181  tgradU.clear();
182 
184  (
185  fvm::ddt(alpha, rho, k_)
186  + fvm::div(alphaRhoPhi, k_)
187  - fvm::laplacian(alpha*rho*DkEff(), k_)
188  ==
189  alpha*rho*G
190  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
191  - fvm::Sp(this->Ce_*alpha*rho*sqrt(k_)/this->delta(), k_)
192  + kSource()
193  );
194 
195  kEqn().relax();
196  solve(kEqn);
197  bound(k_, this->kMin_);
198 
199  correctNut();
200 }
201 
202 
203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204 
205 } // End namespace LESModels
206 } // End namespace Foam
207 
208 // ************************************************************************* //
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
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
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::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::LESModels::kEqn::correctNut
virtual void correctNut()
Definition: kEqn.C:38
Foam::LESModels::kEqn::epsilon
virtual tmp< volScalarField > epsilon() const
Return sub-grid disipation rate.
Definition: kEqn.C:140
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::LESModels::kEqn::correct
virtual void correct()
Correct eddy-Viscosity and related properties.
Definition: kEqn.C:161
Foam::LESModels::kEqn::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: kEqn.H:108
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
U
U
Definition: pEqn.H:46
Foam::LESModels::kEqn
One equation eddy-viscosity model.
Definition: kEqn.H:74
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
nut
nut
Definition: createTDFields.H:71
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::fvm::SuSp
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::LESModels::kEqn::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: kEqn.H:109
Foam::LESModels::LESeddyViscosity
Eddy viscosity LES SGS model base class.
Definition: LESeddyViscosity.H:55
kEqn.H
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::LESModels::kEqn::kSource
virtual tmp< fvScalarMatrix > kSource() const
Definition: kEqn.C:48
rho
rho
Definition: pEqn.H:3
Foam::LESModels::kEqn::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: kEqn.H:107
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
divU
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::eddyViscosity< LESModel< BasicTurbulenceModel > >::correct
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
Definition: eddyViscosity.C:132
Foam::LESModels::kEqn::kEqn
kEqn(const kEqn &)
Foam::LESModels::kEqn::read
virtual bool read()
Read model coefficients if they have changed.
Definition: kEqn.C:124
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
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:104