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