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-2023 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  A0_
162  (
163  dimensioned<scalar>::lookupOrAddToDict
164  (
165  "A0",
166  this->coeffDict_,
167  4.0
168  )
169  ),
170  C2_
171  (
172  dimensioned<scalar>::lookupOrAddToDict
173  (
174  "C2",
175  this->coeffDict_,
176  1.9
177  )
178  ),
179  sigmak_
180  (
181  dimensioned<scalar>::lookupOrAddToDict
182  (
183  "sigmak",
184  this->coeffDict_,
185  1.0
186  )
187  ),
188  sigmaEps_
189  (
190  dimensioned<scalar>::lookupOrAddToDict
191  (
192  "sigmaEps",
193  this->coeffDict_,
194  1.2
195  )
196  ),
197 
198  k_
199  (
200  IOobject
201  (
202  this->groupName("k"),
203  this->runTime_.name(),
204  this->mesh_,
205  IOobject::MUST_READ,
206  IOobject::AUTO_WRITE
207  ),
208  this->mesh_
209  ),
210  epsilon_
211  (
212  IOobject
213  (
214  this->groupName("epsilon"),
215  this->runTime_.name(),
216  this->mesh_,
217  IOobject::MUST_READ,
218  IOobject::AUTO_WRITE
219  ),
220  this->mesh_
221  )
222 {
223  bound(k_, this->kMin_);
224  boundEpsilon();
225 
226  if (type == typeName)
227  {
228  this->printCoeffs(type);
229  }
230 }
231 
232 
233 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
234 
235 template<class BasicMomentumTransportModel>
237 {
239  {
240  A0_.readIfPresent(this->coeffDict());
241  C2_.readIfPresent(this->coeffDict());
242  sigmak_.readIfPresent(this->coeffDict());
243  sigmaEps_.readIfPresent(this->coeffDict());
244 
245  return true;
246  }
247  else
248  {
249  return false;
250  }
251 }
252 
253 
254 template<class BasicMomentumTransportModel>
256 {
257  if (!this->turbulence_)
258  {
259  return;
260  }
261 
262  // Local references
263  const alphaField& alpha = this->alpha_;
264  const rhoField& rho = this->rho_;
265  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
266  const volVectorField& U = this->U_;
267  volScalarField& nut = this->nut_;
268  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
270  (
271  Foam::fvConstraints::New(this->mesh_)
272  );
273 
275 
277  (
278  typedName("divU"),
279  fvc::div(fvc::absolute(this->phi(), U))()
280  );
281 
282  const volTensorField gradU(fvc::grad(U));
283  const volScalarField S2(typedName("S2"), 2*magSqr(dev(symm(gradU))));
284  const volScalarField magS(typedName("magS"), sqrt(S2));
285 
286  const volScalarField::Internal eta
287  (
288  typedName("eta"), magS()*k_()/epsilon_()
289  );
290  const volScalarField::Internal C1
291  (
292  typedName("C1"),
293  max(eta/(scalar(5) + eta), scalar(0.43))
294  );
295 
297  (
298  this->GName(),
299  nut*(gradU.v() && dev(twoSymm(gradU.v())))
300  );
301 
302  // Update epsilon and G at the wall
303  epsilon_.boundaryFieldRef().updateCoeffs();
304 
305  // Dissipation equation
306  tmp<fvScalarMatrix> epsEqn
307  (
308  fvm::ddt(alpha, rho, epsilon_)
309  + fvm::div(alphaRhoPhi, epsilon_)
310  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
311  ==
312  C1*alpha()*rho()*magS()*epsilon_()
313  - fvm::Sp
314  (
315  C2_*alpha()*rho()*epsilon_()/(k_() + sqrt(this->nu()()*epsilon_())),
316  epsilon_
317  )
318  + epsilonSource()
319  + fvModels.source(alpha, rho, epsilon_)
320  );
321 
322  epsEqn.ref().relax();
323  fvConstraints.constrain(epsEqn.ref());
324  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
325  solve(epsEqn);
326  fvConstraints.constrain(epsilon_);
327  boundEpsilon();
328 
329 
330  // Turbulent kinetic energy equation
331 
333  (
334  fvm::ddt(alpha, rho, k_)
335  + fvm::div(alphaRhoPhi, k_)
336  - fvm::laplacian(alpha*rho*DkEff(), k_)
337  ==
338  alpha()*rho()*G
339  - fvm::SuSp(2.0/3.0*alpha()*rho()*divU, k_)
340  - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
341  + kSource()
342  + fvModels.source(alpha, rho, k_)
343  );
344 
345  kEqn.ref().relax();
346  fvConstraints.constrain(kEqn.ref());
347  solve(kEqn);
349  bound(k_, this->kMin_);
350 
351  correctNut(gradU, S2, magS);
352 }
353 
354 
355 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
356 
357 } // End namespace RASModels
358 } // End namespace Foam
359 
360 // ************************************************************************* //
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:31
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:255
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:236
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:100
Dimension set for the base types.
Definition: dimensionSet.H:125
Generic dimensioned Type class.
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:237
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
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 SuType &Su)
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.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedSymmTensor sqr(const dimensionedVector &dv)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
scalarList W(const fluidMulticomponentThermo &thermo)
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimVolume
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:176
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
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedTensor skew(const dimensionedTensor &dt)
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.