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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "kEqn.H"
29 #include "fvOptions.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace LESModels
36 {
37 
38 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
39 
40 template<class BasicTurbulenceModel>
42 {
43  this->nut_ = Ck_*sqrt(k_)*this->delta();
44  this->nut_.correctBoundaryConditions();
45  fv::options::New(this->mesh_).correct(this->nut_);
46 
47  BasicTurbulenceModel::correctNut();
48 }
49 
50 
51 template<class BasicTurbulenceModel>
53 {
54  return tmp<fvScalarMatrix>
55  (
56  new fvScalarMatrix
57  (
58  k_,
59  dimVolume*this->rho_.dimensions()*k_.dimensions()
60  /dimTime
61  )
62  );
63 }
64 
65 
66 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
67 
68 template<class BasicTurbulenceModel>
70 (
71  const alphaField& alpha,
72  const rhoField& rho,
73  const volVectorField& U,
74  const surfaceScalarField& alphaRhoPhi,
75  const surfaceScalarField& phi,
76  const transportModel& transport,
77  const word& propertiesName,
78  const word& type
79 )
80 :
82  (
83  type,
84  alpha,
85  rho,
86  U,
87  alphaRhoPhi,
88  phi,
89  transport,
90  propertiesName
91  ),
92 
93  k_
94  (
95  IOobject
96  (
97  IOobject::groupName("k", this->alphaRhoPhi_.group()),
98  this->runTime_.timeName(),
99  this->mesh_,
102  ),
103  this->mesh_
104  ),
105 
106  Ck_
107  (
109  (
110  "Ck",
111  this->coeffDict_,
112  0.094
113  )
114  )
115 {
116  bound(k_, this->kMin_);
117 
118  if (type == typeName)
119  {
120  this->printCoeffs(type);
121  }
122 }
123 
124 
125 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
126 
127 template<class BasicTurbulenceModel>
129 {
131  {
132  Ck_.readIfPresent(this->coeffDict());
133 
134  return true;
135  }
136 
137  return false;
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->alphaRhoPhi_.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 
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 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
fvOptions.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fv::optionList::correct
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
Definition: fvOptionListTemplates.C:321
Foam::tmp::clear
void clear() const noexcept
Definition: tmpI.H:325
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:104
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:35
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
Foam::LESModels::kEqn::correctNut
virtual void correctNut()
Definition: kEqn.C:41
Foam::LESModels::kEqn::epsilon
virtual tmp< volScalarField > epsilon() const
Return sub-grid dissipation rate.
Definition: kEqn.C:142
Foam::LESModels::kEqn::correct
virtual void correct()
Correct eddy-Viscosity and related properties.
Definition: kEqn.C:163
Foam::LESModels::kEqn::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: kEqn.H:111
rho
rho
Definition: readInitialConditions.H:96
Foam::LESModels::kEqn
One equation eddy-viscosity model.
Definition: kEqn.H:77
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::fvm::SuSp
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::solve
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
divU
zeroField divU
Definition: alphaSuSp.H:3
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::LESModels::kEqn::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: kEqn.H:112
Foam::fv::options
Finite-volume options.
Definition: fvOptions.H:55
Foam::LESModels::LESeddyViscosity
Eddy viscosity LES SGS model base class.
Definition: LESeddyViscosity.H:58
kEqn.H
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:43
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::LESModels::kEqn::kSource
virtual tmp< fvScalarMatrix > kSource() const
Definition: kEqn.C:52
Foam::LESModels::kEqn::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: kEqn.H:110
U
U
Definition: pEqn.H:72
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::IOobject::groupName
static word groupName(StringType name, const word &group)
Create dot-delimited name.group.
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
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:76
Foam::eddyViscosity< LESModel< BasicTurbulenceModel > >::correct
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
Definition: eddyViscosity.C:134
Foam::LESModels::kEqn::kEqn
kEqn(const kEqn &)=delete
No copy construct.
Foam::LESModels::kEqn::read
virtual bool read()
Read model coefficients if they have changed.
Definition: kEqn.C:128
nut
scalar nut
Definition: evaluateNearWall.H:5
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:95
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:190
Foam::IOobject::MUST_READ
Definition: IOobject.H:120
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106