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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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 "SpalartAllmarasDES.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 namespace LESModels
33 {
34 
35 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
36 
37 template<class BasicTurbulenceModel>
39 {
40  return nuTilda_/this->nu();
41 }
42 
43 
44 template<class BasicTurbulenceModel>
46 (
47  const volScalarField& chi
48 ) const
49 {
50  const volScalarField chi3("chi3", pow3(chi));
51  return chi3/(chi3 + pow3(Cv1_));
52 }
53 
54 
55 template<class BasicTurbulenceModel>
57 (
58  const volScalarField& chi,
59  const volScalarField& fv1
60 ) const
61 {
62  return 1.0 - chi/(1.0 + chi*fv1);
63 }
64 
65 
66 template<class BasicTurbulenceModel>
68 (
69  const volScalarField& chi
70 ) const
71 {
72  return Ct3_*exp(-Ct4_*sqr(chi));
73 }
74 
75 
76 template<class BasicTurbulenceModel>
78 (
79  const volTensorField& gradU
80 ) const
81 {
82  return sqrt(2.0)*mag(skew(gradU));
83 }
84 
85 
86 template<class BasicTurbulenceModel>
88 (
89  const volScalarField& chi,
90  const volScalarField& fv1,
91  const volScalarField& Omega,
92  const volScalarField& dTilda
93 ) const
94 {
95  return
96  (
97  max
98  (
99  Omega
100  + fv2(chi, fv1)*nuTilda_/sqr(kappa_*dTilda),
101  Cs_*Omega
102  )
103  );
104 }
105 
106 
107 template<class BasicTurbulenceModel>
109 (
110  const volScalarField& nur,
111  const volScalarField& Omega,
112  const volScalarField& dTilda
113 ) const
114 {
116  (
117  min
118  (
119  nur
120  /(
121  max
122  (
123  Omega,
124  dimensionedScalar("SMALL", Omega.dimensions(), SMALL)
125  )
126  *sqr(kappa_*dTilda)
127  ),
128  scalar(10)
129  )
130  );
131  tr().boundaryField() == 0.0;
132 
133  return tr;
134 }
135 
136 
137 template<class BasicTurbulenceModel>
139 (
140  const volScalarField& Omega,
141  const volScalarField& dTilda
142 ) const
143 {
144  const volScalarField r(this->r(nuTilda_, Omega, dTilda));
145  const volScalarField g(r + Cw2_*(pow6(r) - r));
146 
147  return g*pow((1 + pow6(Cw3_))/(pow6(g) + pow6(Cw3_)), 1.0/6.0);
148 }
149 
150 
151 template<class BasicTurbulenceModel>
153 (
154  const volScalarField& chi,
155  const volScalarField& fv1
156 ) const
157 {
159  (
160  new volScalarField
161  (
162  IOobject
163  (
164  type() + ":psi",
165  this->time().timeName(),
166  this->mesh(),
169  ),
170  this->mesh(),
171  dimensionedScalar("one", dimless, 1)
172  )
173  );
174 
175  if (lowReCorrection_)
176  {
177  volScalarField& psi = tpsi();
178 
179  const volScalarField fv2(this->fv2(chi, fv1));
180  const volScalarField ft2(this->ft2(chi));
181 
182  psi =
183  sqrt
184  (
185  min
186  (
187  scalar(100),
188  (1 - Cb1_/(Cw1_*sqr(kappa_)*fwStar_)*(ft2 + (1 - ft2)*fv2))
189  /max(SMALL, (fv1*max(scalar(1e-10), 1 - ft2)))
190  )
191  );
192  }
193 
194  return tpsi;
195 }
196 
197 
198 template<class BasicTurbulenceModel>
200 (
201  const volScalarField& chi,
202  const volScalarField& fv1,
203  const volTensorField& gradU
204 ) const
205 {
206  tmp<volScalarField> tdTilda(psi(chi, fv1)*CDES_*this->delta());
207  min(tdTilda().dimensionedInternalField(), tdTilda(), y_);
208  return tdTilda;
209 }
210 
211 
212 template<class BasicTurbulenceModel>
214 (
215  const volScalarField& fv1
216 )
217 {
218  this->nut_ = nuTilda_*fv1;
219  this->nut_.correctBoundaryConditions();
220 
221  BasicTurbulenceModel::correctNut();
222 }
223 
224 
225 template<class BasicTurbulenceModel>
227 {
228  correctNut(fv1(this->chi()));
229 }
230 
231 
232 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
233 
234 template<class BasicTurbulenceModel>
236 (
237  const alphaField& alpha,
238  const rhoField& rho,
239  const volVectorField& U,
240  const surfaceScalarField& alphaRhoPhi,
241  const surfaceScalarField& phi,
242  const transportModel& transport,
243  const word& propertiesName,
244  const word& type
245 )
246 :
248  (
249  type,
250  alpha,
251  rho,
252  U,
253  alphaRhoPhi,
254  phi,
255  transport,
256  propertiesName
257  ),
258 
259  sigmaNut_
260  (
262  (
263  "sigmaNut",
264  this->coeffDict_,
265  0.66666
266  )
267  ),
268  kappa_
269  (
271  (
272  "kappa",
273  this->coeffDict_,
274  0.41
275  )
276  ),
277  Cb1_
278  (
280  (
281  "Cb1",
282  this->coeffDict_,
283  0.1355
284  )
285  ),
286  Cb2_
287  (
289  (
290  "Cb2",
291  this->coeffDict_,
292  0.622
293  )
294  ),
295  Cw1_(Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_),
296  Cw2_
297  (
299  (
300  "Cw2",
301  this->coeffDict_,
302  0.3
303  )
304  ),
305  Cw3_
306  (
308  (
309  "Cw3",
310  this->coeffDict_,
311  2.0
312  )
313  ),
314  Cv1_
315  (
317  (
318  "Cv1",
319  this->coeffDict_,
320  7.1
321  )
322  ),
323  Cs_
324  (
326  (
327  "Cs",
328  this->coeffDict_,
329  0.3
330  )
331  ),
332  CDES_
333  (
335  (
336  "CDES",
337  this->coeffDict_,
338  0.65
339  )
340  ),
341  ck_
342  (
344  (
345  "ck",
346  this->coeffDict_,
347  0.07
348  )
349  ),
350  lowReCorrection_
351  (
353  (
354  "lowReCorrection",
355  this->coeffDict_,
356  true
357  )
358  ),
359  Ct3_
360  (
362  (
363  "Ct3",
364  this->coeffDict_,
365  1.2
366  )
367  ),
368  Ct4_
369  (
371  (
372  "Ct4",
373  this->coeffDict_,
374  0.5
375  )
376  ),
377  fwStar_
378  (
380  (
381  "fwStar",
382  this->coeffDict_,
383  0.424
384  )
385  ),
386 
387  nuTilda_
388  (
389  IOobject
390  (
391  "nuTilda",
392  this->runTime_.timeName(),
393  this->mesh_,
396  ),
397  this->mesh_
398  ),
399 
400  y_(wallDist::New(this->mesh_).y())
401 {
402  if (type == typeName)
403  {
404  this->printCoeffs(type);
405  }
406 }
407 
408 
409 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
410 
411 template<class BasicTurbulenceModel>
413 {
415  {
416  sigmaNut_.readIfPresent(this->coeffDict());
417  kappa_.readIfPresent(*this);
418 
419  Cb1_.readIfPresent(this->coeffDict());
420  Cb2_.readIfPresent(this->coeffDict());
421  Cw1_ = Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_;
422  Cw2_.readIfPresent(this->coeffDict());
423  Cw3_.readIfPresent(this->coeffDict());
424  Cv1_.readIfPresent(this->coeffDict());
425  Cs_.readIfPresent(this->coeffDict());
426 
427  CDES_.readIfPresent(this->coeffDict());
428  ck_.readIfPresent(this->coeffDict());
429 
430  lowReCorrection_.readIfPresent("lowReCorrection", this->coeffDict());
431  Ct3_.readIfPresent(this->coeffDict());
432  Ct4_.readIfPresent(this->coeffDict());
433  fwStar_.readIfPresent(this->coeffDict());
434 
435  return true;
436  }
437  else
438  {
439  return false;
440  }
441 }
442 
443 
444 template<class BasicTurbulenceModel>
446 DnuTildaEff() const
447 {
448  return tmp<volScalarField>
449  (
450  new volScalarField("DnuTildaEff", (nuTilda_ + this->nu())/sigmaNut_)
451  );
452 }
453 
454 
455 template<class BasicTurbulenceModel>
457 {
458  const volScalarField chi(this->chi());
459  const volScalarField fv1(this->fv1(chi));
460  return sqr(this->nut()/ck_/dTilda(chi, fv1, fvc::grad(this->U_)));
461 }
462 
463 
464 template<class BasicTurbulenceModel>
466 {
467  const volScalarField chi(this->chi());
468  const volScalarField fv1(this->fv1(chi));
469 
470  tmp<volScalarField> tLESRegion
471  (
472  new volScalarField
473  (
474  IOobject
475  (
476  "DES::LESRegion",
477  this->mesh_.time().timeName(),
478  this->mesh_,
481  ),
482  neg(dTilda(chi, fv1, fvc::grad(this->U_)) - y_)
483  )
484  );
485 
486  return tLESRegion;
487 }
488 
489 
490 template<class BasicTurbulenceModel>
492 {
493  if (!this->turbulence_)
494  {
495  return;
496  }
497 
498  // Local references
499  const alphaField& alpha = this->alpha_;
500  const rhoField& rho = this->rho_;
501  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
502  const volVectorField& U = this->U_;
503 
505 
506  const volScalarField chi(this->chi());
507  const volScalarField fv1(this->fv1(chi));
508  const volScalarField ft2(this->ft2(chi));
509 
510  tmp<volTensorField> tgradU = fvc::grad(U);
511  const volScalarField Omega(this->Omega(tgradU()));
512  const volScalarField dTilda(this->dTilda(chi, fv1, tgradU()));
513  const volScalarField Stilda(this->Stilda(chi, fv1, Omega, dTilda));
514 
515  tmp<fvScalarMatrix> nuTildaEqn
516  (
517  fvm::ddt(alpha, rho, nuTilda_)
518  + fvm::div(alphaRhoPhi, nuTilda_)
519  - fvm::laplacian(alpha*rho*DnuTildaEff(), nuTilda_)
520  - Cb2_/sigmaNut_*alpha*rho*magSqr(fvc::grad(nuTilda_))
521  ==
522  Cb1_*alpha*rho*Stilda*nuTilda_*(scalar(1) - ft2)
523  - fvm::Sp
524  (
525  (Cw1_*fw(Stilda, dTilda) - Cb1_/sqr(kappa_)*ft2)
526  *alpha*rho*nuTilda_/sqr(dTilda),
527  nuTilda_
528  )
529  );
530 
531  nuTildaEqn().relax();
532  solve(nuTildaEqn);
533  bound(nuTilda_, dimensionedScalar("zero", nuTilda_.dimensions(), 0.0));
534  nuTilda_.correctBoundaryConditions();
535 
536  correctNut(fv1);
537 }
538 
539 
540 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
541 
542 } // End namespace LESModels
543 } // End namespace Foam
544 
545 // ************************************************************************* //
Foam::LESModels::SpalartAllmarasDES::Stilda
tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1, const volScalarField &Omega, const volScalarField &dTilda) const
Definition: SpalartAllmarasDES.C:88
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
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::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::LESModels::SpalartAllmarasDES::DnuTildaEff
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
Definition: SpalartAllmarasDES.C:446
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::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
Foam::LESModels::SpalartAllmarasDES::fv2
tmp< volScalarField > fv2(const volScalarField &chi, const volScalarField &fv1) const
Definition: SpalartAllmarasDES.C:57
dimensionedInternalField
rDeltaT dimensionedInternalField()
Foam::LESModels::SpalartAllmarasDES::r
tmp< volScalarField > r(const volScalarField &nur, const volScalarField &Omega, const volScalarField &dTilda) const
Definition: SpalartAllmarasDES.C:109
Foam::LESModels::DESModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: DESModel.H:71
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
g
const dimensionedVector & g
Definition: setRegionFluidFields.H:33
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
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::LESModels::SpalartAllmarasDES::chi
tmp< volScalarField > chi() const
Definition: SpalartAllmarasDES.C:38
nu
volScalarField & nu
Definition: readMechanicalProperties.H:179
Foam::LESModels::SpalartAllmarasDES::LESRegion
virtual tmp< volScalarField > LESRegion() const
Return the LES field indicator.
Definition: SpalartAllmarasDES.C:465
Foam::pow6
dimensionedScalar pow6(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:120
Foam::LESModels::SpalartAllmarasDES::k
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition: SpalartAllmarasDES.C:456
nut
nut
Definition: createTDFields.H:71
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:87
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
DESModel< BasicTurbulenceModel >
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Foam::LESModels::SpalartAllmarasDES::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: SpalartAllmarasDES.C:412
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
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
SpalartAllmarasDES.H
Foam::LESModels::DESModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: DESModel.H:72
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
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::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::LESModels::SpalartAllmarasDES::psi
tmp< volScalarField > psi(const volScalarField &chi, const volScalarField &fv1) const
Definition: SpalartAllmarasDES.C:153
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
psi
const volScalarField & psi
Definition: setRegionFluidFields.H:13
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::LESModels::SpalartAllmarasDES::ft2
tmp< volScalarField > ft2(const volScalarField &chi) const
Definition: SpalartAllmarasDES.C:68
Foam::LESModels::SpalartAllmarasDES::dTilda
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
Length scale.
Definition: SpalartAllmarasDES.C:200
Foam::eddyViscosity< LESModel< BasicTurbulenceModel > >::correct
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
Definition: eddyViscosity.C:132
Foam::tr
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:49
Foam::LESModels::DESModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: DESModel.H:70
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::LESModels::SpalartAllmarasDES::SpalartAllmarasDES
SpalartAllmarasDES(const SpalartAllmarasDES &)
Foam::LESModels::SpalartAllmarasDES::correct
virtual void correct()
Correct nuTilda and related properties.
Definition: SpalartAllmarasDES.C:491
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::LESModels::SpalartAllmarasDES::fv1
tmp< volScalarField > fv1(const volScalarField &chi) const
Definition: SpalartAllmarasDES.C:46
Foam::neg
dimensionedScalar neg(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:201
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Foam::LESModels::SpalartAllmarasDES::fw
tmp< volScalarField > fw(const volScalarField &Omega, const volScalarField &dTilda) const
Definition: SpalartAllmarasDES.C:139
Foam::LESModels::SpalartAllmarasDES::Omega
tmp< volScalarField > Omega(const volTensorField &gradU) const
Definition: SpalartAllmarasDES.C:78
Foam::LESModels::SpalartAllmarasDES::correctNut
virtual void correctNut()
Definition: SpalartAllmarasDES.C:226
Foam::LESModels::DESModel
Definition: DESModel.H:51
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::Switch::lookupOrAddToDict
static Switch lookupOrAddToDict(const word &, dictionary &, const Switch &defaultValue=false)
Construct from dictionary, supplying default value so that if the.
Definition: Switch.C:107