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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2012-2017 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "v2f.H"
30 #include "fvOptions.H"
31 #include "bound.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace RASModels
38 {
39 
40 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41 
42 template<class BasicTurbulenceModel>
44 {
45  // SAF: limiting thermo->nu(). If psiThermo is used rho might be < 0
46  // temporarily when p < 0 then nu < 0 which needs limiting
47  return
48  max
49  (
50  k_/epsilon_,
51  6.0*sqrt
52  (
53  max
54  (
55  this->nu(),
56  dimensionedScalar(this->nu()().dimensions(), Zero)
57  )
58  / epsilon_
59  )
60  );
61 }
62 
63 
64 template<class BasicTurbulenceModel>
66 {
67  // SAF: limiting thermo->nu(). If psiThermo is used rho might be < 0
68  // temporarily when p < 0 then nu < 0 which needs limiting
69  return
70  CL_
71  * max
72  (
73  pow(k_, 1.5)/epsilon_,
74  Ceta_*pow025
75  (
76  pow3
77  (
78  max
79  (
80  this->nu(),
81  dimensionedScalar(this->nu()().dimensions(), Zero)
82  )
83  )/epsilon_
84  )
85  );
86 }
87 
88 
89 template<class BasicTurbulenceModel>
91 {
92  this->nut_ = min(CmuKEps_*sqr(k_)/epsilon_, this->Cmu_*v2_*Ts());
93  this->nut_.correctBoundaryConditions();
94  fv::options::New(this->mesh_).correct(this->nut_);
95 
96  BasicTurbulenceModel::correctNut();
97 }
98 
99 
100 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
101 
102 template<class BasicTurbulenceModel>
104 (
105  const alphaField& alpha,
106  const rhoField& rho,
107  const volVectorField& U,
108  const surfaceScalarField& alphaRhoPhi,
109  const surfaceScalarField& phi,
110  const transportModel& transport,
111  const word& propertiesName,
112  const word& type
113 )
114 :
116  (
117  type,
118  alpha,
119  rho,
120  U,
121  alphaRhoPhi,
122  phi,
123  transport,
124  propertiesName
125  ),
126  v2fBase(),
127 
128  Cmu_
129  (
131  (
132  "Cmu",
133  this->coeffDict_,
134  0.22
135  )
136  ),
137  CmuKEps_
138  (
140  (
141  "CmuKEps",
142  this->coeffDict_,
143  0.09
144  )
145  ),
146  C1_
147  (
149  (
150  "C1",
151  this->coeffDict_,
152  1.4
153  )
154  ),
155  C2_
156  (
158  (
159  "C2",
160  this->coeffDict_,
161  0.3
162  )
163  ),
164  CL_
165  (
167  (
168  "CL",
169  this->coeffDict_,
170  0.23
171  )
172  ),
173  Ceta_
174  (
176  (
177  "Ceta",
178  this->coeffDict_,
179  70.0
180  )
181  ),
182  Ceps2_
183  (
185  (
186  "Ceps2",
187  this->coeffDict_,
188  1.9
189  )
190  ),
191  Ceps3_
192  (
194  (
195  "Ceps3",
196  this->coeffDict_,
197  -0.33
198  )
199  ),
200  sigmaK_
201  (
203  (
204  "sigmaK",
205  this->coeffDict_,
206  1.0
207  )
208  ),
209  sigmaEps_
210  (
212  (
213  "sigmaEps",
214  this->coeffDict_,
215  1.3
216  )
217  ),
218 
219  k_
220  (
221  IOobject
222  (
223  IOobject::groupName("k", alphaRhoPhi.group()),
224  this->runTime_.timeName(),
225  this->mesh_,
228  ),
229  this->mesh_
230  ),
231  epsilon_
232  (
233  IOobject
234  (
235  IOobject::groupName("epsilon", alphaRhoPhi.group()),
236  this->runTime_.timeName(),
237  this->mesh_,
240  ),
241  this->mesh_
242  ),
243  v2_
244  (
245  IOobject
246  (
247  IOobject::groupName("v2", alphaRhoPhi.group()),
248  this->runTime_.timeName(),
249  this->mesh_,
252  ),
253  this->mesh_
254  ),
255  f_
256  (
257  IOobject
258  (
259  IOobject::groupName("f", alphaRhoPhi.group()),
260  this->runTime_.timeName(),
261  this->mesh_,
264  ),
265  this->mesh_
266  ),
267  v2Min_(dimensionedScalar("v2Min", v2_.dimensions(), SMALL)),
268  fMin_(dimensionedScalar("fMin", f_.dimensions(), Zero))
269 {
270  bound(k_, this->kMin_);
271  bound(epsilon_, this->epsilonMin_);
272  bound(v2_, v2Min_);
273  bound(f_, fMin_);
274 
275  if (type == typeName)
276  {
277  this->printCoeffs(type);
278  }
279 }
280 
281 
282 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
283 
284 template<class BasicTurbulenceModel>
286 {
288  {
289  Cmu_.readIfPresent(this->coeffDict());
290  CmuKEps_.readIfPresent(this->coeffDict());
291  C1_.readIfPresent(this->coeffDict());
292  C2_.readIfPresent(this->coeffDict());
293  CL_.readIfPresent(this->coeffDict());
294  Ceta_.readIfPresent(this->coeffDict());
295  Ceps2_.readIfPresent(this->coeffDict());
296  Ceps3_.readIfPresent(this->coeffDict());
297  sigmaK_.readIfPresent(this->coeffDict());
298  sigmaEps_.readIfPresent(this->coeffDict());
299 
300  return true;
301  }
302 
303  return false;
304 }
305 
306 
307 template<class BasicTurbulenceModel>
309 {
310  if (!this->turbulence_)
311  {
312  return;
313  }
314 
315  // Local references
316  const alphaField& alpha = this->alpha_;
317  const rhoField& rho = this->rho_;
318  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
319  const volVectorField& U = this->U_;
320  volScalarField& nut = this->nut_;
321  fv::options& fvOptions(fv::options::New(this->mesh_));
322 
324 
326 
327  // Use N=6 so that f=0 at walls
328  const dimensionedScalar N("N", dimless, 6.0);
329 
330  const volTensorField gradU(fvc::grad(U));
331  const volScalarField S2(2*magSqr(dev(symm(gradU))));
332 
333  const volScalarField G(this->GName(), nut*S2);
334  const volScalarField Ts(this->Ts());
335  const volScalarField L2(type() + ":L2", sqr(Ls()));
336  const volScalarField v2fAlpha
337  (
338  type() + ":alpha",
339  1.0/Ts*((C1_ - N)*v2_ - 2.0/3.0*k_*(C1_ - 1.0))
340  );
341 
342  const volScalarField Ceps1
343  (
344  "Ceps1",
345  1.4*(1.0 + 0.05*min(sqrt(k_/v2_), scalar(100)))
346  );
347 
348  // Update epsilon (and possibly G) at the wall
349  epsilon_.boundaryFieldRef().updateCoeffs();
350 
351  // Dissipation equation
352  tmp<fvScalarMatrix> epsEqn
353  (
354  fvm::ddt(alpha, rho, epsilon_)
355  + fvm::div(alphaRhoPhi, epsilon_)
356  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
357  ==
358  Ceps1*alpha*rho*G/Ts
359  - fvm::SuSp(((2.0/3.0)*Ceps1 + Ceps3_)*alpha*rho*divU, epsilon_)
360  - fvm::Sp(Ceps2_*alpha*rho/Ts, epsilon_)
361  + fvOptions(alpha, rho, epsilon_)
362  );
363 
364  epsEqn.ref().relax();
365  fvOptions.constrain(epsEqn.ref());
366  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
367  solve(epsEqn);
368  fvOptions.correct(epsilon_);
369  bound(epsilon_, this->epsilonMin_);
370 
371 
372  // Turbulent kinetic energy equation
374  (
375  fvm::ddt(alpha, rho, k_)
376  + fvm::div(alphaRhoPhi, k_)
377  - fvm::laplacian(alpha*rho*DkEff(), k_)
378  ==
379  alpha*rho*G
380  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
381  - fvm::Sp(alpha*rho*epsilon_/k_, k_)
382  + fvOptions(alpha, rho, k_)
383  );
384 
385  kEqn.ref().relax();
386  fvOptions.constrain(kEqn.ref());
387  solve(kEqn);
388  fvOptions.correct(k_);
389  bound(k_, this->kMin_);
390 
391 
392  // Relaxation function equation
394  (
395  - fvm::laplacian(f_)
396  ==
397  - fvm::Sp(1.0/L2, f_)
398  - 1.0/L2/k_*(v2fAlpha - C2_*G)
399  );
400 
401  fEqn.ref().relax();
402  fvOptions.constrain(fEqn.ref());
403  solve(fEqn);
404  fvOptions.correct(f_);
405  bound(f_, fMin_);
406 
407 
408  // Turbulence stress normal to streamlines equation
409  tmp<fvScalarMatrix> v2Eqn
410  (
411  fvm::ddt(alpha, rho, v2_)
412  + fvm::div(alphaRhoPhi, v2_)
413  - fvm::laplacian(alpha*rho*DkEff(), v2_)
414  ==
415  alpha*rho*min(k_*f_, C2_*G - v2fAlpha)
416  - fvm::Sp(N*alpha*rho*epsilon_/k_, v2_)
417  + fvOptions(alpha, rho, v2_)
418  );
419 
420  v2Eqn.ref().relax();
421  fvOptions.constrain(v2Eqn.ref());
422  solve(v2Eqn);
423  fvOptions.correct(v2_);
424  bound(v2_, v2Min_);
425 
426  correctNut();
427 }
428 
429 
430 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
431 
432 } // End namespace RASModels
433 } // End namespace Foam
434 
435 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:84
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
fvOptions.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fv::optionList::correct
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
Definition: fvOptionListTemplates.C:321
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp< volScalarField >
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::RASModel
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:51
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:104
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:35
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
rho
rho
Definition: readInitialConditions.H:96
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::RASModels::v2f::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: v2f.C:285
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::RASModels::v2f::Ts
tmp< volScalarField > Ts() const
Return time scale, Ts.
Definition: v2f.C:43
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:133
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::fvm::SuSp
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
bound.H
Bound the given scalar field if it has gone unbounded.
correct
fvOptions correct(rho)
Foam::solve
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
divU
zeroField divU
Definition: alphaSuSp.H:3
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::RASModels::v2f::v2f
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:104
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::fv::options
Finite-volume options.
Definition: fvOptions.H:55
Foam::RASModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:97
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:43
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::RASModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:96
Foam::RASModels::v2f::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: v2f.C:308
Foam::RASModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:98
U
U
Definition: pEqn.H:72
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::IOobject::groupName
static word groupName(StringType name, const word &group)
Create dot-delimited name.group.
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::RASModels::v2f::correctNut
virtual void correctNut()
Definition: v2f.C:90
v2f.H
Foam::eddyViscosity
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:55
N
const Vector< label > N(dict.get< Vector< label >>("N"))
nut
scalar nut
Definition: evaluateNearWall.H:5
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
Foam::RASModels::v2fBase
Abstract base-class for v2-f models to provide BCs access to the v2 and f fields.
Definition: v2fBase.H:60
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:190
Foam::RASModels::v2f::Ls
tmp< volScalarField > Ls() const
Return length scale, Ls.
Definition: v2f.C:65
Foam::IOobject::MUST_READ
Definition: IOobject.H:120
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106