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-2018 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 "fvOptions.H"
28 #include "bound.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace RASModels
35 {
36 
37 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
38 
39 template<class BasicTurbulenceModel>
41 (
42  const volTensorField& gradU,
43  const volScalarField& S2,
44  const volScalarField& magS
45 )
46 {
47  tmp<volSymmTensorField> tS = dev(symm(gradU));
48  const volSymmTensorField& S = tS();
49 
51  (
52  (2*sqrt(2.0))*((S&S)&&S)
53  /(
54  magS*S2
55  + dimensionedScalar("small", dimensionSet(0, 0, -3, 0, 0), small)
56  )
57  );
58 
59  tS.clear();
60 
61  volScalarField phis
62  (
63  (1.0/3.0)*acos(min(max(sqrt(6.0)*W, -scalar(1)), scalar(1)))
64  );
65  volScalarField As(sqrt(6.0)*cos(phis));
66  volScalarField Us(sqrt(S2/2.0 + magSqr(skew(gradU))));
67 
68  return 1.0/(A0_ + As*Us*k_/epsilon_);
69 }
70 
71 
72 template<class BasicTurbulenceModel>
74 (
75  const volTensorField& gradU,
76  const volScalarField& S2,
77  const volScalarField& magS
78 )
79 {
80  this->nut_ = rCmu(gradU, S2, magS)*sqr(k_)/epsilon_;
81  this->nut_.correctBoundaryConditions();
82  fv::options::New(this->mesh_).correct(this->nut_);
83 
84  BasicTurbulenceModel::correctNut();
85 }
86 
87 
88 template<class BasicTurbulenceModel>
90 {
91  tmp<volTensorField> tgradU = fvc::grad(this->U_);
92  volScalarField S2(2*magSqr(dev(symm(tgradU()))));
93  volScalarField magS(sqrt(S2));
94  correctNut(tgradU(), S2, magS);
95 }
96 
97 
98 template<class BasicTurbulenceModel>
100 {
101  return tmp<fvScalarMatrix>
102  (
103  new fvScalarMatrix
104  (
105  k_,
106  dimVolume*this->rho_.dimensions()*k_.dimensions()
107  /dimTime
108  )
109  );
110 }
111 
112 
113 template<class BasicTurbulenceModel>
115 {
116  return tmp<fvScalarMatrix>
117  (
118  new fvScalarMatrix
119  (
120  epsilon_,
121  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
122  /dimTime
123  )
124  );
125 }
126 
127 
128 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
129 
130 template<class BasicTurbulenceModel>
132 (
133  const alphaField& alpha,
134  const rhoField& rho,
135  const volVectorField& U,
136  const surfaceScalarField& alphaRhoPhi,
137  const surfaceScalarField& phi,
138  const transportModel& transport,
139  const word& propertiesName,
140  const word& type
141 )
142 :
144  (
145  type,
146  alpha,
147  rho,
148  U,
149  alphaRhoPhi,
150  phi,
151  transport,
152  propertiesName
153  ),
154  A0_
155  (
157  (
158  "A0",
159  this->coeffDict_,
160  4.0
161  )
162  ),
163  C2_
164  (
166  (
167  "C2",
168  this->coeffDict_,
169  1.9
170  )
171  ),
172  sigmak_
173  (
175  (
176  "sigmak",
177  this->coeffDict_,
178  1.0
179  )
180  ),
181  sigmaEps_
182  (
184  (
185  "sigmaEps",
186  this->coeffDict_,
187  1.2
188  )
189  ),
190 
191  k_
192  (
193  IOobject
194  (
195  "k",
196  this->runTime_.timeName(),
197  this->mesh_,
200  ),
201  this->mesh_
202  ),
203  epsilon_
204  (
205  IOobject
206  (
207  "epsilon",
208  this->runTime_.timeName(),
209  this->mesh_,
212  ),
213  this->mesh_
214  )
215 {
216  bound(k_, this->kMin_);
217  bound(epsilon_, this->epsilonMin_);
218 
219  if (type == typeName)
220  {
221  this->printCoeffs(type);
222  }
223 }
224 
225 
226 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
227 
228 template<class BasicTurbulenceModel>
230 {
232  {
233  A0_.readIfPresent(this->coeffDict());
234  C2_.readIfPresent(this->coeffDict());
235  sigmak_.readIfPresent(this->coeffDict());
236  sigmaEps_.readIfPresent(this->coeffDict());
237 
238  return true;
239  }
240  else
241  {
242  return false;
243  }
244 }
245 
246 
247 template<class BasicTurbulenceModel>
249 {
250  if (!this->turbulence_)
251  {
252  return;
253  }
254 
255  // Local references
256  const alphaField& alpha = this->alpha_;
257  const rhoField& rho = this->rho_;
258  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
259  const volVectorField& U = this->U_;
260  volScalarField& nut = this->nut_;
261  fv::options& fvOptions(fv::options::New(this->mesh_));
262 
264 
266 
267  tmp<volTensorField> tgradU = fvc::grad(U);
268  volScalarField S2(2*magSqr(dev(symm(tgradU()))));
269  volScalarField magS(sqrt(S2));
270 
271  volScalarField eta(magS*k_/epsilon_);
272  volScalarField C1(max(eta/(scalar(5) + eta), scalar(0.43)));
273 
274  volScalarField G(this->GName(), nut*(tgradU() && dev(twoSymm(tgradU()))));
275 
276  // Update epsilon and G at the wall
277  epsilon_.boundaryFieldRef().updateCoeffs();
278 
279  // Dissipation equation
280  tmp<fvScalarMatrix> epsEqn
281  (
282  fvm::ddt(alpha, rho, epsilon_)
283  + fvm::div(alphaRhoPhi, epsilon_)
284  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
285  ==
286  C1*alpha*rho*magS*epsilon_
287  - fvm::Sp
288  (
289  C2_*alpha*rho*epsilon_/(k_ + sqrt(this->nu()*epsilon_)),
290  epsilon_
291  )
292  + epsilonSource()
293  + fvOptions(alpha, rho, epsilon_)
294  );
295 
296  epsEqn.ref().relax();
297  fvOptions.constrain(epsEqn.ref());
298  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
299  solve(epsEqn);
300  fvOptions.correct(epsilon_);
301  bound(epsilon_, this->epsilonMin_);
302 
303 
304  // Turbulent kinetic energy equation
305 
307  (
308  fvm::ddt(alpha, rho, k_)
309  + fvm::div(alphaRhoPhi, k_)
310  - fvm::laplacian(alpha*rho*DkEff(), k_)
311  ==
312  alpha*rho*G
313  - fvm::SuSp(2.0/3.0*alpha*rho*divU, k_)
314  - fvm::Sp(alpha*rho*epsilon_/k_, k_)
315  + kSource()
316  + fvOptions(alpha, rho, k_)
317  );
318 
319  kEqn.ref().relax();
320  fvOptions.constrain(kEqn.ref());
321  solve(kEqn);
322  fvOptions.correct(k_);
323  bound(k_, this->kMin_);
324 
325  correctNut(tgradU(), S2, magS);
326 }
327 
328 
329 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
330 
331 } // End namespace RASModels
332 } // End namespace Foam
333 
334 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
dimensionedScalar acos(const dimensionedScalar &ds)
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
surfaceScalarField & phi
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedTensor skew(const dimensionedTensor &dt)
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
const dimensionedScalar G
Newtonian constant of gravitation.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:100
dimensionedScalar sqrt(const dimensionedScalar &ds)
Generic dimensioned Type class.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
virtual tmp< fvScalarMatrix > kSource() const
Definition: realizableKE.C:99
Finite-volume options.
Definition: fvOptions.H:52
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:101
virtual bool read()
Re-read model coefficients if they have changed.
Definition: realizableKE.C:229
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
tmp< volScalarField > rCmu(const volTensorField &gradU, const volScalarField &S2, const volScalarField &magS)
Definition: realizableKE.C:41
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Dimension set for the base types.
Definition: dimensionSet.H:120
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar cos(const dimensionedScalar &ds)
A class for handling words, derived from string.
Definition: word.H:59
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Bound the given scalar field if it has gone unbounded.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/ubuntu/OpenFOAM-6/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:72
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: realizableKE.C:114
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
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
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: realizableKE.C:248
const scalarList W(::W(thermo))
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
realizableKE(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
Definition: realizableKE.C:132
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
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
volScalarField & nu
scalar nut
Namespace for OpenFOAM.
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:99