v2f.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) 2012-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 "v2f.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  return max(k_/epsilon_, 6.0*sqrt(this->nu()/epsilon_));
43 }
44 
45 
46 template<class BasicTurbulenceModel>
48 {
49  return
50  CL_*max(pow(k_, 1.5)
51  /epsilon_, Ceta_*pow025(pow3(this->nu())/epsilon_));
52 }
53 
54 
55 template<class BasicTurbulenceModel>
57 {
58  this->nut_ = min(CmuKEps_*sqr(k_)/epsilon_, this->Cmu_*v2_*Ts());
59  this->nut_.correctBoundaryConditions();
60  fv::options::New(this->mesh_).correct(this->nut_);
61 
62  BasicTurbulenceModel::correctNut();
63 }
64 
65 
66 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
67 
68 template<class BasicTurbulenceModel>
70 (
71  const alphaField& alpha,
72  const rhoField& rho,
73  const volVectorField& U,
74  const surfaceScalarField& alphaRhoPhi,
75  const surfaceScalarField& phi,
76  const transportModel& transport,
77  const word& propertiesName,
78  const word& type
79 )
80 :
82  (
83  type,
84  alpha,
85  rho,
86  U,
87  alphaRhoPhi,
88  phi,
89  transport,
90  propertiesName
91  ),
92  v2fBase(),
93 
94  Cmu_
95  (
97  (
98  "Cmu",
99  this->coeffDict_,
100  0.22
101  )
102  ),
103  CmuKEps_
104  (
106  (
107  "CmuKEps",
108  this->coeffDict_,
109  0.09
110  )
111  ),
112  C1_
113  (
115  (
116  "C1",
117  this->coeffDict_,
118  1.4
119  )
120  ),
121  C2_
122  (
124  (
125  "C2",
126  this->coeffDict_,
127  0.3
128  )
129  ),
130  CL_
131  (
133  (
134  "CL",
135  this->coeffDict_,
136  0.23
137  )
138  ),
139  Ceta_
140  (
142  (
143  "Ceta",
144  this->coeffDict_,
145  70.0
146  )
147  ),
148  Ceps2_
149  (
151  (
152  "Ceps2",
153  this->coeffDict_,
154  1.9
155  )
156  ),
157  Ceps3_
158  (
160  (
161  "Ceps3",
162  this->coeffDict_,
163  -0.33
164  )
165  ),
166  sigmaK_
167  (
169  (
170  "sigmaK",
171  this->coeffDict_,
172  1.0
173  )
174  ),
175  sigmaEps_
176  (
178  (
179  "sigmaEps",
180  this->coeffDict_,
181  1.3
182  )
183  ),
184 
185  k_
186  (
187  IOobject
188  (
189  IOobject::groupName("k", U.group()),
190  this->runTime_.timeName(),
191  this->mesh_,
194  ),
195  this->mesh_
196  ),
197  epsilon_
198  (
199  IOobject
200  (
201  IOobject::groupName("epsilon", U.group()),
202  this->runTime_.timeName(),
203  this->mesh_,
206  ),
207  this->mesh_
208  ),
209  v2_
210  (
211  IOobject
212  (
213  IOobject::groupName("v2", U.group()),
214  this->runTime_.timeName(),
215  this->mesh_,
218  ),
219  this->mesh_
220  ),
221  f_
222  (
223  IOobject
224  (
225  IOobject::groupName("f", U.group()),
226  this->runTime_.timeName(),
227  this->mesh_,
230  ),
231  this->mesh_
232  ),
233  v2Min_(dimensionedScalar("v2Min", v2_.dimensions(), SMALL)),
234  fMin_(dimensionedScalar("fMin", f_.dimensions(), 0.0))
235 {
236  bound(k_, this->kMin_);
237  bound(epsilon_, this->epsilonMin_);
238  bound(v2_, v2Min_);
239  bound(f_, fMin_);
240 
241  if (type == typeName)
242  {
243  this->printCoeffs(type);
244  }
245 }
246 
247 
248 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
249 
250 template<class BasicTurbulenceModel>
252 {
254  {
255  Cmu_.readIfPresent(this->coeffDict());
256  CmuKEps_.readIfPresent(this->coeffDict());
257  C1_.readIfPresent(this->coeffDict());
258  C2_.readIfPresent(this->coeffDict());
259  CL_.readIfPresent(this->coeffDict());
260  Ceta_.readIfPresent(this->coeffDict());
261  Ceps2_.readIfPresent(this->coeffDict());
262  Ceps3_.readIfPresent(this->coeffDict());
263  sigmaK_.readIfPresent(this->coeffDict());
264  sigmaEps_.readIfPresent(this->coeffDict());
265 
266  return true;
267  }
268  else
269  {
270  return false;
271  }
272 }
273 
274 
275 template<class BasicTurbulenceModel>
277 {
278  if (!this->turbulence_)
279  {
280  return;
281  }
282 
283  // Local references
284  const alphaField& alpha = this->alpha_;
285  const rhoField& rho = this->rho_;
286  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
287  const volVectorField& U = this->U_;
288  volScalarField& nut = this->nut_;
289  fv::options& fvOptions(fv::options::New(this->mesh_));
290 
292 
294 
295  // Use N=6 so that f=0 at walls
296  const dimensionedScalar N("N", dimless, 6.0);
297 
298  const volTensorField gradU(fvc::grad(U));
299  const volScalarField S2(2*magSqr(dev(symm(gradU))));
300 
301  const volScalarField G(this->GName(), nut*S2);
302  const volScalarField Ts(this->Ts());
303  const volScalarField L2(type() + ":L2", sqr(Ls()));
304  const volScalarField v2fAlpha
305  (
306  type() + ":alpha",
307  1.0/Ts*((C1_ - N)*v2_ - 2.0/3.0*k_*(C1_ - 1.0))
308  );
309 
310  const volScalarField Ceps1
311  (
312  "Ceps1",
313  1.4*(1.0 + 0.05*min(sqrt(k_/v2_), scalar(100.0)))
314  );
315 
316  // Update epsilon (and possibly G) at the wall
317  epsilon_.boundaryFieldRef().updateCoeffs();
318 
319  // Dissipation equation
320  tmp<fvScalarMatrix> epsEqn
321  (
322  fvm::ddt(alpha, rho, epsilon_)
323  + fvm::div(alphaRhoPhi, epsilon_)
324  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
325  ==
326  Ceps1*alpha*rho*G/Ts
327  - fvm::SuSp(((2.0/3.0)*Ceps1 + Ceps3_)*alpha*rho*divU, epsilon_)
328  - fvm::Sp(Ceps2_*alpha*rho/Ts, epsilon_)
329  + fvOptions(alpha, rho, epsilon_)
330  );
331 
332  epsEqn.ref().relax();
333  fvOptions.constrain(epsEqn.ref());
334  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
335  solve(epsEqn);
336  fvOptions.correct(epsilon_);
337  bound(epsilon_, this->epsilonMin_);
338 
339 
340  // Turbulent kinetic energy equation
342  (
343  fvm::ddt(alpha, rho, k_)
344  + fvm::div(alphaRhoPhi, k_)
345  - fvm::laplacian(alpha*rho*DkEff(), k_)
346  ==
347  alpha*rho*G
348  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
349  - fvm::Sp(alpha*rho*epsilon_/k_, k_)
350  + fvOptions(alpha, rho, k_)
351  );
352 
353  kEqn.ref().relax();
354  fvOptions.constrain(kEqn.ref());
355  solve(kEqn);
356  fvOptions.correct(k_);
357  bound(k_, this->kMin_);
358 
359 
360  // Relaxation function equation
362  (
363  - fvm::laplacian(f_)
364  ==
365  - fvm::Sp(1.0/L2, f_)
366  - 1.0/L2/k_*(v2fAlpha - C2_*G)
367  );
368 
369  fEqn.ref().relax();
370  fvOptions.constrain(fEqn.ref());
371  solve(fEqn);
372  fvOptions.correct(f_);
373  bound(f_, fMin_);
374 
375 
376  // Turbulence stress normal to streamlines equation
377  tmp<fvScalarMatrix> v2Eqn
378  (
379  fvm::ddt(alpha, rho, v2_)
380  + fvm::div(alphaRhoPhi, v2_)
381  - fvm::laplacian(alpha*rho*DkEff(), v2_)
382  ==
383  alpha*rho*min(k_*f_, C2_*G - v2fAlpha)
384  - fvm::Sp(N*alpha*rho*epsilon_/k_, v2_)
385  + fvOptions(alpha, rho, v2_)
386  );
387 
388  v2Eqn.ref().relax();
389  fvOptions.constrain(v2Eqn.ref());
390  solve(v2Eqn);
391  fvOptions.correct(v2_);
392  bound(v2_, v2Min_);
393 
394  correctNut();
395 }
396 
397 
398 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
399 
400 } // End namespace RASModels
401 } // End namespace Foam
402 
403 // ************************************************************************* //
surfaceScalarField & phi
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
U
Definition: pEqn.H:83
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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)
dimensionedScalar pow025(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:52
Finite-volume options.
Definition: fvOptions.H:52
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:101
fv::options & fvOptions
tmp< volScalarField > Ls() const
Return length scale, Ls.
Definition: v2f.C:47
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
Abstract base-class for v2-f models to provide BCs access to the v2 and f fields. ...
Definition: v2fBase.H:57
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
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.
static word groupName(Name name, const word &group)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: v2f.C:276
v2f(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: v2f.C:70
virtual void correctNut()
Definition: v2f.C:56
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-5.0/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 > &)
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:485
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:188
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
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 bool read()
Read RASProperties dictionary.
Definition: v2f.C:251
label N
Definition: createFields.H:22
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:103
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.
tmp< volScalarField > Ts() const
Return time scale, Ts.
Definition: v2f.C:40
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:99