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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "SpalartAllmaras.H"
30 #include "fvOptions.H"
31 #include "bound.H"
32 #include "wallDist.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace RASModels
39 {
40 
41 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
42 
43 template<class BasicTurbulenceModel>
45 {
46  return nuTilda_/this->nu();
47 }
48 
49 
50 template<class BasicTurbulenceModel>
52 (
53  const volScalarField& chi
54 ) const
55 {
56  const volScalarField chi3(pow3(chi));
57  return chi3/(chi3 + pow3(Cv1_));
58 }
59 
60 
61 template<class BasicTurbulenceModel>
63 (
64  const volScalarField& chi,
65  const volScalarField& fv1
66 ) const
67 {
68  return 1.0 - chi/(1.0 + chi*fv1);
69 }
70 
71 
72 template<class BasicTurbulenceModel>
74 (
75  const volScalarField& chi,
76  const volScalarField& fv1
77 ) const
78 {
79  volScalarField Omega(::sqrt(2.0)*mag(skew(fvc::grad(this->U_))));
80 
81  return
82  (
83  max
84  (
85  Omega
86  + fv2(chi, fv1)*nuTilda_/sqr(kappa_*y_),
87  Cs_*Omega
88  )
89  );
90 }
91 
92 
93 template<class BasicTurbulenceModel>
95 (
96  const volScalarField& Stilda
97 ) const
98 {
100  (
101  min
102  (
103  nuTilda_
104  /(
105  max
106  (
107  Stilda,
108  dimensionedScalar("SMALL", Stilda.dimensions(), SMALL)
109  )
110  *sqr(kappa_*y_)
111  ),
112  scalar(10)
113  )
114  );
115  r.boundaryFieldRef() == 0.0;
116 
117  const volScalarField g(r + Cw2_*(pow6(r) - r));
118 
119  return g*pow((1.0 + pow6(Cw3_))/(pow6(g) + pow6(Cw3_)), 1.0/6.0);
120 }
121 
122 
123 template<class BasicTurbulenceModel>
125 (
126  const volScalarField& fv1
127 )
128 {
129  this->nut_ = nuTilda_*fv1;
130  this->nut_.correctBoundaryConditions();
131  fv::options::New(this->mesh_).correct(this->nut_);
132 
133  BasicTurbulenceModel::correctNut();
134 }
135 
136 
137 template<class BasicTurbulenceModel>
139 {
140  correctNut(fv1(this->chi()));
141 }
142 
143 
144 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
145 
146 template<class BasicTurbulenceModel>
148 (
149  const alphaField& alpha,
150  const rhoField& rho,
151  const volVectorField& U,
152  const surfaceScalarField& alphaRhoPhi,
153  const surfaceScalarField& phi,
154  const transportModel& transport,
155  const word& propertiesName,
156  const word& type
157 )
158 :
160  (
161  type,
162  alpha,
163  rho,
164  U,
165  alphaRhoPhi,
166  phi,
167  transport,
168  propertiesName
169  ),
170 
171  sigmaNut_
172  (
174  (
175  "sigmaNut",
176  this->coeffDict_,
177  0.66666
178  )
179  ),
180  kappa_
181  (
183  (
184  "kappa",
185  this->coeffDict_,
186  0.41
187  )
188  ),
189 
190  Cb1_
191  (
193  (
194  "Cb1",
195  this->coeffDict_,
196  0.1355
197  )
198  ),
199  Cb2_
200  (
202  (
203  "Cb2",
204  this->coeffDict_,
205  0.622
206  )
207  ),
208  Cw1_(Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_),
209  Cw2_
210  (
212  (
213  "Cw2",
214  this->coeffDict_,
215  0.3
216  )
217  ),
218  Cw3_
219  (
221  (
222  "Cw3",
223  this->coeffDict_,
224  2.0
225  )
226  ),
227  Cv1_
228  (
230  (
231  "Cv1",
232  this->coeffDict_,
233  7.1
234  )
235  ),
236  Cs_
237  (
239  (
240  "Cs",
241  this->coeffDict_,
242  0.3
243  )
244  ),
245 
246  nuTilda_
247  (
248  IOobject
249  (
250  "nuTilda",
251  this->runTime_.timeName(),
252  this->mesh_,
255  ),
256  this->mesh_
257  ),
258 
259  y_(wallDist::New(this->mesh_).y())
260 {
261  if (type == typeName)
262  {
263  this->printCoeffs(type);
264  }
265 }
266 
267 
268 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
269 
270 template<class BasicTurbulenceModel>
272 {
274  {
275  sigmaNut_.readIfPresent(this->coeffDict());
276  kappa_.readIfPresent(this->coeffDict());
277 
278  Cb1_.readIfPresent(this->coeffDict());
279  Cb2_.readIfPresent(this->coeffDict());
280  Cw1_ = Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_;
281  Cw2_.readIfPresent(this->coeffDict());
282  Cw3_.readIfPresent(this->coeffDict());
283  Cv1_.readIfPresent(this->coeffDict());
284  Cs_.readIfPresent(this->coeffDict());
285 
286  return true;
287  }
288 
289  return false;
290 }
291 
292 
293 template<class BasicTurbulenceModel>
295 {
296  return tmp<volScalarField>
297  (
298  new volScalarField("DnuTildaEff", (nuTilda_ + this->nu())/sigmaNut_)
299  );
300 }
301 
302 
303 template<class BasicTurbulenceModel>
305 {
307  << "Turbulence kinetic energy not defined for "
308  << "Spalart-Allmaras model. Returning zero field"
309  << endl;
310 
312  (
313  IOobject
314  (
315  "k",
316  this->runTime_.timeName(),
317  this->mesh_
318  ),
319  this->mesh_,
322  );
323 }
324 
325 
326 template<class BasicTurbulenceModel>
328 {
330  << "Turbulence kinetic energy dissipation rate not defined for "
331  << "Spalart-Allmaras model. Returning zero field"
332  << endl;
333 
335  (
336  IOobject
337  (
338  "epsilon",
339  this->runTime_.timeName(),
340  this->mesh_
341  ),
342  this->mesh_,
345  );
346 }
347 
348 
349 template<class BasicTurbulenceModel>
351 {
353  << "Specific dissipation rate not defined for "
354  << "Spalart-Allmaras model. Returning zero field"
355  << endl;
356 
358  (
359  IOobject
360  (
361  IOobject::groupName("omega", this->alphaRhoPhi_.group()),
362  this->runTime_.timeName(),
363  this->mesh_
364  ),
365  this->mesh_,
367  );
368 }
369 
370 
371 template<class BasicTurbulenceModel>
373 {
374  if (!this->turbulence_)
375  {
376  return;
377  }
378 
379  {
380  // Local references
381  const alphaField& alpha = this->alpha_;
382  const rhoField& rho = this->rho_;
383  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
384  fv::options& fvOptions(fv::options::New(this->mesh_));
385 
387 
388  const volScalarField chi(this->chi());
389  const volScalarField fv1(this->fv1(chi));
390 
391  const volScalarField Stilda(this->Stilda(chi, fv1));
392 
393  tmp<fvScalarMatrix> nuTildaEqn
394  (
395  fvm::ddt(alpha, rho, nuTilda_)
396  + fvm::div(alphaRhoPhi, nuTilda_)
397  - fvm::laplacian(alpha*rho*DnuTildaEff(), nuTilda_)
398  - Cb2_/sigmaNut_*alpha*rho*magSqr(fvc::grad(nuTilda_))
399  ==
400  Cb1_*alpha*rho*Stilda*nuTilda_
401  - fvm::Sp(Cw1_*alpha*rho*fw(Stilda)*nuTilda_/sqr(y_), nuTilda_)
402  + fvOptions(alpha, rho, nuTilda_)
403  );
404 
405  nuTildaEqn.ref().relax();
406  fvOptions.constrain(nuTildaEqn.ref());
407  solve(nuTildaEqn);
408  fvOptions.correct(nuTilda_);
409  bound(nuTilda_, dimensionedScalar(nuTilda_.dimensions(), Zero));
410  nuTilda_.correctBoundaryConditions();
411  }
412 
413  // Update nut with latest available k,epsilon
414  correctNut();
415 }
416 
417 
418 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
419 
420 } // End namespace RASModels
421 } // End namespace Foam
422 
423 // ************************************************************************* //
wallDist.H
Foam::RASModels::SpalartAllmaras::DnuTildaEff
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
Definition: SpalartAllmaras.C:294
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
fvOptions.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::RASModels::SpalartAllmaras::correctNut
virtual void correctNut()
Definition: SpalartAllmaras.C:138
Foam::fv::optionList::correct
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
Definition: fvOptionListTemplates.C:323
Foam::skew
dimensionedTensor skew(const dimensionedTensor &dt)
Definition: dimensionedTensor.C:138
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
SpalartAllmaras.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::RASModels::SpalartAllmaras::fv1
tmp< volScalarField > fv1(const volScalarField &chi) const
Definition: SpalartAllmaras.C:52
Foam::RASModel
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:51
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:104
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:35
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::MeshObject< fvMesh, UpdateableMeshObject, wallDist >::New
static const wallDist & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::RASModels::SpalartAllmaras::k
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: SpalartAllmaras.C:304
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
rho
rho
Definition: readInitialConditions.H:88
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::RASModels::SpalartAllmaras::omega
virtual tmp< volScalarField > omega() const
Return the (estimated) specific dissipation rate.
Definition: SpalartAllmaras.C:350
Foam::zeroGradientFvPatchField
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
Definition: zeroGradientFvPatchField.H:64
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::RASModels::SpalartAllmaras::fv2
tmp< volScalarField > fv2(const volScalarField &chi, const volScalarField &fv1) const
Definition: SpalartAllmaras.C:63
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::RASModels::SpalartAllmaras::Stilda
tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1) const
Definition: SpalartAllmaras.C:74
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::RASModels::SpalartAllmaras::fw
tmp< volScalarField > fw(const volScalarField &Stilda) const
Definition: SpalartAllmaras.C:95
Foam::pow6
dimensionedScalar pow6(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:122
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::dimensioned::readIfPresent
bool readIfPresent(const dictionary &dict)
Definition: dimensionedType.C:483
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::RASModels::SpalartAllmaras::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: SpalartAllmaras.C:271
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
bound.H
Bound the given scalar field if it has gone unbounded.
correct
fvOptions correct(rho)
Foam::solve
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::RASModels::SpalartAllmaras::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: SpalartAllmaras.C:372
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::fv::options
Finite-volume options.
Definition: fvOptions.H:55
Foam::RASModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:97
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:43
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:24
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::RASModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:96
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:940
Foam::RASModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:98
U
U
Definition: pEqn.H:72
Foam::RASModels::SpalartAllmaras
Spalart-Allmaras one-eqn mixing-length model for incompressible and compressible external flows.
Definition: SpalartAllmaras.H:91
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::RASModels::SpalartAllmaras::chi
tmp< volScalarField > chi() const
Definition: SpalartAllmaras.C:44
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::eddyViscosity
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:55
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:303
Foam::RASModels::SpalartAllmaras::epsilon
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: SpalartAllmaras.C:327
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::IOobject::MUST_READ
Definition: IOobject.H:120