All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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  tmp<volScalarField> tCmuk2(CmuKEps_*sqr(k_));
44  epsilon_ = max(epsilon_, tCmuk2()/(this->nutMaxCoeff_*this->nu()));
45  return tCmuk2;
46 }
47 
48 
49 template<class BasicMomentumTransportModel>
51 {
52  return max(k_/epsilon_, 6.0*sqrt(this->nu()/epsilon_));
53 }
54 
55 
56 template<class BasicMomentumTransportModel>
58 {
59  return
60  CL_*max(pow(k_, 1.5)
61  /epsilon_, Ceta_*pow025(pow3(this->nu())/epsilon_));
62 }
63 
64 
65 template<class BasicMomentumTransportModel>
67 {
68  this->nut_ = min(boundEpsilon()/epsilon_, this->Cmu_*v2_*Ts());
69  this->nut_.correctBoundaryConditions();
70  fvConstraints::New(this->mesh_).constrain(this->nut_);
71 }
72 
73 
74 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
75 
76 template<class BasicMomentumTransportModel>
78 (
79  const alphaField& alpha,
80  const rhoField& rho,
81  const volVectorField& U,
82  const surfaceScalarField& alphaRhoPhi,
83  const surfaceScalarField& phi,
84  const viscosity& viscosity,
85  const word& type
86 )
87 :
88  eddyViscosity<RASModel<BasicMomentumTransportModel>>
89  (
90  type,
91  alpha,
92  rho,
93  U,
94  alphaRhoPhi,
95  phi,
96  viscosity
97  ),
98  v2fBase(),
99 
100  Cmu_
101  (
102  dimensioned<scalar>::lookupOrAddToDict
103  (
104  "Cmu",
105  this->coeffDict_,
106  0.22
107  )
108  ),
109  CmuKEps_
110  (
111  dimensioned<scalar>::lookupOrAddToDict
112  (
113  "CmuKEps",
114  this->coeffDict_,
115  0.09
116  )
117  ),
118  C1_
119  (
120  dimensioned<scalar>::lookupOrAddToDict
121  (
122  "C1",
123  this->coeffDict_,
124  1.4
125  )
126  ),
127  C2_
128  (
129  dimensioned<scalar>::lookupOrAddToDict
130  (
131  "C2",
132  this->coeffDict_,
133  0.3
134  )
135  ),
136  CL_
137  (
138  dimensioned<scalar>::lookupOrAddToDict
139  (
140  "CL",
141  this->coeffDict_,
142  0.23
143  )
144  ),
145  Ceta_
146  (
147  dimensioned<scalar>::lookupOrAddToDict
148  (
149  "Ceta",
150  this->coeffDict_,
151  70.0
152  )
153  ),
154  Ceps2_
155  (
156  dimensioned<scalar>::lookupOrAddToDict
157  (
158  "Ceps2",
159  this->coeffDict_,
160  1.9
161  )
162  ),
163  Ceps3_
164  (
165  dimensioned<scalar>::lookupOrAddToDict
166  (
167  "Ceps3",
168  this->coeffDict_,
169  -0.33
170  )
171  ),
172  sigmaK_
173  (
174  dimensioned<scalar>::lookupOrAddToDict
175  (
176  "sigmaK",
177  this->coeffDict_,
178  1.0
179  )
180  ),
181  sigmaEps_
182  (
183  dimensioned<scalar>::lookupOrAddToDict
184  (
185  "sigmaEps",
186  this->coeffDict_,
187  1.3
188  )
189  ),
190 
191  k_
192  (
193  IOobject
194  (
195  this->groupName("k"),
196  this->runTime_.name(),
197  this->mesh_,
198  IOobject::MUST_READ,
199  IOobject::AUTO_WRITE
200  ),
201  this->mesh_
202  ),
203  epsilon_
204  (
205  IOobject
206  (
207  this->groupName("epsilon"),
208  this->runTime_.name(),
209  this->mesh_,
210  IOobject::MUST_READ,
211  IOobject::AUTO_WRITE
212  ),
213  this->mesh_
214  ),
215  v2_
216  (
217  IOobject
218  (
219  this->groupName("v2"),
220  this->runTime_.name(),
221  this->mesh_,
222  IOobject::MUST_READ,
223  IOobject::AUTO_WRITE
224  ),
225  this->mesh_
226  ),
227  f_
228  (
229  IOobject
230  (
231  this->groupName("f"),
232  this->runTime_.name(),
233  this->mesh_,
234  IOobject::MUST_READ,
235  IOobject::AUTO_WRITE
236  ),
237  this->mesh_
238  ),
239  v2Min_(dimensionedScalar(v2_.dimensions(), small)),
240  fMin_(dimensionedScalar(f_.dimensions(), 0))
241 {
242  bound(k_, this->kMin_);
243  boundEpsilon();
244  bound(v2_, v2Min_);
245  bound(f_, fMin_);
246 
247  if (type == typeName)
248  {
249  this->printCoeffs(type);
250  }
251 }
252 
253 
254 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
255 
256 template<class BasicMomentumTransportModel>
258 {
260  {
261  Cmu_.readIfPresent(this->coeffDict());
262  CmuKEps_.readIfPresent(this->coeffDict());
263  C1_.readIfPresent(this->coeffDict());
264  C2_.readIfPresent(this->coeffDict());
265  CL_.readIfPresent(this->coeffDict());
266  Ceta_.readIfPresent(this->coeffDict());
267  Ceps2_.readIfPresent(this->coeffDict());
268  Ceps3_.readIfPresent(this->coeffDict());
269  sigmaK_.readIfPresent(this->coeffDict());
270  sigmaEps_.readIfPresent(this->coeffDict());
271 
272  return true;
273  }
274  else
275  {
276  return false;
277  }
278 }
279 
280 
281 template<class BasicMomentumTransportModel>
283 {
284  if (!this->turbulence_)
285  {
286  return;
287  }
288 
289  // Local references
290  const alphaField& alpha = this->alpha_;
291  const rhoField& rho = this->rho_;
292  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
293  const volVectorField& U = this->U_;
294  volScalarField& nut = this->nut_;
295  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
297  (
298  Foam::fvConstraints::New(this->mesh_)
299  );
300 
302 
303  volScalarField divU(fvc::div(fvc::absolute(this->phi(), U)));
304 
305  // Use N=6 so that f=0 at walls
306  const dimensionedScalar N("N", dimless, 6.0);
307 
308  const volTensorField gradU(fvc::grad(U));
309  const volScalarField S2(2*magSqr(dev(symm(gradU))));
310 
311  const volScalarField G(this->GName(), nut*S2);
312  const volScalarField Ts(this->Ts());
313  const volScalarField L2(typedName("L2"), sqr(Ls()));
314  const volScalarField v2fAlpha
315  (
316  typedName("alpha"),
317  1.0/Ts*((C1_ - N)*v2_ - 2.0/3.0*k_*(C1_ - 1.0))
318  );
319 
320  const volScalarField Ceps1
321  (
322  typedName("Ceps1"),
323  1.4*(1.0 + 0.05*min(sqrt(k_/v2_), scalar(100.0)))
324  );
325 
326  // Update epsilon (and possibly G) at the wall
327  epsilon_.boundaryFieldRef().updateCoeffs();
328 
329  // Dissipation equation
330  tmp<fvScalarMatrix> epsEqn
331  (
332  fvm::ddt(alpha, rho, epsilon_)
333  + fvm::div(alphaRhoPhi, epsilon_)
334  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
335  ==
336  Ceps1*alpha*rho*G/Ts
337  - fvm::SuSp(((2.0/3.0)*Ceps1 + Ceps3_)*alpha*rho*divU, epsilon_)
338  - fvm::Sp(Ceps2_*alpha*rho/Ts, epsilon_)
339  + fvModels.source(alpha, rho, epsilon_)
340  );
341 
342  epsEqn.ref().relax();
343  fvConstraints.constrain(epsEqn.ref());
344  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
345  solve(epsEqn);
346  fvConstraints.constrain(epsilon_);
347  boundEpsilon();
348 
349 
350  // Turbulent kinetic energy equation
352  (
353  fvm::ddt(alpha, rho, k_)
354  + fvm::div(alphaRhoPhi, k_)
355  - fvm::laplacian(alpha*rho*DkEff(), k_)
356  ==
357  alpha*rho*G
358  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
359  - fvm::Sp(alpha*rho*epsilon_/k_, k_)
360  + fvModels.source(alpha, rho, k_)
361  );
362 
363  kEqn.ref().relax();
364  fvConstraints.constrain(kEqn.ref());
365  solve(kEqn);
367  bound(k_, this->kMin_);
368 
369 
370  // Relaxation function equation
372  (
373  - fvm::laplacian(f_)
374  ==
375  - fvm::Sp(1.0/L2, f_)
376  - 1.0/L2/k_*(v2fAlpha - C2_*G)
377  );
378 
379  fEqn.ref().relax();
380  fvConstraints.constrain(fEqn.ref());
381  solve(fEqn);
383  bound(f_, fMin_);
384 
385 
386  // Turbulence stress normal to streamlines equation
387  tmp<fvScalarMatrix> v2Eqn
388  (
389  fvm::ddt(alpha, rho, v2_)
390  + fvm::div(alphaRhoPhi, v2_)
391  - fvm::laplacian(alpha*rho*DkEff(), v2_)
392  ==
393  alpha*rho*min(k_*f_, C2_*G - v2fAlpha)
394  - fvm::Sp(N*alpha*rho*epsilon_/k_, v2_)
395  + fvModels.source(alpha, rho, v2_)
396  );
397 
398  v2Eqn.ref().relax();
399  fvConstraints.constrain(v2Eqn.ref());
400  solve(v2Eqn);
402  bound(v2_, v2Min_);
403 
404  correctNut();
405 }
406 
407 
408 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
409 
410 } // End namespace RASModels
411 } // End namespace Foam
412 
413 // ************************************************************************* //
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 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:78
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: v2f.C:282
tmp< volScalarField > boundEpsilon()
Bound epsilon and return Cmu*sqr(k) for nut.
Definition: v2f.C:41
dimensionedScalar fMin_
Definition: v2f.H:159
dimensionedScalar v2Min_
Definition: v2f.H:158
tmp< volScalarField > Ts() const
Return time scale, Ts.
Definition: v2f.C:50
virtual void correctNut()
Correct the eddy-viscosity nut.
Definition: v2f.C:66
volScalarField f_
Damping function.
Definition: v2f.H:153
tmp< volScalarField > Ls() const
Return length scale, Ls.
Definition: v2f.C:57
virtual bool read()
Read RASProperties dictionary.
Definition: v2f.C:257
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:67
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Finite volume models.
Definition: fvModels.H:65
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
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 > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
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)
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
dimensionedScalar pow3(const dimensionedScalar &ds)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
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:176
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
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)