SpalartAllmaras.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) 2011-2023 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 "fvModels.H"
28 #include "fvConstraints.H"
29 #include "bound.H"
30 #include "wallDist.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace RASModels
37 {
38 
39 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
40 
41 template<class BasicMomentumTransportModel>
43 {
44  return volScalarField::New(typedName("chi"), nuTilda_/this->nu());
45 }
46 
47 
48 template<class BasicMomentumTransportModel>
50 (
51  const volScalarField& chi
52 ) const
53 {
54  const volScalarField chi3(typedName("chi3"), pow3(chi));
55  return volScalarField::New(typedName("fv1"), chi3/(chi3 + pow3(Cv1_)));
56 }
57 
58 
59 template<class BasicMomentumTransportModel>
61 (
62  const volScalarField::Internal& chi,
63  const volScalarField::Internal& fv1
64 ) const
65 {
67  (
68  typedName("fv2"),
69  1.0 - chi/(1.0 + chi*fv1)
70  );
71 }
72 
73 
74 template<class BasicMomentumTransportModel>
77 (
78  const volScalarField::Internal& chi,
79  const volScalarField::Internal& fv1
80 ) const
81 {
82  const volScalarField::Internal Omega
83  (
84  typedName("Omega"),
85  ::sqrt(2.0)*mag(skew(fvc::grad(this->U_)().v()))
86  );
87 
89  (
90  typedName("Stilda"),
91  (
92  max
93  (
94  Omega
95  + fv2(chi, fv1)*nuTilda_/sqr(kappa_*this->y()()),
96  Cs_*Omega
97  )
98  )
99  );
100 }
101 
102 
103 template<class BasicMomentumTransportModel>
105 (
106  const volScalarField::Internal& Stilda
107 ) const
108 {
110  (
111  typedName("r"),
112  min
113  (
114  nuTilda_()
115  /(
116  max
117  (
118  Stilda,
119  dimensionedScalar(Stilda.dimensions(), small)
120  )
121  *sqr(kappa_*this->y()())
122  ),
123  scalar(10.0)
124  )
125  );
126 
127  const volScalarField::Internal g(typedName("g"), r + Cw2_*(pow6(r) - r));
128 
130  (
131  typedName("fw"),
132  g*pow((1.0 + pow6(Cw3_))/(pow6(g) + pow6(Cw3_)), 1.0/6.0)
133  );
134 }
135 
136 
137 template<class BasicMomentumTransportModel>
139 (
140  const volScalarField& fv1
141 )
142 {
143  this->nut_ = nuTilda_*fv1;
144  this->nut_.correctBoundaryConditions();
145  fvConstraints::New(this->mesh_).constrain(this->nut_);
146 }
147 
148 
149 template<class BasicMomentumTransportModel>
151 {
152  correctNut(fv1(this->chi()));
153 }
154 
155 
156 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
157 
158 template<class BasicMomentumTransportModel>
160 (
161  const alphaField& alpha,
162  const rhoField& rho,
163  const volVectorField& U,
164  const surfaceScalarField& alphaRhoPhi,
165  const surfaceScalarField& phi,
166  const viscosity& viscosity,
167  const word& type
168 )
169 :
170  eddyViscosity<RASModel<BasicMomentumTransportModel>>
171  (
172  type,
173  alpha,
174  rho,
175  U,
176  alphaRhoPhi,
177  phi,
178  viscosity
179  ),
180 
181  sigmaNut_
182  (
183  dimensioned<scalar>::lookupOrAddToDict
184  (
185  "sigmaNut",
186  this->coeffDict_,
187  0.66666
188  )
189  ),
190  kappa_
191  (
192  dimensioned<scalar>::lookupOrAddToDict
193  (
194  "kappa",
195  this->coeffDict_,
196  0.41
197  )
198  ),
199 
200  Cb1_
201  (
202  dimensioned<scalar>::lookupOrAddToDict
203  (
204  "Cb1",
205  this->coeffDict_,
206  0.1355
207  )
208  ),
209  Cb2_
210  (
211  dimensioned<scalar>::lookupOrAddToDict
212  (
213  "Cb2",
214  this->coeffDict_,
215  0.622
216  )
217  ),
218  Cw1_(Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_),
219  Cw2_
220  (
221  dimensioned<scalar>::lookupOrAddToDict
222  (
223  "Cw2",
224  this->coeffDict_,
225  0.3
226  )
227  ),
228  Cw3_
229  (
230  dimensioned<scalar>::lookupOrAddToDict
231  (
232  "Cw3",
233  this->coeffDict_,
234  2.0
235  )
236  ),
237  Cv1_
238  (
239  dimensioned<scalar>::lookupOrAddToDict
240  (
241  "Cv1",
242  this->coeffDict_,
243  7.1
244  )
245  ),
246  Cs_
247  (
248  dimensioned<scalar>::lookupOrAddToDict
249  (
250  "Cs",
251  this->coeffDict_,
252  0.3
253  )
254  ),
255 
256  nuTilda_
257  (
258  IOobject
259  (
260  "nuTilda",
261  this->runTime_.name(),
262  this->mesh_,
263  IOobject::MUST_READ,
264  IOobject::AUTO_WRITE
265  ),
266  this->mesh_
267  )
268 {
269  if (type == typeName)
270  {
271  this->printCoeffs(type);
272  }
273 }
274 
275 
276 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
277 
278 template<class BasicMomentumTransportModel>
280 {
282  {
283  sigmaNut_.readIfPresent(this->coeffDict());
284  kappa_.readIfPresent(this->coeffDict());
285 
286  Cb1_.readIfPresent(this->coeffDict());
287  Cb2_.readIfPresent(this->coeffDict());
288  Cw1_ = Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_;
289  Cw2_.readIfPresent(this->coeffDict());
290  Cw3_.readIfPresent(this->coeffDict());
291  Cv1_.readIfPresent(this->coeffDict());
292  Cs_.readIfPresent(this->coeffDict());
293 
294  return true;
295  }
296  else
297  {
298  return false;
299  }
300 }
301 
302 
303 template<class BasicMomentumTransportModel>
306 {
307  return volScalarField::New
308  (
309  "DnuTildaEff",
310  (nuTilda_ + this->nu())/sigmaNut_
311  );
312 }
313 
314 
315 template<class BasicMomentumTransportModel>
317 {
318  return volScalarField::New
319  (
320  "k",
321  this->mesh_,
322  dimensionedScalar(dimensionSet(0, 2, -2, 0, 0), 0)
323  );
324 }
325 
326 
327 template<class BasicMomentumTransportModel>
330 {
332  << "Turbulence kinetic energy dissipation rate not defined for "
333  << "Spalart-Allmaras model. Returning zero field"
334  << endl;
335 
336  return volScalarField::New
337  (
338  "epsilon",
339  this->mesh_,
340  dimensionedScalar(dimensionSet(0, 2, -3, 0, 0), 0)
341  );
342 }
343 
344 
345 template<class BasicMomentumTransportModel>
348 {
350  << "Turbulence specific dissipation rate not defined for "
351  << "Spalart-Allmaras model. Returning zero field"
352  << endl;
353 
354  return volScalarField::New
355  (
356  "omega",
357  this->mesh_,
359  );
360 }
361 
362 
363 template<class BasicMomentumTransportModel>
365 {
366  if (!this->turbulence_)
367  {
368  return;
369  }
370 
371  // Local references
372  const alphaField& alpha = this->alpha_;
373  const rhoField& rho = this->rho_;
374  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
375  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
377  (
378  Foam::fvConstraints::New(this->mesh_)
379  );
380 
382 
383  const volScalarField chi(this->chi());
384  const volScalarField fv1(this->fv1(chi));
385 
386  const volScalarField::Internal Stilda(this->Stilda(chi, fv1));
387 
388  tmp<fvScalarMatrix> nuTildaEqn
389  (
390  fvm::ddt(alpha, rho, nuTilda_)
391  + fvm::div(alphaRhoPhi, nuTilda_)
392  - fvm::laplacian(alpha*rho*DnuTildaEff(), nuTilda_)
393  - Cb2_/sigmaNut_*alpha*rho*magSqr(fvc::grad(nuTilda_))
394  ==
395  Cb1_*alpha()*rho()*Stilda*nuTilda_()
396  - fvm::Sp
397  (
398  Cw1_*alpha()*rho()*fw(Stilda)*nuTilda_()/sqr(this->y()()),
399  nuTilda_
400  )
401  + fvModels.source(alpha, rho, nuTilda_)
402  );
403 
404  nuTildaEqn.ref().relax();
405  fvConstraints.constrain(nuTildaEqn.ref());
406  solve(nuTildaEqn);
407  fvConstraints.constrain(nuTilda_);
408  bound(nuTilda_, dimensionedScalar(nuTilda_.dimensions(), 0));
409  nuTilda_.correctBoundaryConditions();
410 
411  correctNut(fv1);
412 }
413 
414 
415 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
416 
417 } // End namespace RASModels
418 } // End namespace Foam
419 
420 // ************************************************************************* //
scalar y
Bound the given scalar field where it is below the specified minimum.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &, const Field< Type > &)
Return a temporary field constructed from name, mesh,.
Generic GeometricField class.
void correctBoundaryConditions()
Correct boundary field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
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
tmp< volScalarField > chi() const
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
tmp< volScalarField > fv1(const volScalarField &chi) const
tmp< volScalarField::Internal > fv2(const volScalarField::Internal &chi, const volScalarField::Internal &fv1) const
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
tmp< volScalarField::Internal > fw(const volScalarField::Internal &Stilda) const
SpalartAllmaras(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.
tmp< volScalarField::Internal > Stilda(const volScalarField::Internal &chi, const volScalarField::Internal &fv1) const
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
virtual bool read()
Read RASProperties dictionary.
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:100
Dimension set for the base types.
Definition: dimensionSet.H:125
Generic dimensioned Type class.
bool readIfPresent(const dictionary &)
Update the value of dimensioned<Type> if found in the dictionary.
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:181
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))
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
#define WarningInFunction
Report a warning using Foam::Warning.
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
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 > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
dimensionedScalar pow6(const dimensionedScalar &ds)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
dimensionedScalar pow3(const dimensionedScalar &ds)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionSet dimTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:176
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
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dimensionedTensor skew(const dimensionedTensor &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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.