kEpsilon.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-2024 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 "kEpsilon.H"
27 #include "fvModels.H"
28 #include "fvConstraints.H"
29 #include "bound.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace RASModels
36 {
37 
38 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
39 
40 template<class BasicMomentumTransportModel>
42 {
43  tmp<volScalarField> tCmuk2(Cmu_*sqr(k_));
44  epsilon_ = max(epsilon_, tCmuk2()/(this->nutMaxCoeff_*this->nu()));
45  return tCmuk2;
46 }
47 
48 
49 template<class BasicMomentumTransportModel>
51 {
52  this->nut_ = boundEpsilon()/epsilon_;
53  this->nut_.correctBoundaryConditions();
54  fvConstraints::New(this->mesh_).constrain(this->nut_);
55 }
56 
57 
58 template<class BasicMomentumTransportModel>
60 {
61  return tmp<fvScalarMatrix>
62  (
63  new fvScalarMatrix
64  (
65  k_,
66  dimVolume*this->rho_.dimensions()*k_.dimensions()
67  /dimTime
68  )
69  );
70 }
71 
72 
73 template<class BasicMomentumTransportModel>
75 {
76  return tmp<fvScalarMatrix>
77  (
78  new fvScalarMatrix
79  (
80  epsilon_,
81  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
82  /dimTime
83  )
84  );
85 }
86 
87 
88 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
89 
90 template<class BasicMomentumTransportModel>
92 (
93  const alphaField& alpha,
94  const rhoField& rho,
95  const volVectorField& U,
96  const surfaceScalarField& alphaRhoPhi,
97  const surfaceScalarField& phi,
98  const viscosity& viscosity,
99  const word& type
100 )
101 :
102  eddyViscosity<RASModel<BasicMomentumTransportModel>>
103  (
104  type,
105  alpha,
106  rho,
107  U,
108  alphaRhoPhi,
109  phi,
110  viscosity
111  ),
112 
113  Cmu_("Cmu", this->coeffDict(), 0.09),
114  C1_("C1", this->coeffDict(), 1.44),
115  C2_("C2", this->coeffDict(), 1.92),
116  C3_("C3", this->coeffDict(), 0),
117  sigmak_("sigmak", this->coeffDict(), 1.0),
118  sigmaEps_("sigmaEps", this->coeffDict(), 1.3),
119 
120  k_
121  (
122  IOobject
123  (
124  this->groupName("k"),
125  this->runTime_.name(),
126  this->mesh_,
127  IOobject::MUST_READ,
128  IOobject::AUTO_WRITE
129  ),
130  this->mesh_
131  ),
132  epsilon_
133  (
135  (
136  this->groupName("epsilon"),
137  this->runTime_.name(),
138  this->mesh_,
139  IOobject::MUST_READ,
140  IOobject::AUTO_WRITE
141  ),
142  this->mesh_
143  )
144 {
145  bound(k_, this->kMin_);
146  boundEpsilon();
147 }
148 
149 
150 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151 
152 template<class BasicMomentumTransportModel>
154 {
156  {
157  Cmu_.readIfPresent(this->coeffDict());
158  C1_.readIfPresent(this->coeffDict());
159  C2_.readIfPresent(this->coeffDict());
160  C3_.readIfPresent(this->coeffDict());
161  sigmak_.readIfPresent(this->coeffDict());
162  sigmaEps_.readIfPresent(this->coeffDict());
163 
164  return true;
165  }
166  else
167  {
168  return false;
169  }
170 }
171 
172 
173 template<class BasicMomentumTransportModel>
175 {
176  if (!this->turbulence_)
177  {
178  return;
179  }
180 
181  // Local references
182  const alphaField& alpha = this->alpha_;
183  const rhoField& rho = this->rho_;
184  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
185  const volVectorField& U = this->U_;
186  volScalarField& nut = this->nut_;
187  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
189  (
190  Foam::fvConstraints::New(this->mesh_)
191  );
192 
194 
196  (
197  fvc::div(fvc::absolute(this->phi(), U))()
198  );
199 
200  tmp<volTensorField> tgradU = fvc::grad(U);
202  (
203  this->GName(),
204  nut()*(dev(twoSymm(tgradU().v())) && tgradU().v())
205  );
206  tgradU.clear();
207 
208  // Update epsilon and G at the wall
209  epsilon_.boundaryFieldRef().updateCoeffs();
210 
211  // Dissipation equation
212  tmp<fvScalarMatrix> epsEqn
213  (
214  fvm::ddt(alpha, rho, epsilon_)
215  + fvm::div(alphaRhoPhi, epsilon_)
216  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
217  ==
218  C1_*alpha()*rho()*G*epsilon_()/k_()
219  - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*alpha()*rho()*divU, epsilon_)
220  - fvm::Sp(C2_*alpha()*rho()*epsilon_()/k_(), epsilon_)
221  + epsilonSource()
222  + fvModels.source(alpha, rho, epsilon_)
223  );
224 
225  epsEqn.ref().relax();
226  fvConstraints.constrain(epsEqn.ref());
227  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
228  solve(epsEqn);
229  fvConstraints.constrain(epsilon_);
230  boundEpsilon();
231 
232  // Turbulent kinetic energy equation
234  (
235  fvm::ddt(alpha, rho, k_)
236  + fvm::div(alphaRhoPhi, k_)
237  - fvm::laplacian(alpha*rho*DkEff(), k_)
238  ==
239  alpha()*rho()*G
240  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
241  - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
242  + kSource()
243  + fvModels.source(alpha, rho, k_)
244  );
245 
246  kEqn.ref().relax();
247  fvConstraints.constrain(kEqn.ref());
248  solve(kEqn);
250  bound(k_, this->kMin_);
251 
252  correctNut();
253 }
254 
255 
256 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257 
258 } // End namespace RASModels
259 } // End namespace Foam
260 
261 // ************************************************************************* //
Bound the given scalar field where it is below the specified minimum.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:56
virtual tmp< fvScalarMatrix > epsilonSource() const
Source term for the epsilon equation.
Definition: kEpsilon.C:74
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kEpsilon.C:174
tmp< volScalarField > boundEpsilon()
Bound epsilon and return Cmu*sqr(k) for nut.
Definition: kEpsilon.C:41
kEpsilon(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
Definition: kEpsilon.C:92
virtual void correctNut()
Correct the eddy-viscosity nut.
Definition: kEpsilon.C:50
virtual tmp< fvScalarMatrix > kSource() const
Source term for the k equation.
Definition: kEpsilon.C:59
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kEpsilon.C:153
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:103
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
Finite volume constraints.
Definition: fvConstraints.H:67
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Finite volume models.
Definition: fvModels.H:65
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
A class for managing temporary objects.
Definition: tmp.H:55
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:253
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
const scalar nut
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp)
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
void dev(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void twoSymm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
const dimensionSet dimTime
const dimensionSet dimVolume
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
bool bound(volScalarField &, const dimensionedScalar &min)
Bound the given scalar field where it is below the specified min value.
Definition: bound.C:31
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.