realizableKE.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 "realizableKE.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  epsilon_ = max(epsilon_, 0.09*sqr(k_)/(this->nutMaxCoeff_*this->nu()));
44 }
45 
46 
47 template<class BasicMomentumTransportModel>
49 (
50  const volTensorField& gradU,
51  const volScalarField& S2,
52  const volScalarField& magS
53 )
54 {
55  tmp<volSymmTensorField> tS = dev(symm(gradU));
56  const volSymmTensorField& S = tS();
57 
58  const volScalarField W
59  (
60  (2*sqrt(2.0))*((S&S)&&S)
61  /(
62  magS*S2
63  + dimensionedScalar(dimensionSet(0, 0, -3, 0, 0), small)
64  )
65  );
66 
67  tS.clear();
68 
69  const volScalarField phis
70  (
71  (1.0/3.0)*acos(min(max(sqrt(6.0)*W, -scalar(1)), scalar(1)))
72  );
73  const volScalarField As(sqrt(6.0)*cos(phis));
74  const volScalarField Us(sqrt(S2/2.0 + magSqr(skew(gradU))));
75 
76  return 1.0/(A0_ + As*Us*k_/epsilon_);
77 }
78 
79 
80 template<class BasicMomentumTransportModel>
82 (
83  const volTensorField& gradU,
84  const volScalarField& S2,
85  const volScalarField& magS
86 )
87 {
88  boundEpsilon();
89  this->nut_ = rCmu(gradU, S2, magS)*sqr(k_)/epsilon_;
90  this->nut_.correctBoundaryConditions();
91  fvConstraints::New(this->mesh_).constrain(this->nut_);
92 }
93 
94 
95 template<class BasicMomentumTransportModel>
97 {
98  const volTensorField gradU(fvc::grad(this->U_));
99  const volScalarField S2(typedName("S2"), 2*magSqr(dev(symm(gradU))));
100  const volScalarField magS(typedName("magS"), sqrt(S2));
101 
102  correctNut(gradU, S2, magS);
103 }
104 
105 
106 template<class BasicMomentumTransportModel>
108 {
109  return tmp<fvScalarMatrix>
110  (
111  new fvScalarMatrix
112  (
113  k_,
114  dimVolume*this->rho_.dimensions()*k_.dimensions()
115  /dimTime
116  )
117  );
118 }
119 
120 
121 template<class BasicMomentumTransportModel>
124 {
125  return tmp<fvScalarMatrix>
126  (
127  new fvScalarMatrix
128  (
129  epsilon_,
130  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
131  /dimTime
132  )
133  );
134 }
135 
136 
137 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
138 
139 template<class BasicMomentumTransportModel>
141 (
142  const alphaField& alpha,
143  const rhoField& rho,
144  const volVectorField& U,
145  const surfaceScalarField& alphaRhoPhi,
146  const surfaceScalarField& phi,
147  const viscosity& viscosity,
148  const word& type
149 )
150 :
151  eddyViscosity<RASModel<BasicMomentumTransportModel>>
152  (
153  type,
154  alpha,
155  rho,
156  U,
157  alphaRhoPhi,
158  phi,
159  viscosity
160  ),
161 
162  A0_("A0", this->coeffDict(), 4.0),
163  C2_
164  ("C2", this->coeffDict(), 1.9),
165  sigmak_("sigmak", this->coeffDict(), 1.0),
166  sigmaEps_("sigmaEps", this->coeffDict(), 1.2),
167 
168  k_
169  (
170  IOobject
171  (
172  this->groupName("k"),
173  this->runTime_.name(),
174  this->mesh_,
175  IOobject::MUST_READ,
176  IOobject::AUTO_WRITE
177  ),
178  this->mesh_
179  ),
180  epsilon_
181  (
182  IOobject
183  (
184  this->groupName("epsilon"),
185  this->runTime_.name(),
186  this->mesh_,
187  IOobject::MUST_READ,
188  IOobject::AUTO_WRITE
189  ),
190  this->mesh_
191  )
192 {
193  bound(k_, this->kMin_);
194  boundEpsilon();
195 }
196 
197 
198 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
199 
200 template<class BasicMomentumTransportModel>
202 {
204  {
205  A0_.readIfPresent(this->coeffDict());
206  C2_.readIfPresent(this->coeffDict());
207  sigmak_.readIfPresent(this->coeffDict());
208  sigmaEps_.readIfPresent(this->coeffDict());
209 
210  return true;
211  }
212  else
213  {
214  return false;
215  }
216 }
217 
218 
219 template<class BasicMomentumTransportModel>
221 {
222  if (!this->turbulence_)
223  {
224  return;
225  }
226 
227  // Local references
228  const alphaField& alpha = this->alpha_;
229  const rhoField& rho = this->rho_;
230  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
231  const volVectorField& U = this->U_;
232  volScalarField& nut = this->nut_;
233  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
235  (
236  Foam::fvConstraints::New(this->mesh_)
237  );
238 
240 
242  (
243  typedName("divU"),
244  fvc::div(fvc::absolute(this->phi(), U))()
245  );
246 
247  const volTensorField gradU(fvc::grad(U));
248  const volScalarField S2(typedName("S2"), 2*magSqr(dev(symm(gradU))));
249  const volScalarField magS(typedName("magS"), sqrt(S2));
250 
251  const volScalarField::Internal eta
252  (
253  typedName("eta"), magS()*k_()/epsilon_()
254  );
255  const volScalarField::Internal C1
256  (
257  typedName("C1"),
258  max(eta/(scalar(5) + eta), scalar(0.43))
259  );
260 
262  (
263  this->GName(),
264  nut*(gradU.v() && dev(twoSymm(gradU.v())))
265  );
266 
267  // Update epsilon and G at the wall
268  epsilon_.boundaryFieldRef().updateCoeffs();
269 
270  // Dissipation equation
271  tmp<fvScalarMatrix> epsEqn
272  (
273  fvm::ddt(alpha, rho, epsilon_)
274  + fvm::div(alphaRhoPhi, epsilon_)
275  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
276  ==
277  C1*alpha()*rho()*magS()*epsilon_()
278  - fvm::Sp
279  (
280  C2_*alpha()*rho()*epsilon_()/(k_() + sqrt(this->nu()()*epsilon_())),
281  epsilon_
282  )
283  + epsilonSource()
284  + fvModels.source(alpha, rho, epsilon_)
285  );
286 
287  epsEqn.ref().relax();
288  fvConstraints.constrain(epsEqn.ref());
289  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
290  solve(epsEqn);
291  fvConstraints.constrain(epsilon_);
292  boundEpsilon();
293 
294 
295  // Turbulent kinetic energy equation
296 
298  (
299  fvm::ddt(alpha, rho, k_)
300  + fvm::div(alphaRhoPhi, k_)
301  - fvm::laplacian(alpha*rho*DkEff(), k_)
302  ==
303  alpha()*rho()*G
304  - fvm::SuSp(2.0/3.0*alpha()*rho()*divU, k_)
305  - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
306  + kSource()
307  + fvModels.source(alpha, rho, k_)
308  );
309 
310  kEqn.ref().relax();
311  fvConstraints.constrain(kEqn.ref());
312  solve(kEqn);
314  bound(k_, this->kMin_);
315 
316  correctNut(gradU, S2, magS);
317 }
318 
319 
320 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
321 
322 } // End namespace RASModels
323 } // End namespace Foam
324 
325 // ************************************************************************* //
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.
const Internal & v() const
Return a const-reference to the dimensioned internal field.
Definition: volFieldsI.H:29
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
tmp< volScalarField > rCmu(const volTensorField &gradU, const volScalarField &S2, const volScalarField &magS)
Definition: realizableKE.C:49
virtual tmp< fvScalarMatrix > epsilonSource() const
Source term for the epsilon equation.
Definition: realizableKE.C:123
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: realizableKE.C:220
void boundEpsilon()
Bound epsilon.
Definition: realizableKE.C:41
virtual void correctNut()
Correct the eddy-viscosity nut.
Definition: realizableKE.C:96
virtual tmp< fvScalarMatrix > kSource() const
Source term for the k equation.
Definition: realizableKE.C:107
realizableKE(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: realizableKE.C:141
virtual bool read()
Re-read model coefficients if they have changed.
Definition: realizableKE.C:201
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:103
Dimension set for the base types.
Definition: dimensionSet.H:125
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
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
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 > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
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.
void skew(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
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)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
void symm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
scalarList W(const fluidMulticomponentThermo &thermo)
const dimensionSet dimVolume
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:181
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.
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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.