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