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