SpalartAllmaras.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) 2011-2016 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 "SpalartAllmaras.H"
27 #include "fvOptions.H"
28 #include "bound.H"
29 #include "wallDist.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace RASModels
36 {
37 
38 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
39 
40 template<class BasicTurbulenceModel>
42 {
43  return nuTilda_/this->nu();
44 }
45 
46 
47 template<class BasicTurbulenceModel>
49 (
50  const volScalarField& chi
51 ) const
52 {
53  const volScalarField chi3(pow3(chi));
54  return chi3/(chi3 + pow3(Cv1_));
55 }
56 
57 
58 template<class BasicTurbulenceModel>
60 (
61  const volScalarField& chi,
62  const volScalarField& fv1
63 ) const
64 {
65  return 1.0 - chi/(1.0 + chi*fv1);
66 }
67 
68 
69 template<class BasicTurbulenceModel>
71 (
72  const volScalarField& chi,
73  const volScalarField& fv1
74 ) const
75 {
76  volScalarField Omega(::sqrt(2.0)*mag(skew(fvc::grad(this->U_))));
77 
78  return
79  (
80  max
81  (
82  Omega
83  + fv2(chi, fv1)*nuTilda_/sqr(kappa_*y_),
84  Cs_*Omega
85  )
86  );
87 }
88 
89 
90 template<class BasicTurbulenceModel>
92 (
93  const volScalarField& Stilda
94 ) const
95 {
97  (
98  min
99  (
100  nuTilda_
101  /(
102  max
103  (
104  Stilda,
105  dimensionedScalar("SMALL", Stilda.dimensions(), SMALL)
106  )
107  *sqr(kappa_*y_)
108  ),
109  scalar(10.0)
110  )
111  );
112  r.boundaryFieldRef() == 0.0;
113 
114  const volScalarField g(r + Cw2_*(pow6(r) - r));
115 
116  return g*pow((1.0 + pow6(Cw3_))/(pow6(g) + pow6(Cw3_)), 1.0/6.0);
117 }
118 
119 
120 template<class BasicTurbulenceModel>
122 (
123  const volScalarField& fv1
124 )
125 {
126  this->nut_ = nuTilda_*fv1;
127  this->nut_.correctBoundaryConditions();
128  fv::options::New(this->mesh_).correct(this->nut_);
129 
130  BasicTurbulenceModel::correctNut();
131 }
132 
133 
134 template<class BasicTurbulenceModel>
136 {
137  correctNut(fv1(this->chi()));
138 }
139 
140 
141 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
142 
143 template<class BasicTurbulenceModel>
145 (
146  const alphaField& alpha,
147  const rhoField& rho,
148  const volVectorField& U,
149  const surfaceScalarField& alphaRhoPhi,
150  const surfaceScalarField& phi,
151  const transportModel& transport,
152  const word& propertiesName,
153  const word& type
154 )
155 :
157  (
158  type,
159  alpha,
160  rho,
161  U,
162  alphaRhoPhi,
163  phi,
164  transport,
165  propertiesName
166  ),
167 
168  sigmaNut_
169  (
171  (
172  "sigmaNut",
173  this->coeffDict_,
174  0.66666
175  )
176  ),
177  kappa_
178  (
180  (
181  "kappa",
182  this->coeffDict_,
183  0.41
184  )
185  ),
186 
187  Cb1_
188  (
190  (
191  "Cb1",
192  this->coeffDict_,
193  0.1355
194  )
195  ),
196  Cb2_
197  (
199  (
200  "Cb2",
201  this->coeffDict_,
202  0.622
203  )
204  ),
205  Cw1_(Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_),
206  Cw2_
207  (
209  (
210  "Cw2",
211  this->coeffDict_,
212  0.3
213  )
214  ),
215  Cw3_
216  (
218  (
219  "Cw3",
220  this->coeffDict_,
221  2.0
222  )
223  ),
224  Cv1_
225  (
227  (
228  "Cv1",
229  this->coeffDict_,
230  7.1
231  )
232  ),
233  Cs_
234  (
236  (
237  "Cs",
238  this->coeffDict_,
239  0.3
240  )
241  ),
242 
243  nuTilda_
244  (
245  IOobject
246  (
247  "nuTilda",
248  this->runTime_.timeName(),
249  this->mesh_,
252  ),
253  this->mesh_
254  ),
255 
256  y_(wallDist::New(this->mesh_).y())
257 {
258  if (type == typeName)
259  {
260  this->printCoeffs(type);
261  }
262 }
263 
264 
265 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
266 
267 template<class BasicTurbulenceModel>
269 {
271  {
272  sigmaNut_.readIfPresent(this->coeffDict());
273  kappa_.readIfPresent(this->coeffDict());
274 
275  Cb1_.readIfPresent(this->coeffDict());
276  Cb2_.readIfPresent(this->coeffDict());
277  Cw1_ = Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_;
278  Cw2_.readIfPresent(this->coeffDict());
279  Cw3_.readIfPresent(this->coeffDict());
280  Cv1_.readIfPresent(this->coeffDict());
281  Cs_.readIfPresent(this->coeffDict());
282 
283  return true;
284  }
285  else
286  {
287  return false;
288  }
289 }
290 
291 
292 template<class BasicTurbulenceModel>
294 {
295  return tmp<volScalarField>
296  (
297  new volScalarField("DnuTildaEff", (nuTilda_ + this->nu())/sigmaNut_)
298  );
299 }
300 
301 
302 template<class BasicTurbulenceModel>
304 {
305  return tmp<volScalarField>
306  (
307  new volScalarField
308  (
309  IOobject
310  (
311  "k",
312  this->runTime_.timeName(),
313  this->mesh_
314  ),
315  this->mesh_,
316  dimensionedScalar("0", dimensionSet(0, 2, -2, 0, 0), 0)
317  )
318  );
319 }
320 
321 
322 template<class BasicTurbulenceModel>
324 {
326  << "Turbulence kinetic energy dissipation rate not defined for "
327  << "Spalart-Allmaras model. Returning zero field"
328  << endl;
329 
330  return tmp<volScalarField>
331  (
332  new volScalarField
333  (
334  IOobject
335  (
336  "epsilon",
337  this->runTime_.timeName(),
338  this->mesh_
339  ),
340  this->mesh_,
341  dimensionedScalar("0", dimensionSet(0, 2, -3, 0, 0), 0)
342  )
343  );
344 }
345 
346 
347 template<class BasicTurbulenceModel>
349 {
350  if (!this->turbulence_)
351  {
352  return;
353  }
354 
355  // Local references
356  const alphaField& alpha = this->alpha_;
357  const rhoField& rho = this->rho_;
358  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
359  fv::options& fvOptions(fv::options::New(this->mesh_));
360 
362 
363  const volScalarField chi(this->chi());
364  const volScalarField fv1(this->fv1(chi));
365 
366  const volScalarField Stilda(this->Stilda(chi, fv1));
367 
368  tmp<fvScalarMatrix> nuTildaEqn
369  (
370  fvm::ddt(alpha, rho, nuTilda_)
371  + fvm::div(alphaRhoPhi, nuTilda_)
372  - fvm::laplacian(alpha*rho*DnuTildaEff(), nuTilda_)
373  - Cb2_/sigmaNut_*alpha*rho*magSqr(fvc::grad(nuTilda_))
374  ==
375  Cb1_*alpha*rho*Stilda*nuTilda_
376  - fvm::Sp(Cw1_*alpha*rho*fw(Stilda)*nuTilda_/sqr(y_), nuTilda_)
377  + fvOptions(alpha, rho, nuTilda_)
378  );
379 
380  nuTildaEqn.ref().relax();
381  fvOptions.constrain(nuTildaEqn.ref());
382  solve(nuTildaEqn);
383  fvOptions.correct(nuTilda_);
384  bound(nuTilda_, dimensionedScalar("0", nuTilda_.dimensions(), 0.0));
385  nuTilda_.correctBoundaryConditions();
386 
387  correctNut(fv1);
388 }
389 
390 
391 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
392 
393 } // End namespace RASModels
394 } // End namespace Foam
395 
396 // ************************************************************************* //
surfaceScalarField & phi
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
U
Definition: pEqn.H:83
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedTensor skew(const dimensionedTensor &dt)
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1) const
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:100
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
bool readIfPresent(const dictionary &)
Update the value of dimensioned<Type> if found in the dictionary.
Generic dimensioned Type class.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
Finite-volume options.
Definition: fvOptions.H:52
tmp< volScalarField > chi() const
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:101
fv::options & fvOptions
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
static const wallDist & New(const fvMesh &mesh)
Definition: MeshObject.C:44
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
scalar y
const dimensionSet & dimensions() const
Return dimensions.
Dimension set for the base types.
Definition: dimensionSet.H:120
bool read(const char *, int32_t &)
Definition: int32IO.C:85
A class for handling words, derived from string.
Definition: word.H:59
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
virtual bool read()
Read RASProperties dictionary.
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Bound the given scalar field if it has gone unbounded.
const dimensionedVector & g
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/ubuntu/OpenFOAM-5.0/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:72
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:485
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
Spalart-Allmaras one-eqn mixing-length model for incompressible and compressible external flows...
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
#define WarningInFunction
Report a warning using Foam::Warning.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
tmp< volScalarField > fv1(const volScalarField &chi) const
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
void correctBoundaryConditions()
Correct boundary field.
dimensionedScalar pow6(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:103
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
tmp< volScalarField > fv2(const volScalarField &chi, const volScalarField &fv1) const
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
volScalarField & nu
tmp< volScalarField > fw(const volScalarField &Stilda) const
Namespace for OpenFOAM.
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:99