v2f.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) 2012-2020 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 BasicMomentumTransportModel>
41 {
42  return max(k_/epsilon_, 6.0*sqrt(this->nu()/epsilon_));
43 }
44 
45 
46 template<class BasicMomentumTransportModel>
48 {
49  return
50  CL_*max(pow(k_, 1.5)
51  /epsilon_, Ceta_*pow025(pow3(this->nu())/epsilon_));
52 }
53 
54 
55 template<class BasicMomentumTransportModel>
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 
63 
64 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
65 
66 template<class BasicMomentumTransportModel>
68 (
69  const alphaField& alpha,
70  const rhoField& rho,
71  const volVectorField& U,
72  const surfaceScalarField& alphaRhoPhi,
73  const surfaceScalarField& phi,
74  const transportModel& transport,
75  const word& type
76 )
77 :
79  (
80  type,
81  alpha,
82  rho,
83  U,
84  alphaRhoPhi,
85  phi,
86  transport
87  ),
88  v2fBase(),
89 
90  Cmu_
91  (
93  (
94  "Cmu",
95  this->coeffDict_,
96  0.22
97  )
98  ),
99  CmuKEps_
100  (
102  (
103  "CmuKEps",
104  this->coeffDict_,
105  0.09
106  )
107  ),
108  C1_
109  (
111  (
112  "C1",
113  this->coeffDict_,
114  1.4
115  )
116  ),
117  C2_
118  (
120  (
121  "C2",
122  this->coeffDict_,
123  0.3
124  )
125  ),
126  CL_
127  (
129  (
130  "CL",
131  this->coeffDict_,
132  0.23
133  )
134  ),
135  Ceta_
136  (
138  (
139  "Ceta",
140  this->coeffDict_,
141  70.0
142  )
143  ),
144  Ceps2_
145  (
147  (
148  "Ceps2",
149  this->coeffDict_,
150  1.9
151  )
152  ),
153  Ceps3_
154  (
156  (
157  "Ceps3",
158  this->coeffDict_,
159  -0.33
160  )
161  ),
162  sigmaK_
163  (
165  (
166  "sigmaK",
167  this->coeffDict_,
168  1.0
169  )
170  ),
171  sigmaEps_
172  (
174  (
175  "sigmaEps",
176  this->coeffDict_,
177  1.3
178  )
179  ),
180 
181  k_
182  (
183  IOobject
184  (
185  IOobject::groupName("k", alphaRhoPhi.group()),
186  this->runTime_.timeName(),
187  this->mesh_,
190  ),
191  this->mesh_
192  ),
193  epsilon_
194  (
195  IOobject
196  (
197  IOobject::groupName("epsilon", alphaRhoPhi.group()),
198  this->runTime_.timeName(),
199  this->mesh_,
202  ),
203  this->mesh_
204  ),
205  v2_
206  (
207  IOobject
208  (
209  IOobject::groupName("v2", alphaRhoPhi.group()),
210  this->runTime_.timeName(),
211  this->mesh_,
214  ),
215  this->mesh_
216  ),
217  f_
218  (
219  IOobject
220  (
221  IOobject::groupName("f", alphaRhoPhi.group()),
222  this->runTime_.timeName(),
223  this->mesh_,
226  ),
227  this->mesh_
228  ),
229  v2Min_(dimensionedScalar(v2_.dimensions(), small)),
230  fMin_(dimensionedScalar(f_.dimensions(), 0))
231 {
232  bound(k_, this->kMin_);
233  bound(epsilon_, this->epsilonMin_);
234  bound(v2_, v2Min_);
235  bound(f_, fMin_);
236 
237  if (type == typeName)
238  {
239  this->printCoeffs(type);
240  }
241 }
242 
243 
244 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
245 
246 template<class BasicMomentumTransportModel>
248 {
250  {
251  Cmu_.readIfPresent(this->coeffDict());
252  CmuKEps_.readIfPresent(this->coeffDict());
253  C1_.readIfPresent(this->coeffDict());
254  C2_.readIfPresent(this->coeffDict());
255  CL_.readIfPresent(this->coeffDict());
256  Ceta_.readIfPresent(this->coeffDict());
257  Ceps2_.readIfPresent(this->coeffDict());
258  Ceps3_.readIfPresent(this->coeffDict());
259  sigmaK_.readIfPresent(this->coeffDict());
260  sigmaEps_.readIfPresent(this->coeffDict());
261 
262  return true;
263  }
264  else
265  {
266  return false;
267  }
268 }
269 
270 
271 template<class BasicMomentumTransportModel>
273 {
274  if (!this->turbulence_)
275  {
276  return;
277  }
278 
279  // Local references
280  const alphaField& alpha = this->alpha_;
281  const rhoField& rho = this->rho_;
282  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
283  const volVectorField& U = this->U_;
284  volScalarField& nut = this->nut_;
285  fv::options& fvOptions(fv::options::New(this->mesh_));
286 
288 
290 
291  // Use N=6 so that f=0 at walls
292  const dimensionedScalar N("N", dimless, 6.0);
293 
294  const volTensorField gradU(fvc::grad(U));
295  const volScalarField S2(2*magSqr(dev(symm(gradU))));
296 
297  const volScalarField G(this->GName(), nut*S2);
298  const volScalarField Ts(this->Ts());
299  const volScalarField L2(type() + ":L2", sqr(Ls()));
300  const volScalarField v2fAlpha
301  (
302  type() + ":alpha",
303  1.0/Ts*((C1_ - N)*v2_ - 2.0/3.0*k_*(C1_ - 1.0))
304  );
305 
306  const volScalarField Ceps1
307  (
308  "Ceps1",
309  1.4*(1.0 + 0.05*min(sqrt(k_/v2_), scalar(100.0)))
310  );
311 
312  // Update epsilon (and possibly G) at the wall
313  epsilon_.boundaryFieldRef().updateCoeffs();
314 
315  // Dissipation equation
316  tmp<fvScalarMatrix> epsEqn
317  (
318  fvm::ddt(alpha, rho, epsilon_)
319  + fvm::div(alphaRhoPhi, epsilon_)
320  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
321  ==
322  Ceps1*alpha*rho*G/Ts
323  - fvm::SuSp(((2.0/3.0)*Ceps1 + Ceps3_)*alpha*rho*divU, epsilon_)
324  - fvm::Sp(Ceps2_*alpha*rho/Ts, epsilon_)
325  + fvOptions(alpha, rho, epsilon_)
326  );
327 
328  epsEqn.ref().relax();
329  fvOptions.constrain(epsEqn.ref());
330  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
331  solve(epsEqn);
332  fvOptions.correct(epsilon_);
333  bound(epsilon_, this->epsilonMin_);
334 
335 
336  // Turbulent kinetic energy equation
338  (
339  fvm::ddt(alpha, rho, k_)
340  + fvm::div(alphaRhoPhi, k_)
341  - fvm::laplacian(alpha*rho*DkEff(), k_)
342  ==
343  alpha*rho*G
344  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
345  - fvm::Sp(alpha*rho*epsilon_/k_, k_)
346  + fvOptions(alpha, rho, k_)
347  );
348 
349  kEqn.ref().relax();
350  fvOptions.constrain(kEqn.ref());
351  solve(kEqn);
352  fvOptions.correct(k_);
353  bound(k_, this->kMin_);
354 
355 
356  // Relaxation function equation
358  (
359  - fvm::laplacian(f_)
360  ==
361  - fvm::Sp(1.0/L2, f_)
362  - 1.0/L2/k_*(v2fAlpha - C2_*G)
363  );
364 
365  fEqn.ref().relax();
366  fvOptions.constrain(fEqn.ref());
367  solve(fEqn);
368  fvOptions.correct(f_);
369  bound(f_, fMin_);
370 
371 
372  // Turbulence stress normal to streamlines equation
373  tmp<fvScalarMatrix> v2Eqn
374  (
375  fvm::ddt(alpha, rho, v2_)
376  + fvm::div(alphaRhoPhi, v2_)
377  - fvm::laplacian(alpha*rho*DkEff(), v2_)
378  ==
379  alpha*rho*min(k_*f_, C2_*G - v2fAlpha)
380  - fvm::Sp(N*alpha*rho*epsilon_/k_, v2_)
381  + fvOptions(alpha, rho, v2_)
382  );
383 
384  v2Eqn.ref().relax();
385  fvOptions.constrain(v2Eqn.ref());
386  solve(v2Eqn);
387  fvOptions.correct(v2_);
388  bound(v2_, v2Min_);
389 
390  correctNut();
391 }
392 
393 
394 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
395 
396 } // End namespace RASModels
397 } // End namespace Foam
398 
399 // ************************************************************************* //
v2f(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
Definition: v2f.C:68
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
fv::options & fvOptions
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
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
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:49
Finite-volume options.
Definition: fvOptions.H:52
BasicMomentumTransportModel::alphaField alphaField
Definition: RASModel.H:88
phi
Definition: pEqn.H:104
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:54
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
BasicMomentumTransportModel::transportModel transportModel
Definition: RASModel.H:90
A class for handling words, derived from string.
Definition: word.H:59
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:272
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(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-8/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:68
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
BasicMomentumTransportModel::rhoField rhoField
Definition: RASModel.H:89
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
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:247
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.
label N
Definition: createFields.H:22
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
const dimensionedScalar & G
Newtonian constant of gravitation.
scalar nut
Namespace for OpenFOAM.
tmp< volScalarField > Ts() const
Return time scale, Ts.
Definition: v2f.C:40