All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
kEqn.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 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 BasicMomentumTransportModel>
40 {
41  this->nut_ = Ck_*sqrt(k_)*this->delta();
42  this->nut_.correctBoundaryConditions();
43  fv::options::New(this->mesh_).correct(this->nut_);
44 }
45 
46 
47 template<class BasicMomentumTransportModel>
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 BasicMomentumTransportModel>
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& type
74 )
75 :
77  (
78  type,
79  alpha,
80  rho,
81  U,
82  alphaRhoPhi,
83  phi,
84  transport
85  ),
86 
87  k_
88  (
89  IOobject
90  (
91  IOobject::groupName("k", this->alphaRhoPhi_.group()),
92  this->runTime_.timeName(),
93  this->mesh_,
96  ),
97  this->mesh_
98  ),
99 
100  Ck_
101  (
103  (
104  "Ck",
105  this->coeffDict_,
106  0.094
107  )
108  )
109 {
110  bound(k_, this->kMin_);
111 
112  if (type == typeName)
113  {
114  this->printCoeffs(type);
115  }
116 }
117 
118 
119 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120 
121 template<class BasicMomentumTransportModel>
123 {
125  {
126  Ck_.readIfPresent(this->coeffDict());
127 
128  return true;
129  }
130  else
131  {
132  return false;
133  }
134 }
135 
136 
137 template<class BasicMomentumTransportModel>
139 {
140  return volScalarField::New
141  (
142  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
143  this->Ce_*k()*sqrt(k())/this->delta()
144  );
145 }
146 
147 
148 template<class BasicMomentumTransportModel>
150 {
151  if (!this->turbulence_)
152  {
153  return;
154  }
155 
156  // Local references
157  const alphaField& alpha = this->alpha_;
158  const rhoField& rho = this->rho_;
159  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
160  const volVectorField& U = this->U_;
161  volScalarField& nut = this->nut_;
162  fv::options& fvOptions(fv::options::New(this->mesh_));
163 
165 
167 
168  tmp<volTensorField> tgradU(fvc::grad(U));
169  volScalarField G(this->GName(), nut*(tgradU() && dev(twoSymm(tgradU()))));
170  tgradU.clear();
171 
173  (
174  fvm::ddt(alpha, rho, k_)
175  + fvm::div(alphaRhoPhi, k_)
176  - fvm::laplacian(alpha*rho*DkEff(), k_)
177  ==
178  alpha*rho*G
179  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
180  - fvm::Sp(this->Ce_*alpha*rho*sqrt(k_)/this->delta(), k_)
181  + kSource()
182  + fvOptions(alpha, rho, k_)
183  );
184 
185  kEqn.ref().relax();
186  fvOptions.constrain(kEqn.ref());
187  solve(kEqn);
188  fvOptions.correct(k_);
189  bound(k_, this->kMin_);
190 
191  correctNut();
192 }
193 
194 
195 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196 
197 } // End namespace LESModels
198 } // End namespace Foam
199 
200 // ************************************************************************* //
Eddy viscosity LES SGS model base class.
virtual tmp< volScalarField > epsilon() const
Return sub-grid disipation rate.
Definition: kEqn.C:138
scalar delta
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
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:230
fv::options & fvOptions
kEqn(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Constructor from components.
Definition: kEqn.C:66
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
dimensionedScalar sqrt(const dimensionedScalar &ds)
label k
Boltzmann constant.
Generic dimensioned Type class.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Finite-volume options.
Definition: fvOptions.H:52
One equation eddy-viscosity model.
Definition: kEqn.H:71
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
phi
Definition: pEqn.H:104
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)
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
virtual void correct()
Correct eddy-Viscosity and related properties.
Definition: kEqn.C:149
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
virtual tmp< fvScalarMatrix > kSource() const
Definition: kEqn.C:48
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:122
BasicMomentumTransportModel::rhoField rhoField
Definition: kEqn.H:98
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
U
Definition: pEqn.H:72
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
zeroField divU
Definition: alphaSuSp.H:3
A class for managing temporary objects.
Definition: PtrList.H:53
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:101
BasicMomentumTransportModel::transportModel transportModel
Definition: kEqn.H:99
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
const dimensionedScalar & G
Newtonian constant of gravitation.
BasicMomentumTransportModel::alphaField alphaField
Definition: kEqn.H:97
scalar nut
Namespace for OpenFOAM.