realizableKE.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 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 // * * * * * * * * * * * * Private 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 
155  Cmu_
156  (
158  (
159  "Cmu",
160  this->coeffDict_,
161  0.09
162  )
163  ),
164  A0_
165  (
167  (
168  "A0",
169  this->coeffDict_,
170  4.0
171  )
172  ),
173  C2_
174  (
176  (
177  "C2",
178  this->coeffDict_,
179  1.9
180  )
181  ),
182  sigmak_
183  (
185  (
186  "sigmak",
187  this->coeffDict_,
188  1.0
189  )
190  ),
191  sigmaEps_
192  (
194  (
195  "sigmaEps",
196  this->coeffDict_,
197  1.2
198  )
199  ),
200 
201  k_
202  (
203  IOobject
204  (
205  "k",
206  this->runTime_.timeName(),
207  this->mesh_,
210  ),
211  this->mesh_
212  ),
213  epsilon_
214  (
215  IOobject
216  (
217  "epsilon",
218  this->runTime_.timeName(),
219  this->mesh_,
222  ),
223  this->mesh_
224  )
225 {
226  bound(k_, this->kMin_);
227  bound(epsilon_, this->epsilonMin_);
228 
229  if (type == typeName)
230  {
231  this->printCoeffs(type);
232  }
233 }
234 
235 
236 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
237 
238 template<class BasicTurbulenceModel>
240 {
242  {
243  Cmu_.readIfPresent(this->coeffDict());
244  A0_.readIfPresent(this->coeffDict());
245  C2_.readIfPresent(this->coeffDict());
246  sigmak_.readIfPresent(this->coeffDict());
247  sigmaEps_.readIfPresent(this->coeffDict());
248 
249  return true;
250  }
251  else
252  {
253  return false;
254  }
255 }
256 
257 
258 template<class BasicTurbulenceModel>
260 {
261  if (!this->turbulence_)
262  {
263  return;
264  }
265 
266  // Local references
267  const alphaField& alpha = this->alpha_;
268  const rhoField& rho = this->rho_;
269  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
270  const volVectorField& U = this->U_;
271  volScalarField& nut = this->nut_;
272  fv::options& fvOptions(fv::options::New(this->mesh_));
273 
275 
277 
278  tmp<volTensorField> tgradU = fvc::grad(U);
279  volScalarField S2(2*magSqr(dev(symm(tgradU()))));
280  volScalarField magS(sqrt(S2));
281 
282  volScalarField eta(magS*k_/epsilon_);
283  volScalarField C1(max(eta/(scalar(5) + eta), scalar(0.43)));
284 
285  volScalarField G(this->GName(), nut*(tgradU() && dev(twoSymm(tgradU()))));
286 
287  // Update epsilon and G at the wall
288  epsilon_.boundaryFieldRef().updateCoeffs();
289 
290  // Dissipation equation
291  tmp<fvScalarMatrix> epsEqn
292  (
293  fvm::ddt(alpha, rho, epsilon_)
294  + fvm::div(alphaRhoPhi, epsilon_)
295  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
296  ==
297  C1*alpha*rho*magS*epsilon_
298  - fvm::Sp
299  (
300  C2_*alpha*rho*epsilon_/(k_ + sqrt(this->nu()*epsilon_)),
301  epsilon_
302  )
303  + epsilonSource()
304  + fvOptions(alpha, rho, epsilon_)
305  );
306 
307  epsEqn.ref().relax();
308  fvOptions.constrain(epsEqn.ref());
309  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
310  solve(epsEqn);
311  fvOptions.correct(epsilon_);
312  bound(epsilon_, this->epsilonMin_);
313 
314 
315  // Turbulent kinetic energy equation
316 
318  (
319  fvm::ddt(alpha, rho, k_)
320  + fvm::div(alphaRhoPhi, k_)
321  - fvm::laplacian(alpha*rho*DkEff(), k_)
322  ==
323  alpha*rho*G
324  - fvm::SuSp(2.0/3.0*alpha*rho*divU, k_)
325  - fvm::Sp(alpha*rho*epsilon_/k_, k_)
326  + kSource()
327  + fvOptions(alpha, rho, k_)
328  );
329 
330  kEqn.ref().relax();
331  fvOptions.constrain(kEqn.ref());
332  solve(kEqn);
333  fvOptions.correct(k_);
334  bound(k_, this->kMin_);
335 
336  correctNut(tgradU(), S2, magS);
337 }
338 
339 
340 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
341 
342 } // End namespace RASModels
343 } // End namespace Foam
344 
345 // ************************************************************************* //
surfaceScalarField & phi
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: realizableKE.C:114
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
U
Definition: pEqn.H:83
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedTensor skew(const dimensionedTensor &dt)
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
const dimensionedScalar G
Newtonian constant of gravitation.
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.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
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:239
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
fv::options & fvOptions
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
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:230
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Dimension set for the base types.
Definition: dimensionSet.H:118
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar cos(const dimensionedScalar &ds)
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
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
virtual tmp< fvScalarMatrix > kSource() const
Definition: realizableKE.C:99
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-4.1/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 > &)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:188
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
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:259
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
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
A class for managing temporary objects.
Definition: PtrList.H:54
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:103
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
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:91
volScalarField & nu
scalar nut
Namespace for OpenFOAM.
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:99