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-2024 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_("Cmu", this->coeffDict(), 0.22),
101  CmuKEps_("CmuKEps", this->coeffDict(), 0.09),
102  C1_("C1", this->coeffDict(), 1.4),
103  C2_("C2", this->coeffDict(), 0.3),
104  CL_("CL", this->coeffDict(), 0.23),
105  Ceta_("Ceta", this->coeffDict(), 70.0),
106  Ceps2_("Ceps2", this->coeffDict(), 1.9),
107  Ceps3_("Ceps3", this->coeffDict(), -0.33),
108  sigmaK_("sigmaK", this->coeffDict(), 1.0),
109  sigmaEps_("sigmaEps", this->coeffDict(), 1.3),
110 
111  k_
112  (
113  IOobject
114  (
115  this->groupName("k"),
116  this->runTime_.name(),
117  this->mesh_,
118  IOobject::MUST_READ,
119  IOobject::AUTO_WRITE
120  ),
121  this->mesh_
122  ),
123  epsilon_
124  (
125  IOobject
126  (
127  this->groupName("epsilon"),
128  this->runTime_.name(),
129  this->mesh_,
130  IOobject::MUST_READ,
131  IOobject::AUTO_WRITE
132  ),
133  this->mesh_
134  ),
135  v2_
136  (
137  IOobject
138  (
139  this->groupName("v2"),
140  this->runTime_.name(),
141  this->mesh_,
142  IOobject::MUST_READ,
143  IOobject::AUTO_WRITE
144  ),
145  this->mesh_
146  ),
147  f_
148  (
149  IOobject
150  (
151  this->groupName("f"),
152  this->runTime_.name(),
153  this->mesh_,
154  IOobject::MUST_READ,
155  IOobject::AUTO_WRITE
156  ),
157  this->mesh_
158  ),
159  v2Min_(dimensionedScalar(v2_.dimensions(), small)),
160  fMin_(dimensionedScalar(f_.dimensions(), 0))
161 {
162  bound(k_, this->kMin_);
163  boundEpsilon();
164  bound(v2_, v2Min_);
165  bound(f_, fMin_);
166 }
167 
168 
169 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
170 
171 template<class BasicMomentumTransportModel>
173 {
175  {
176  Cmu_.readIfPresent(this->coeffDict());
177  CmuKEps_.readIfPresent(this->coeffDict());
178  C1_.readIfPresent(this->coeffDict());
179  C2_.readIfPresent(this->coeffDict());
180  CL_.readIfPresent(this->coeffDict());
181  Ceta_.readIfPresent(this->coeffDict());
182  Ceps2_.readIfPresent(this->coeffDict());
183  Ceps3_.readIfPresent(this->coeffDict());
184  sigmaK_.readIfPresent(this->coeffDict());
185  sigmaEps_.readIfPresent(this->coeffDict());
186 
187  return true;
188  }
189  else
190  {
191  return false;
192  }
193 }
194 
195 
196 template<class BasicMomentumTransportModel>
198 {
199  if (!this->turbulence_)
200  {
201  return;
202  }
203 
204  // Local references
205  const alphaField& alpha = this->alpha_;
206  const rhoField& rho = this->rho_;
207  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
208  const volVectorField& U = this->U_;
209  volScalarField& nut = this->nut_;
210  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
212  (
213  Foam::fvConstraints::New(this->mesh_)
214  );
215 
217 
218  volScalarField divU(fvc::div(fvc::absolute(this->phi(), U)));
219 
220  // Use N=6 so that f=0 at walls
221  const dimensionedScalar N("N", dimless, 6.0);
222 
223  const volTensorField gradU(fvc::grad(U));
224  const volScalarField S2(2*magSqr(dev(symm(gradU))));
225 
226  const volScalarField G(this->GName(), nut*S2);
227  const volScalarField Ts(this->Ts());
228  const volScalarField L2(typedName("L2"), sqr(Ls()));
229  const volScalarField v2fAlpha
230  (
231  typedName("alpha"),
232  1.0/Ts*((C1_ - N)*v2_ - 2.0/3.0*k_*(C1_ - 1.0))
233  );
234 
235  const volScalarField Ceps1
236  (
237  typedName("Ceps1"),
238  1.4*(1.0 + 0.05*min(sqrt(k_/v2_), scalar(100.0)))
239  );
240 
241  // Update epsilon (and possibly G) at the wall
242  epsilon_.boundaryFieldRef().updateCoeffs();
243 
244  // Dissipation equation
245  tmp<fvScalarMatrix> epsEqn
246  (
247  fvm::ddt(alpha, rho, epsilon_)
248  + fvm::div(alphaRhoPhi, epsilon_)
249  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
250  ==
251  Ceps1*alpha*rho*G/Ts
252  - fvm::SuSp(((2.0/3.0)*Ceps1 + Ceps3_)*alpha*rho*divU, epsilon_)
253  - fvm::Sp(Ceps2_*alpha*rho/Ts, epsilon_)
254  + fvModels.source(alpha, rho, epsilon_)
255  );
256 
257  epsEqn.ref().relax();
258  fvConstraints.constrain(epsEqn.ref());
259  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
260  solve(epsEqn);
261  fvConstraints.constrain(epsilon_);
262  boundEpsilon();
263 
264 
265  // Turbulent kinetic energy equation
267  (
268  fvm::ddt(alpha, rho, k_)
269  + fvm::div(alphaRhoPhi, k_)
270  - fvm::laplacian(alpha*rho*DkEff(), k_)
271  ==
272  alpha*rho*G
273  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
274  - fvm::Sp(alpha*rho*epsilon_/k_, k_)
275  + fvModels.source(alpha, rho, k_)
276  );
277 
278  kEqn.ref().relax();
279  fvConstraints.constrain(kEqn.ref());
280  solve(kEqn);
282  bound(k_, this->kMin_);
283 
284 
285  // Relaxation function equation
287  (
288  - fvm::laplacian(f_)
289  ==
290  - fvm::Sp(1.0/L2, f_)
291  - 1.0/L2/k_*(v2fAlpha - C2_*G)
292  );
293 
294  fEqn.ref().relax();
295  fvConstraints.constrain(fEqn.ref());
296  solve(fEqn);
298  bound(f_, fMin_);
299 
300 
301  // Turbulence stress normal to streamlines equation
302  tmp<fvScalarMatrix> v2Eqn
303  (
304  fvm::ddt(alpha, rho, v2_)
305  + fvm::div(alphaRhoPhi, v2_)
306  - fvm::laplacian(alpha*rho*DkEff(), v2_)
307  ==
308  alpha*rho*min(k_*f_, C2_*G - v2fAlpha)
309  - fvm::Sp(N*alpha*rho*epsilon_/k_, v2_)
310  + fvModels.source(alpha, rho, v2_)
311  );
312 
313  v2Eqn.ref().relax();
314  fvConstraints.constrain(v2Eqn.ref());
315  solve(v2Eqn);
317  bound(v2_, v2Min_);
318 
319  correctNut();
320 }
321 
322 
323 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
324 
325 } // End namespace RASModels
326 } // End namespace Foam
327 
328 // ************************************************************************* //
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:197
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:172
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:103
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:197
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 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.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
void dev(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void pow025(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
const dimensionSet dimless
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void symm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:181
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
void pow3(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
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.