Loading...
Searching...
No Matches
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 Copyright (C) 2020-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "kEqn.H"
30#include "fvOptions.H"
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34namespace Foam
36namespace LESModels
37{
38
39// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
40
41template<class BasicTurbulenceModel>
43{
44 this->nut_ = Ck_*sqrt(k_)*this->delta();
45 this->nut_.correctBoundaryConditions();
46 fv::options::New(this->mesh_).correct(this->nut_);
47
48 BasicTurbulenceModel::correctNut();
49}
50
51
52template<class BasicTurbulenceModel>
54{
56 (
57 k_,
58 dimVolume*this->rho_.dimensions()*k_.dimensions()/dimTime
59 );
60}
61
62
63// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
64
65template<class BasicTurbulenceModel>
67(
68 const alphaField& alpha,
69 const rhoField& rho,
70 const volVectorField& U,
71 const surfaceScalarField& alphaRhoPhi,
73 const transportModel& transport,
74 const word& propertiesName,
75 const word& type
76)
77:
78 LESeddyViscosity<BasicTurbulenceModel>
79 (
80 type,
81 alpha,
82 rho,
83 U,
84 alphaRhoPhi,
85 phi,
86 transport,
87 propertiesName
88 ),
89
90 k_
91 (
93 (
94 IOobject::groupName("k", this->alphaRhoPhi_.group()),
95 this->runTime_.timeName(),
96 this->mesh_,
97 IOobject::MUST_READ,
98 IOobject::AUTO_WRITE
99 ),
100 this->mesh_
101 ),
102
103 Ck_
104 (
105 dimensioned<scalar>::getOrAddToDict
106 (
107 "Ck",
108 this->coeffDict_,
109 0.094
110 )
111 )
112{
113 bound(k_, this->kMin_);
114
115 if (type == typeName)
116 {
117 this->printCoeffs(type);
119}
120
121
122// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
123
124template<class BasicTurbulenceModel>
126{
128 {
129 Ck_.readIfPresent(this->coeffDict());
130
131 return true;
133
134 return false;
135}
136
137
138template<class BasicTurbulenceModel>
140{
141 if (!this->turbulence_)
142 {
143 return;
144 }
145
146 // Local references
147 const alphaField& alpha = this->alpha_;
148 const rhoField& rho = this->rho_;
149 const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
150 const volVectorField& U = this->U_;
151 volScalarField& nut = this->nut_;
152 fv::options& fvOptions(fv::options::New(this->mesh_));
153
155
157
158 tmp<volTensorField> tgradU(fvc::grad(U));
159 volScalarField G(this->GName(), nut*(tgradU() && devTwoSymm(tgradU())));
160 tgradU.clear();
161
162 tmp<fvScalarMatrix> kEqn
163 (
164 fvm::ddt(alpha, rho, k_)
165 + fvm::div(alphaRhoPhi, k_)
166 - fvm::laplacian(alpha*rho*DkEff(), k_)
167 ==
168 alpha*rho*G
169 - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
170 - fvm::Sp(this->Ce_*alpha*rho*sqrt(k_)/this->delta(), k_)
171 + kSource()
172 + fvOptions(alpha, rho, k_)
173 );
174
175 kEqn.ref().relax();
176 fvOptions.constrain(kEqn.ref());
177 solve(kEqn);
178 fvOptions.correct(k_);
179 bound(k_, this->kMin_);
180
181 correctNut();
182}
183
184
185// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
186
187} // End namespace LESModels
188} // End namespace Foam
189
190// ************************************************************************* //
scalar delta
fv::options & fvOptions
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Eddy viscosity LES SGS model base class.
virtual bool read()
Read model coefficients if they have changed.
BasicTurbulenceModel::alphaField alphaField
Definition kEqn.H:110
kEqn(const kEqn &)=delete
No copy construct.
BasicTurbulenceModel::rhoField rhoField
Definition kEqn.H:111
volScalarField k_
Definition kEqn.H:83
dimensionedScalar Ck_
Definition kEqn.H:88
virtual void correct()
Correct eddy-Viscosity and related properties.
Definition kEqn.C:132
virtual void correctNut()
Definition kEqn.C:35
virtual tmp< fvScalarMatrix > kSource() const
Definition kEqn.C:46
BasicTurbulenceModel::transportModel transportModel
Definition kEqn.H:112
virtual bool read()
Read model coefficients if they have changed.
Definition kEqn.C:118
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition kEqn.H:163
Generic dimensioned Type class.
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
Finite-volume options, which is an IOdictionary of values and a fv::optionList.
Definition fvOptions.H:69
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
Definition fvOptions.C:116
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
A class for managing temporary objects.
Definition tmp.H:75
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition tmpI.H:289
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
Base-class for all transport models used by the incompressible turbulence models.
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
scalar nut
zeroField divU
Definition alphaSuSp.H:3
word timeName
Definition getTimeIndex.H:3
Namespace for LES SGS models.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition fvcMeshPhi.C:183
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvmDiv.C:41
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvmDdt.C:41
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition bound.C:29
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimVolume(pow3(dimLength))
volScalarField & alpha
CEqn solve()