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 -------------------------------------------------------------------------------
11 License
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 
34 namespace Foam
35 {
36 namespace LESModels
37 {
38 
39 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
40 
41 template<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 
52 template<class BasicTurbulenceModel>
54 {
55  return tmp<fvScalarMatrix>
56  (
57  new fvScalarMatrix
58  (
59  k_,
60  dimVolume*this->rho_.dimensions()*k_.dimensions()
61  /dimTime
62  )
63  );
64 }
65 
66 
67 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
68 
69 template<class BasicTurbulenceModel>
71 (
72  const alphaField& alpha,
73  const rhoField& rho,
74  const volVectorField& U,
75  const surfaceScalarField& alphaRhoPhi,
76  const surfaceScalarField& phi,
77  const transportModel& transport,
78  const word& propertiesName,
79  const word& type
80 )
81 :
82  LESeddyViscosity<BasicTurbulenceModel>
83  (
84  type,
85  alpha,
86  rho,
87  U,
88  alphaRhoPhi,
89  phi,
90  transport,
91  propertiesName
92  ),
93 
94  k_
95  (
96  IOobject
97  (
98  IOobject::groupName("k", this->alphaRhoPhi_.group()),
99  this->runTime_.timeName(),
100  this->mesh_,
101  IOobject::MUST_READ,
102  IOobject::AUTO_WRITE
103  ),
104  this->mesh_
105  ),
106 
107  Ck_
108  (
109  dimensioned<scalar>::getOrAddToDict
110  (
111  "Ck",
112  this->coeffDict_,
113  0.094
114  )
115  )
116 {
117  bound(k_, this->kMin_);
118 
119  if (type == typeName)
120  {
121  this->printCoeffs(type);
122  }
123 }
124 
125 
126 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127 
128 template<class BasicTurbulenceModel>
130 {
132  {
133  Ck_.readIfPresent(this->coeffDict());
134 
135  return true;
136  }
137 
138  return false;
139 }
140 
141 
142 template<class BasicTurbulenceModel>
144 {
145  if (!this->turbulence_)
146  {
147  return;
148  }
149 
150  // Local references
151  const alphaField& alpha = this->alpha_;
152  const rhoField& rho = this->rho_;
153  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
154  const volVectorField& U = this->U_;
155  volScalarField& nut = this->nut_;
156  fv::options& fvOptions(fv::options::New(this->mesh_));
157 
159 
161 
162  tmp<volTensorField> tgradU(fvc::grad(U));
163  volScalarField G(this->GName(), nut*(tgradU() && devTwoSymm(tgradU())));
164  tgradU.clear();
165 
166  tmp<fvScalarMatrix> kEqn
167  (
168  fvm::ddt(alpha, rho, k_)
169  + fvm::div(alphaRhoPhi, k_)
170  - fvm::laplacian(alpha*rho*DkEff(), k_)
171  ==
172  alpha*rho*G
173  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
174  - fvm::Sp(this->Ce_*alpha*rho*sqrt(k_)/this->delta(), k_)
175  + kSource()
176  + fvOptions(alpha, rho, k_)
177  );
178 
179  kEqn.ref().relax();
180  fvOptions.constrain(kEqn.ref());
181  solve(kEqn);
182  fvOptions.correct(k_);
183  bound(k_, this->kMin_);
184 
185  correctNut();
186 }
187 
188 
189 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
190 
191 } // End namespace LESModels
192 } // End namespace Foam
193 
194 // ************************************************************************* //
Eddy viscosity LES SGS model base class.
scalar delta
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: LESModel.C:27
Generic dimensioned Type class.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
fv::options & fvOptions
constexpr const char *const group
Group name for atomic constants.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
word timeName
Definition: getTimeIndex.H:3
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
volScalarField k_
Definition: kEqn.H:83
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual void correct()
Correct eddy-Viscosity and related properties.
Definition: kEqn.C:136
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.
virtual tmp< fvScalarMatrix > kSource() const
Definition: kEqn.C:46
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
virtual bool read()
Read model coefficients if they have changed.
Definition: kEqn.C:122
virtual bool read()
Read model coefficients if they have changed.
kEqn(const kEqn &)=delete
No copy construct.
virtual void correctNut()
Definition: kEqn.C:35
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:183
U
Definition: pEqn.H:72
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:29
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
Base-class for all transport models used by the incompressible turbulence models. ...
dimensionedScalar kMin_
Lower limit of k.
Definition: LESModel.H:92
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
zeroField divU
Definition: alphaSuSp.H:3
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:96
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
scalar nut
Namespace for OpenFOAM.
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
Definition: SymmTensorI.H:491