SpalartAllmarasDES.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) 2016-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 "SpalartAllmarasDES.H"
30 #include "fvOptions.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace LESModels
37 {
38 
39 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
40 
41 template<class BasicTurbulenceModel>
43 {
44  return nuTilda_/this->nu();
45 }
46 
47 
48 template<class BasicTurbulenceModel>
50 (
51  const volScalarField& chi
52 ) const
53 {
54  const volScalarField chi3("chi3", pow3(chi));
55  return chi3/(chi3 + pow3(Cv1_));
56 }
57 
58 
59 template<class BasicTurbulenceModel>
61 (
62  const volScalarField& chi,
63  const volScalarField& fv1
64 ) const
65 {
66  return 1.0 - chi/(1.0 + chi*fv1);
67 }
68 
69 
70 template<class BasicTurbulenceModel>
72 (
73  const volScalarField& chi
74 ) const
75 {
76  return Ct3_*exp(-Ct4_*sqr(chi));
77 }
78 
79 
80 template<class BasicTurbulenceModel>
82 (
83  const volTensorField& gradU
84 ) const
85 {
86  return sqrt(2.0)*mag(skew(gradU));
87 }
88 
89 
90 template<class BasicTurbulenceModel>
92 (
93  const volScalarField& chi,
94  const volScalarField& fv1,
95  const volScalarField& Omega,
96  const volScalarField& dTilda
97 ) const
98 {
99  return
100  (
101  max
102  (
103  Omega
104  + fv2(chi, fv1)*nuTilda_/sqr(kappa_*dTilda),
105  Cs_*Omega
106  )
107  );
108 }
109 
110 
111 template<class BasicTurbulenceModel>
113 (
114  const volScalarField& nur,
115  const volScalarField& Omega,
116  const volScalarField& dTilda
117 ) const
118 {
120  (
121  min
122  (
123  nur
124  /(
125  max
126  (
127  Omega,
128  dimensionedScalar("SMALL", Omega.dimensions(), SMALL)
129  )
130  *sqr(kappa_*dTilda)
131  ),
132  scalar(10)
133  )
134  );
135  tr.ref().boundaryFieldRef() == 0.0;
136 
137  return tr;
138 }
139 
140 
141 template<class BasicTurbulenceModel>
143 (
144  const volScalarField& Omega,
145  const volScalarField& dTilda
146 ) const
147 {
148  const volScalarField r(this->r(nuTilda_, Omega, dTilda));
149  const volScalarField g(r + Cw2_*(pow6(r) - r));
150 
151  return g*pow((1 + pow6(Cw3_))/(pow6(g) + pow6(Cw3_)), 1.0/6.0);
152 }
153 
154 
155 template<class BasicTurbulenceModel>
157 (
158  const volScalarField& chi,
159  const volScalarField& fv1
160 ) const
161 {
163  (
164  new volScalarField
165  (
166  IOobject
167  (
168  type() + ":psi",
169  this->time().timeName(),
170  this->mesh(),
173  ),
174  this->mesh(),
175  dimensionedScalar("one", dimless, 1)
176  )
177  );
178 
179  if (lowReCorrection_)
180  {
181  volScalarField& psi = tpsi.ref();
182 
183  const volScalarField fv2(this->fv2(chi, fv1));
184  const volScalarField ft2(this->ft2(chi));
185 
186  psi =
187  sqrt
188  (
189  min
190  (
191  scalar(100),
192  (1 - Cb1_/(Cw1_*sqr(kappa_)*fwStar_)*(ft2 + (1 - ft2)*fv2))
193  /max(SMALL, (fv1*max(scalar(1e-10), 1 - ft2)))
194  )
195  );
196  }
197 
198  return tpsi;
199 }
200 
201 
202 template<class BasicTurbulenceModel>
204 (
205  const volScalarField& chi,
206  const volScalarField& fv1,
207  const volTensorField& gradU
208 ) const
209 {
210  tmp<volScalarField> tdTilda(psi(chi, fv1)*CDES_*this->delta());
211  min(tdTilda.ref().ref(), tdTilda(), y_);
212  return tdTilda;
213 }
214 
215 
216 template<class BasicTurbulenceModel>
218 (
219  const volScalarField& fv1
220 )
221 {
222  this->nut_ = nuTilda_*fv1;
223  this->nut_.correctBoundaryConditions();
224  fv::options::New(this->mesh_).correct(this->nut_);
225 
226  BasicTurbulenceModel::correctNut();
227 }
228 
229 
230 template<class BasicTurbulenceModel>
232 {
233  correctNut(fv1(this->chi()));
234 }
235 
236 
237 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
238 
239 template<class BasicTurbulenceModel>
241 (
242  const alphaField& alpha,
243  const rhoField& rho,
244  const volVectorField& U,
245  const surfaceScalarField& alphaRhoPhi,
246  const surfaceScalarField& phi,
247  const transportModel& transport,
248  const word& propertiesName,
249  const word& type
250 )
251 :
253  (
254  type,
255  alpha,
256  rho,
257  U,
258  alphaRhoPhi,
259  phi,
260  transport,
261  propertiesName
262  ),
263 
264  sigmaNut_
265  (
267  (
268  "sigmaNut",
269  this->coeffDict_,
270  0.66666
271  )
272  ),
273  kappa_
274  (
276  (
277  "kappa",
278  this->coeffDict_,
279  0.41
280  )
281  ),
282  Cb1_
283  (
285  (
286  "Cb1",
287  this->coeffDict_,
288  0.1355
289  )
290  ),
291  Cb2_
292  (
294  (
295  "Cb2",
296  this->coeffDict_,
297  0.622
298  )
299  ),
300  Cw1_(Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_),
301  Cw2_
302  (
304  (
305  "Cw2",
306  this->coeffDict_,
307  0.3
308  )
309  ),
310  Cw3_
311  (
313  (
314  "Cw3",
315  this->coeffDict_,
316  2.0
317  )
318  ),
319  Cv1_
320  (
322  (
323  "Cv1",
324  this->coeffDict_,
325  7.1
326  )
327  ),
328  Cs_
329  (
331  (
332  "Cs",
333  this->coeffDict_,
334  0.3
335  )
336  ),
337  CDES_
338  (
340  (
341  "CDES",
342  this->coeffDict_,
343  0.65
344  )
345  ),
346  ck_
347  (
349  (
350  "ck",
351  this->coeffDict_,
352  0.07
353  )
354  ),
355  lowReCorrection_
356  (
358  (
359  "lowReCorrection",
360  this->coeffDict_,
361  true
362  )
363  ),
364  Ct3_
365  (
367  (
368  "Ct3",
369  this->coeffDict_,
370  1.2
371  )
372  ),
373  Ct4_
374  (
376  (
377  "Ct4",
378  this->coeffDict_,
379  0.5
380  )
381  ),
382  fwStar_
383  (
385  (
386  "fwStar",
387  this->coeffDict_,
388  0.424
389  )
390  ),
391 
392  nuTilda_
393  (
394  IOobject
395  (
396  "nuTilda",
397  this->runTime_.timeName(),
398  this->mesh_,
401  ),
402  this->mesh_
403  ),
404 
405  y_(wallDist::New(this->mesh_).y())
406 {
407  if (type == typeName)
408  {
409  this->printCoeffs(type);
410  }
411 }
412 
413 
414 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
415 
416 template<class BasicTurbulenceModel>
418 {
420  {
421  sigmaNut_.readIfPresent(this->coeffDict());
422  kappa_.readIfPresent(*this);
423 
424  Cb1_.readIfPresent(this->coeffDict());
425  Cb2_.readIfPresent(this->coeffDict());
426  Cw1_ = Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_;
427  Cw2_.readIfPresent(this->coeffDict());
428  Cw3_.readIfPresent(this->coeffDict());
429  Cv1_.readIfPresent(this->coeffDict());
430  Cs_.readIfPresent(this->coeffDict());
431 
432  CDES_.readIfPresent(this->coeffDict());
433  ck_.readIfPresent(this->coeffDict());
434 
435  lowReCorrection_.readIfPresent("lowReCorrection", this->coeffDict());
436  Ct3_.readIfPresent(this->coeffDict());
437  Ct4_.readIfPresent(this->coeffDict());
438  fwStar_.readIfPresent(this->coeffDict());
439 
440  return true;
441  }
442 
443  return false;
444 }
445 
446 
447 template<class BasicTurbulenceModel>
449 DnuTildaEff() const
450 {
451  return tmp<volScalarField>
452  (
453  new volScalarField("DnuTildaEff", (nuTilda_ + this->nu())/sigmaNut_)
454  );
455 }
456 
457 
458 template<class BasicTurbulenceModel>
460 {
461  const volScalarField chi(this->chi());
462  const volScalarField fv1(this->fv1(chi));
463 
464  tmp<volScalarField> tdTilda
465  (
466  new volScalarField
467  (
468  IOobject
469  (
470  "dTilda",
471  this->mesh_.time().timeName(),
472  this->mesh_,
475  ),
476  this->mesh_,
479  )
480  );
481  volScalarField& dTildaF = tdTilda.ref();
482  dTildaF = dTilda(chi, fv1, fvc::grad(this->U_));
483  dTildaF.correctBoundaryConditions();
484 
485  return sqr(this->nut()/ck_/dTildaF);
486 }
487 
488 
489 template<class BasicTurbulenceModel>
491 {
492  const volScalarField chi(this->chi());
493  const volScalarField fv1(this->fv1(chi));
494 
495  tmp<volScalarField> tLESRegion
496  (
497  new volScalarField
498  (
499  IOobject
500  (
501  "DES::LESRegion",
502  this->mesh_.time().timeName(),
503  this->mesh_,
506  ),
507  neg(dTilda(chi, fv1, fvc::grad(this->U_)) - y_)
508  )
509  );
510 
511  return tLESRegion;
512 }
513 
514 
515 template<class BasicTurbulenceModel>
517 {
518  if (!this->turbulence_)
519  {
520  return;
521  }
522 
523  // Local references
524  const alphaField& alpha = this->alpha_;
525  const rhoField& rho = this->rho_;
526  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
527  const volVectorField& U = this->U_;
528  fv::options& fvOptions(fv::options::New(this->mesh_));
529 
531 
532  const volScalarField chi(this->chi());
533  const volScalarField fv1(this->fv1(chi));
534  const volScalarField ft2(this->ft2(chi));
535 
536  tmp<volTensorField> tgradU = fvc::grad(U);
537  const volScalarField Omega(this->Omega(tgradU()));
538  const volScalarField dTilda(this->dTilda(chi, fv1, tgradU()));
539  const volScalarField Stilda(this->Stilda(chi, fv1, Omega, dTilda));
540 
541  tmp<fvScalarMatrix> nuTildaEqn
542  (
543  fvm::ddt(alpha, rho, nuTilda_)
544  + fvm::div(alphaRhoPhi, nuTilda_)
545  - fvm::laplacian(alpha*rho*DnuTildaEff(), nuTilda_)
546  - Cb2_/sigmaNut_*alpha*rho*magSqr(fvc::grad(nuTilda_))
547  ==
548  Cb1_*alpha*rho*Stilda*nuTilda_*(scalar(1) - ft2)
549  - fvm::Sp
550  (
551  (Cw1_*fw(Stilda, dTilda) - Cb1_/sqr(kappa_)*ft2)
552  *alpha*rho*nuTilda_/sqr(dTilda),
553  nuTilda_
554  )
555  + fvOptions(alpha, rho, nuTilda_)
556  );
557 
558  nuTildaEqn.ref().relax();
559  fvOptions.constrain(nuTildaEqn.ref());
560  solve(nuTildaEqn);
561  fvOptions.correct(nuTilda_);
562  bound(nuTilda_, dimensionedScalar(nuTilda_.dimensions(), Zero));
563  nuTilda_.correctBoundaryConditions();
564 
565  correctNut();
566 }
567 
568 
569 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
570 
571 } // End namespace LESModels
572 } // End namespace Foam
573 
574 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::LESModels::SpalartAllmarasDES::Stilda
tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1, const volScalarField &Omega, const volScalarField &dTilda) const
Definition: SpalartAllmarasDES.C:92
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
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:65
Foam::fv::optionList::correct
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
Definition: fvOptionListTemplates.C:355
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:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::LESModels::SpalartAllmarasDES::DnuTildaEff
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
Definition: SpalartAllmarasDES.C:449
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:103
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::LESModels::SpalartAllmarasDES::fv2
tmp< volScalarField > fv2(const volScalarField &chi, const volScalarField &fv1) const
Definition: SpalartAllmarasDES.C:61
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::LESModels::SpalartAllmarasDES
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
Definition: SpalartAllmarasDES.H:80
Foam::LESModels::SpalartAllmarasDES::r
tmp< volScalarField > r(const volScalarField &nur, const volScalarField &Omega, const volScalarField &dTilda) const
Definition: SpalartAllmarasDES.C:113
Foam::LESModels::DESModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: DESModel.H:74
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
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::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::LESModels::SpalartAllmarasDES::chi
tmp< volScalarField > chi() const
Definition: SpalartAllmarasDES.C:42
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::LESModels::SpalartAllmarasDES::LESRegion
virtual tmp< volScalarField > LESRegion() const
Return the LES field indicator.
Definition: SpalartAllmarasDES.C:490
Foam::pow6
dimensionedScalar pow6(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:122
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::LESModels::SpalartAllmarasDES::k
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition: SpalartAllmarasDES.C:459
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::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
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:42
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::LESModels::SpalartAllmarasDES::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: SpalartAllmarasDES.C:417
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
timeName
word timeName
Definition: getTimeIndex.H:3
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
SpalartAllmarasDES.H
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:26
Foam::LESModels::DESModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: DESModel.H:75
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:940
U
U
Definition: pEqn.H:72
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Foam::LESModels::SpalartAllmarasDES::psi
tmp< volScalarField > psi(const volScalarField &chi, const volScalarField &fv1) const
Definition: SpalartAllmarasDES.C:157
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
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::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::LESModels::SpalartAllmarasDES::ft2
tmp< volScalarField > ft2(const volScalarField &chi) const
Definition: SpalartAllmarasDES.C:72
Foam::LESModels::SpalartAllmarasDES::dTilda
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
Length scale.
Definition: SpalartAllmarasDES.C:204
Foam::eddyViscosity< LESModel< BasicTurbulenceModel > >::correct
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
Definition: eddyViscosity.C:134
Foam::tr
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:51
nut
scalar nut
Definition: evaluateNearWall.H:5
Foam::LESModels::DESModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: DESModel.H:73
Foam::GeometricField< scalar, fvPatchField, volMesh >
psi
const volScalarField & psi
Definition: createFieldRefs.H:1
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
Foam::LESModels::SpalartAllmarasDES::correct
virtual void correct()
Correct nuTilda and related properties.
Definition: SpalartAllmarasDES.C:516
Foam::LESModels::SpalartAllmarasDES::fv1
tmp< volScalarField > fv1(const volScalarField &chi) const
Definition: SpalartAllmarasDES.C:50
Foam::neg
dimensionedScalar neg(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:199
Foam::LESModels::SpalartAllmarasDES::fw
tmp< volScalarField > fw(const volScalarField &Omega, const volScalarField &dTilda) const
Definition: SpalartAllmarasDES.C:143
Foam::LESModels::SpalartAllmarasDES::Omega
tmp< volScalarField > Omega(const volTensorField &gradU) const
Definition: SpalartAllmarasDES.C:82
Foam::LESModels::SpalartAllmarasDES::correctNut
virtual void correctNut()
Definition: SpalartAllmarasDES.C:231
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::LESModels::DESModel
Templated abstract base class for DES models.
Definition: DESModel.H:54
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::Switch::getOrAddToDict
static Switch getOrAddToDict(const word &key, dictionary &dict, const Switch deflt=switchType::FALSE)
Definition: Switch.C:164
Foam::IOobject::MUST_READ
Definition: IOobject.H:185