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-2023 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 "fvModels.H"
28 #include "fvConstraints.H"
29 #include "bound.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace RASModels
36 {
37 
38 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
39 
40 template<class BasicMomentumTransportModel>
42 {
43  return max(k_/epsilon_, 6.0*sqrt(this->nu()/epsilon_));
44 }
45 
46 
47 template<class BasicMomentumTransportModel>
49 {
50  return
51  CL_*max(pow(k_, 1.5)
52  /epsilon_, Ceta_*pow025(pow3(this->nu())/epsilon_));
53 }
54 
55 
56 template<class BasicMomentumTransportModel>
58 {
59  this->nut_ = min(CmuKEps_*sqr(k_)/epsilon_, this->Cmu_*v2_*Ts());
60  this->nut_.correctBoundaryConditions();
61  fvConstraints::New(this->mesh_).constrain(this->nut_);
62 }
63 
64 
65 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
66 
67 template<class BasicMomentumTransportModel>
69 (
70  const alphaField& alpha,
71  const rhoField& rho,
72  const volVectorField& U,
73  const surfaceScalarField& alphaRhoPhi,
74  const surfaceScalarField& phi,
75  const viscosity& viscosity,
76  const word& type
77 )
78 :
79  eddyViscosity<RASModel<BasicMomentumTransportModel>>
80  (
81  type,
82  alpha,
83  rho,
84  U,
85  alphaRhoPhi,
86  phi,
87  viscosity
88  ),
89  v2fBase(),
90 
91  Cmu_
92  (
93  dimensioned<scalar>::lookupOrAddToDict
94  (
95  "Cmu",
96  this->coeffDict_,
97  0.22
98  )
99  ),
100  CmuKEps_
101  (
102  dimensioned<scalar>::lookupOrAddToDict
103  (
104  "CmuKEps",
105  this->coeffDict_,
106  0.09
107  )
108  ),
109  C1_
110  (
111  dimensioned<scalar>::lookupOrAddToDict
112  (
113  "C1",
114  this->coeffDict_,
115  1.4
116  )
117  ),
118  C2_
119  (
120  dimensioned<scalar>::lookupOrAddToDict
121  (
122  "C2",
123  this->coeffDict_,
124  0.3
125  )
126  ),
127  CL_
128  (
129  dimensioned<scalar>::lookupOrAddToDict
130  (
131  "CL",
132  this->coeffDict_,
133  0.23
134  )
135  ),
136  Ceta_
137  (
138  dimensioned<scalar>::lookupOrAddToDict
139  (
140  "Ceta",
141  this->coeffDict_,
142  70.0
143  )
144  ),
145  Ceps2_
146  (
147  dimensioned<scalar>::lookupOrAddToDict
148  (
149  "Ceps2",
150  this->coeffDict_,
151  1.9
152  )
153  ),
154  Ceps3_
155  (
156  dimensioned<scalar>::lookupOrAddToDict
157  (
158  "Ceps3",
159  this->coeffDict_,
160  -0.33
161  )
162  ),
163  sigmaK_
164  (
165  dimensioned<scalar>::lookupOrAddToDict
166  (
167  "sigmaK",
168  this->coeffDict_,
169  1.0
170  )
171  ),
172  sigmaEps_
173  (
174  dimensioned<scalar>::lookupOrAddToDict
175  (
176  "sigmaEps",
177  this->coeffDict_,
178  1.3
179  )
180  ),
181 
182  k_
183  (
184  IOobject
185  (
186  this->groupName("k"),
187  this->runTime_.name(),
188  this->mesh_,
189  IOobject::MUST_READ,
190  IOobject::AUTO_WRITE
191  ),
192  this->mesh_
193  ),
194  epsilon_
195  (
196  IOobject
197  (
198  this->groupName("epsilon"),
199  this->runTime_.name(),
200  this->mesh_,
201  IOobject::MUST_READ,
202  IOobject::AUTO_WRITE
203  ),
204  this->mesh_
205  ),
206  v2_
207  (
208  IOobject
209  (
210  this->groupName("v2"),
211  this->runTime_.name(),
212  this->mesh_,
213  IOobject::MUST_READ,
214  IOobject::AUTO_WRITE
215  ),
216  this->mesh_
217  ),
218  f_
219  (
220  IOobject
221  (
222  this->groupName("f"),
223  this->runTime_.name(),
224  this->mesh_,
225  IOobject::MUST_READ,
226  IOobject::AUTO_WRITE
227  ),
228  this->mesh_
229  ),
230  v2Min_(dimensionedScalar(v2_.dimensions(), small)),
231  fMin_(dimensionedScalar(f_.dimensions(), 0))
232 {
233  bound(k_, this->kMin_);
234  bound(epsilon_, this->epsilonMin_);
235  bound(v2_, v2Min_);
236  bound(f_, fMin_);
237 
238  if (type == typeName)
239  {
240  this->printCoeffs(type);
241  }
242 }
243 
244 
245 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
246 
247 template<class BasicMomentumTransportModel>
249 {
251  {
252  Cmu_.readIfPresent(this->coeffDict());
253  CmuKEps_.readIfPresent(this->coeffDict());
254  C1_.readIfPresent(this->coeffDict());
255  C2_.readIfPresent(this->coeffDict());
256  CL_.readIfPresent(this->coeffDict());
257  Ceta_.readIfPresent(this->coeffDict());
258  Ceps2_.readIfPresent(this->coeffDict());
259  Ceps3_.readIfPresent(this->coeffDict());
260  sigmaK_.readIfPresent(this->coeffDict());
261  sigmaEps_.readIfPresent(this->coeffDict());
262 
263  return true;
264  }
265  else
266  {
267  return false;
268  }
269 }
270 
271 
272 template<class BasicMomentumTransportModel>
274 {
275  if (!this->turbulence_)
276  {
277  return;
278  }
279 
280  // Local references
281  const alphaField& alpha = this->alpha_;
282  const rhoField& rho = this->rho_;
283  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
284  const volVectorField& U = this->U_;
285  volScalarField& nut = this->nut_;
286  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
288  (
289  Foam::fvConstraints::New(this->mesh_)
290  );
291 
293 
294  volScalarField divU(fvc::div(fvc::absolute(this->phi(), U)));
295 
296  // Use N=6 so that f=0 at walls
297  const dimensionedScalar N("N", dimless, 6.0);
298 
299  const volTensorField gradU(fvc::grad(U));
300  const volScalarField S2(2*magSqr(dev(symm(gradU))));
301 
302  const volScalarField G(this->GName(), nut*S2);
303  const volScalarField Ts(this->Ts());
304  const volScalarField L2(typedName("L2"), sqr(Ls()));
305  const volScalarField v2fAlpha
306  (
307  typedName("alpha"),
308  1.0/Ts*((C1_ - N)*v2_ - 2.0/3.0*k_*(C1_ - 1.0))
309  );
310 
311  const volScalarField Ceps1
312  (
313  typedName("Ceps1"),
314  1.4*(1.0 + 0.05*min(sqrt(k_/v2_), scalar(100.0)))
315  );
316 
317  // Update epsilon (and possibly G) at the wall
318  epsilon_.boundaryFieldRef().updateCoeffs();
319 
320  // Dissipation equation
321  tmp<fvScalarMatrix> epsEqn
322  (
323  fvm::ddt(alpha, rho, epsilon_)
324  + fvm::div(alphaRhoPhi, epsilon_)
325  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
326  ==
327  Ceps1*alpha*rho*G/Ts
328  - fvm::SuSp(((2.0/3.0)*Ceps1 + Ceps3_)*alpha*rho*divU, epsilon_)
329  - fvm::Sp(Ceps2_*alpha*rho/Ts, epsilon_)
330  + fvModels.source(alpha, rho, epsilon_)
331  );
332 
333  epsEqn.ref().relax();
334  fvConstraints.constrain(epsEqn.ref());
335  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
336  solve(epsEqn);
337  fvConstraints.constrain(epsilon_);
338  bound(epsilon_, this->epsilonMin_);
339 
340 
341  // Turbulent kinetic energy equation
343  (
344  fvm::ddt(alpha, rho, k_)
345  + fvm::div(alphaRhoPhi, k_)
346  - fvm::laplacian(alpha*rho*DkEff(), k_)
347  ==
348  alpha*rho*G
349  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
350  - fvm::Sp(alpha*rho*epsilon_/k_, k_)
351  + fvModels.source(alpha, rho, k_)
352  );
353 
354  kEqn.ref().relax();
355  fvConstraints.constrain(kEqn.ref());
356  solve(kEqn);
358  bound(k_, this->kMin_);
359 
360 
361  // Relaxation function equation
363  (
364  - fvm::laplacian(f_)
365  ==
366  - fvm::Sp(1.0/L2, f_)
367  - 1.0/L2/k_*(v2fAlpha - C2_*G)
368  );
369 
370  fEqn.ref().relax();
371  fvConstraints.constrain(fEqn.ref());
372  solve(fEqn);
374  bound(f_, fMin_);
375 
376 
377  // Turbulence stress normal to streamlines equation
378  tmp<fvScalarMatrix> v2Eqn
379  (
380  fvm::ddt(alpha, rho, v2_)
381  + fvm::div(alphaRhoPhi, v2_)
382  - fvm::laplacian(alpha*rho*DkEff(), v2_)
383  ==
384  alpha*rho*min(k_*f_, C2_*G - v2fAlpha)
385  - fvm::Sp(N*alpha*rho*epsilon_/k_, v2_)
386  + fvModels.source(alpha, rho, v2_)
387  );
388 
389  v2Eqn.ref().relax();
390  fvConstraints.constrain(v2Eqn.ref());
391  solve(v2Eqn);
393  bound(v2_, v2Min_);
394 
395  correctNut();
396 }
397 
398 
399 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
400 
401 } // End namespace RASModels
402 } // End namespace Foam
403 
404 // ************************************************************************* //
Bound the given scalar field where it is below the specified minimum.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:56
Abstract base-class for v2-f models to provide BCs access to the v2 and f fields.
Definition: v2fBase.H:55
volScalarField epsilon_
Turbulence dissipation.
Definition: v2f.H:147
volScalarField k_
Turbulence kinetic energy.
Definition: v2f.H:144
v2f(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
Definition: v2f.C:69
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: v2f.C:273
dimensionedScalar fMin_
Definition: v2f.H:159
dimensionedScalar v2Min_
Definition: v2f.H:158
tmp< volScalarField > Ts() const
Return time scale, Ts.
Definition: v2f.C:41
virtual void correctNut()
Definition: v2f.C:57
volScalarField f_
Damping function.
Definition: v2f.H:153
tmp< volScalarField > Ls() const
Return length scale, Ls.
Definition: v2f.C:48
virtual bool read()
Read RASProperties dictionary.
Definition: v2f.C:248
volScalarField v2_
Turbulence stress normal to streamlines.
Definition: v2f.H:150
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:100
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
Finite volume constraints.
Definition: fvConstraints.H:62
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Finite volume models.
Definition: fvModels.H:65
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
const scalar nut
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimless
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:149
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
bool bound(volScalarField &, const dimensionedScalar &min)
Bound the given scalar field where it is below the specified min value.
Definition: bound.C:31
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dimensioned< scalar > magSqr(const dimensioned< Type > &)
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.
dimensionedScalar pow025(const dimensionedScalar &ds)