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-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 "v2f.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  // SAF: limiting thermo->nu(). If psiThermo is used rho might be < 0
42  // temporarily when p < 0 then nu < 0 which needs limiting
43  return
44  max
45  (
46  k_/epsilon_,
47  6.0*sqrt
48  (
49  max
50  (
51  this->nu(),
52  dimensionedScalar("zero", this->nu()().dimensions(), 0.0)
53  )
54  / epsilon_
55  )
56  );
57 }
58 
59 
60 template<class BasicTurbulenceModel>
62 {
63  // SAF: limiting thermo->nu(). If psiThermo is used rho might be < 0
64  // temporarily when p < 0 then nu < 0 which needs limiting
65  return
66  CL_
67  * max
68  (
69  pow(k_, 1.5)/epsilon_,
70  Ceta_*pow025
71  (
72  pow3
73  (
74  max
75  (
76  this->nu(),
77  dimensionedScalar("zero", this->nu()().dimensions(), 0.0)
78  )
79  )/epsilon_
80  )
81  );
82 }
83 
84 
85 template<class BasicTurbulenceModel>
87 {
88  this->nut_ = min(CmuKEps_*sqr(k_)/epsilon_, this->Cmu_*v2_*Ts());
89  this->nut_.correctBoundaryConditions();
90 
91  BasicTurbulenceModel::correctNut();
92 }
93 
94 
95 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
96 
97 template<class BasicTurbulenceModel>
99 (
100  const alphaField& alpha,
101  const rhoField& rho,
102  const volVectorField& U,
103  const surfaceScalarField& alphaRhoPhi,
104  const surfaceScalarField& phi,
105  const transportModel& transport,
106  const word& propertiesName,
107  const word& type
108 )
109 :
111  (
112  type,
113  alpha,
114  rho,
115  U,
116  alphaRhoPhi,
117  phi,
118  transport,
119  propertiesName
120  ),
121  v2fBase(),
122 
123  Cmu_
124  (
126  (
127  "Cmu",
128  this->coeffDict_,
129  0.22
130  )
131  ),
132  CmuKEps_
133  (
135  (
136  "CmuKEps",
137  this->coeffDict_,
138  0.09
139  )
140  ),
141  C1_
142  (
144  (
145  "C1",
146  this->coeffDict_,
147  1.4
148  )
149  ),
150  C2_
151  (
153  (
154  "C2",
155  this->coeffDict_,
156  0.3
157  )
158  ),
159  CL_
160  (
162  (
163  "CL",
164  this->coeffDict_,
165  0.23
166  )
167  ),
168  Ceta_
169  (
171  (
172  "Ceta",
173  this->coeffDict_,
174  70.0
175  )
176  ),
177  Ceps2_
178  (
180  (
181  "Ceps2",
182  this->coeffDict_,
183  1.9
184  )
185  ),
186  Ceps3_
187  (
189  (
190  "Ceps3",
191  this->coeffDict_,
192  -0.33
193  )
194  ),
195  sigmaK_
196  (
198  (
199  "sigmaK",
200  this->coeffDict_,
201  1.0
202  )
203  ),
204  sigmaEps_
205  (
207  (
208  "sigmaEps",
209  this->coeffDict_,
210  1.3
211  )
212  ),
213 
214  k_
215  (
216  IOobject
217  (
218  IOobject::groupName("k", U.group()),
219  this->runTime_.timeName(),
220  this->mesh_,
223  ),
224  this->mesh_
225  ),
226  epsilon_
227  (
228  IOobject
229  (
230  IOobject::groupName("epsilon", U.group()),
231  this->runTime_.timeName(),
232  this->mesh_,
235  ),
236  this->mesh_
237  ),
238  v2_
239  (
240  IOobject
241  (
242  IOobject::groupName("v2", U.group()),
243  this->runTime_.timeName(),
244  this->mesh_,
247  ),
248  this->mesh_
249  ),
250  f_
251  (
252  IOobject
253  (
254  IOobject::groupName("f", U.group()),
255  this->runTime_.timeName(),
256  this->mesh_,
259  ),
260  this->mesh_
261  ),
262  v2Min_(dimensionedScalar("v2Min", v2_.dimensions(), SMALL)),
263  fMin_(dimensionedScalar("fMin", f_.dimensions(), 0.0))
264 {
265  bound(k_, this->kMin_);
266  bound(epsilon_, this->epsilonMin_);
267  bound(v2_, v2Min_);
268  bound(f_, fMin_);
269 
270  if (type == typeName)
271  {
272  this->printCoeffs(type);
273  }
274 }
275 
276 
277 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
278 
279 template<class BasicTurbulenceModel>
281 {
283  {
284  Cmu_.readIfPresent(this->coeffDict());
285  CmuKEps_.readIfPresent(this->coeffDict());
286  C1_.readIfPresent(this->coeffDict());
287  C2_.readIfPresent(this->coeffDict());
288  CL_.readIfPresent(this->coeffDict());
289  Ceta_.readIfPresent(this->coeffDict());
290  Ceps2_.readIfPresent(this->coeffDict());
291  Ceps3_.readIfPresent(this->coeffDict());
292  sigmaK_.readIfPresent(this->coeffDict());
293  sigmaEps_.readIfPresent(this->coeffDict());
294 
295  return true;
296  }
297  else
298  {
299  return false;
300  }
301 }
302 
303 
304 template<class BasicTurbulenceModel>
306 {
307  if (!this->turbulence_)
308  {
309  return;
310  }
311 
312  // Local references
313  const alphaField& alpha = this->alpha_;
314  const rhoField& rho = this->rho_;
315  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
316  const volVectorField& U = this->U_;
317  volScalarField& nut = this->nut_;
318 
320 
322 
323  // Use N=6 so that f=0 at walls
324  const dimensionedScalar N("N", dimless, 6.0);
325 
326  const volTensorField gradU(fvc::grad(U));
327  const volScalarField S2(2*magSqr(dev(symm(gradU))));
328 
329  const volScalarField G(this->GName(), nut*S2);
330  const volScalarField Ts(this->Ts());
331  const volScalarField L2(type() + ":L2", sqr(Ls()));
332  const volScalarField v2fAlpha
333  (
334  type() + ":alpha",
335  1.0/Ts*((C1_ - N)*v2_ - 2.0/3.0*k_*(C1_ - 1.0))
336  );
337 
338  const volScalarField Ceps1
339  (
340  "Ceps1",
341  1.4*(1.0 + 0.05*min(sqrt(k_/v2_), scalar(100.0)))
342  );
343 
344  // Update epsilon (and possibly G) at the wall
345  epsilon_.boundaryField().updateCoeffs();
346 
347  // Dissipation equation
348  tmp<fvScalarMatrix> epsEqn
349  (
350  fvm::ddt(alpha, rho, epsilon_)
351  + fvm::div(alphaRhoPhi, epsilon_)
352  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
353  ==
354  Ceps1*alpha*rho*G/Ts
355  - fvm::SuSp(((2.0/3.0)*Ceps1 + Ceps3_)*alpha*rho*divU, epsilon_)
356  - fvm::Sp(Ceps2_*alpha*rho/Ts, epsilon_)
357  );
358 
359  epsEqn().relax();
360 
361  epsEqn().boundaryManipulate(epsilon_.boundaryField());
362 
363  solve(epsEqn);
364  bound(epsilon_, this->epsilonMin_);
365 
366 
367  // Turbulent kinetic energy equation
369  (
370  fvm::ddt(alpha, rho, k_)
371  + fvm::div(alphaRhoPhi, k_)
372  - fvm::laplacian(alpha*rho*DkEff(), k_)
373  ==
374  alpha*rho*G
375  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
376  - fvm::Sp(alpha*rho*epsilon_/k_, k_)
377  );
378 
379  kEqn().relax();
380  solve(kEqn);
381  bound(k_, this->kMin_);
382 
383 
384  // Relaxation function equation
386  (
387  - fvm::laplacian(f_)
388  ==
389  - fvm::Sp(1.0/L2, f_)
390  - 1.0/L2/k_*(v2fAlpha - C2_*G)
391  );
392 
393  fEqn().relax();
394  solve(fEqn);
395  bound(f_, fMin_);
396 
397 
398  // Turbulence stress normal to streamlines equation
399  tmp<fvScalarMatrix> v2Eqn
400  (
401  fvm::ddt(alpha, rho, v2_)
402  + fvm::div(alphaRhoPhi, v2_)
403  - fvm::laplacian(alpha*rho*DkEff(), v2_)
404  ==
405  alpha*rho*min(k_*f_, C2_*G - v2fAlpha)
406  - fvm::Sp(N*alpha*rho*epsilon_/k_, v2_)
407  );
408 
409  v2Eqn().relax();
410  solve(v2Eqn);
411  bound(v2_, v2Min_);
412 
413  correctNut();
414 }
415 
416 
417 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
418 
419 } // End namespace RASModels
420 } // End namespace Foam
421 
422 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:82
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::RASModel
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::RASModels::v2f::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: v2f.C:280
U
U
Definition: pEqn.H:46
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::RASModels::v2f::Ts
tmp< volScalarField > Ts() const
Return time scale, Ts.
Definition: v2f.C:39
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:131
nu
volScalarField & nu
Definition: readMechanicalProperties.H:179
nut
nut
Definition: createTDFields.H:71
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:87
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
bound.H
Bound the given scalar field if it has gone unbounded.
correct
fvOptions correct(rho)
Foam::fvm::SuSp
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
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:99
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
Foam::RASModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:100
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::RASModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:99
rho
rho
Definition: pEqn.H:3
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::RASModels::v2f::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: v2f.C:305
Foam::RASModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:101
Foam::solve
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
divU
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::RASModels::v2f::correctNut
virtual void correctNut()
Definition: v2f.C:86
v2f.H
Foam::eddyViscosity
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::RASModels::v2fBase
Abstract base-class for v2-f models to provide BCs access to the v2 and f fields.
Definition: v2fBase.H:57
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:187
Foam::RASModels::v2f::Ls
tmp< volScalarField > Ls() const
Return length scale, Ls.
Definition: v2f.C:61
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:104
N
label N
Definition: createFields.H:22