SpalartAllmarasDES.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-2020 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 "SpalartAllmarasDES.H"
27 #include "fvOptions.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace LESModels
34 {
35 
36 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
37 
38 template<class BasicMomentumTransportModel>
40 {
41  return volScalarField::New
42  (
43  modelName("chi"),
44  nuTilda_/this->nu()
45  );
46 }
47 
48 
49 template<class BasicMomentumTransportModel>
51 (
52  const volScalarField& chi
53 ) const
54 {
55  const volScalarField chi3("chi3", pow3(chi));
56  return volScalarField::New
57  (
58  modelName("fv1"),
59  chi3/(chi3 + pow3(Cv1_))
60  );
61 }
62 
63 
64 template<class BasicMomentumTransportModel>
67 (
68  const volScalarField::Internal& chi,
69  const volScalarField::Internal& fv1
70 ) const
71 {
73  (
74  modelName("fv2"),
75  1.0 - chi/(1.0 + chi*fv1)
76  );
77 }
78 
79 
80 template<class BasicMomentumTransportModel>
83 (
84  const volTensorField::Internal& gradU
85 ) const
86 {
88  (
89  modelName("Omega"),
90  sqrt(2.0)*mag(skew(gradU))
91  );
92 }
93 
94 
95 template<class BasicMomentumTransportModel>
98 (
99  const volScalarField::Internal& chi,
100  const volScalarField::Internal& fv1,
101  const volScalarField::Internal& Omega,
102  const volScalarField::Internal& dTilda
103 ) const
104 {
106  (
107  modelName("Stilda"),
108  max
109  (
110  Omega
111  + fv2(chi, fv1)*nuTilda_()/sqr(kappa_*dTilda),
112  Cs_*Omega
113  )
114  );
115 }
116 
117 
118 template<class BasicMomentumTransportModel>
120 (
121  const volScalarField::Internal& nur,
122  const volScalarField::Internal& Omega,
123  const volScalarField::Internal& dTilda
124 ) const
125 {
127  (
128  modelName("r"),
129  min
130  (
131  nur
132  /(
133  max
134  (
135  Omega,
136  dimensionedScalar(Omega.dimensions(), small)
137  )
138  *sqr(kappa_*dTilda)
139  ),
140  scalar(10)
141  )
142  );
143 }
144 
145 
146 template<class BasicMomentumTransportModel>
149 (
150  const volScalarField::Internal& Omega,
151  const volScalarField::Internal& dTilda
152 ) const
153 {
154  const volScalarField::Internal r(this->r(nuTilda_, Omega, dTilda));
155  const volScalarField::Internal g(modelName("g"), r + Cw2_*(pow6(r) - r));
156 
158  (
159  modelName("fw"),
160  g*pow((1 + pow6(Cw3_))/(pow6(g) + pow6(Cw3_)), 1.0/6.0)
161  );
162 }
163 
164 
165 template<class BasicMomentumTransportModel>
168 (
169  const volScalarField::Internal& chi,
170  const volScalarField::Internal& fv1,
171  const volTensorField::Internal& gradU
172 ) const
173 {
175  (
176  modelName("dTilda"),
177  min(CDES_*this->delta()(), y_)
178  );
179 }
180 
181 
182 template<class BasicMomentumTransportModel>
184 (
185  const volScalarField::Internal& dTilda
186 ) const
187 {
188  if (this->mesh_.cacheTemporaryObject(modelName("LESRegion")))
189  {
191  (
192  modelName("LESRegion"),
193  neg(dTilda - y_())
194  );
195  }
196 }
197 
198 
199 template<class BasicMomentumTransportModel>
201 (
202  const volScalarField& fv1
203 )
204 {
205  this->nut_ = nuTilda_*fv1;
206  this->nut_.correctBoundaryConditions();
207  fv::options::New(this->mesh_).correct(this->nut_);
208 }
209 
210 
211 template<class BasicMomentumTransportModel>
213 {
214  correctNut(fv1(this->chi()));
215 }
216 
217 
218 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
219 
220 template<class BasicMomentumTransportModel>
222 (
223  const alphaField& alpha,
224  const rhoField& rho,
225  const volVectorField& U,
226  const surfaceScalarField& alphaRhoPhi,
227  const surfaceScalarField& phi,
228  const transportModel& transport,
229  const word& type
230 )
231 :
233  (
234  type,
235  alpha,
236  rho,
237  U,
238  alphaRhoPhi,
239  phi,
240  transport
241  ),
242 
243  sigmaNut_
244  (
246  (
247  "sigmaNut",
248  this->coeffDict_,
249  0.66666
250  )
251  ),
252  kappa_
253  (
255  (
256  "kappa",
257  this->coeffDict_,
258  0.41
259  )
260  ),
261  Cb1_
262  (
264  (
265  "Cb1",
266  this->coeffDict_,
267  0.1355
268  )
269  ),
270  Cb2_
271  (
273  (
274  "Cb2",
275  this->coeffDict_,
276  0.622
277  )
278  ),
279  Cw1_(Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_),
280  Cw2_
281  (
283  (
284  "Cw2",
285  this->coeffDict_,
286  0.3
287  )
288  ),
289  Cw3_
290  (
292  (
293  "Cw3",
294  this->coeffDict_,
295  2.0
296  )
297  ),
298  Cv1_
299  (
301  (
302  "Cv1",
303  this->coeffDict_,
304  7.1
305  )
306  ),
307  Cs_
308  (
310  (
311  "Cs",
312  this->coeffDict_,
313  0.3
314  )
315  ),
316  CDES_
317  (
319  (
320  "CDES",
321  this->coeffDict_,
322  0.65
323  )
324  ),
325  ck_
326  (
328  (
329  "ck",
330  this->coeffDict_,
331  0.07
332  )
333  ),
334 
335  nuTilda_
336  (
337  IOobject
338  (
339  "nuTilda",
340  this->runTime_.timeName(),
341  this->mesh_,
344  ),
345  this->mesh_
346  ),
347 
348  y_(wallDist::New(this->mesh_).y())
349 {
350  if (type == typeName)
351  {
352  this->printCoeffs(type);
353  }
354 }
355 
356 
357 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
358 
359 template<class BasicMomentumTransportModel>
361 {
363  {
364  sigmaNut_.readIfPresent(this->coeffDict());
365  kappa_.readIfPresent(*this);
366 
367  Cb1_.readIfPresent(this->coeffDict());
368  Cb2_.readIfPresent(this->coeffDict());
369  Cw1_ = Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_;
370  Cw2_.readIfPresent(this->coeffDict());
371  Cw3_.readIfPresent(this->coeffDict());
372  Cv1_.readIfPresent(this->coeffDict());
373  Cs_.readIfPresent(this->coeffDict());
374 
375  CDES_.readIfPresent(this->coeffDict());
376  ck_.readIfPresent(this->coeffDict());
377 
378  return true;
379  }
380  else
381  {
382  return false;
383  }
384 }
385 
386 
387 template<class BasicMomentumTransportModel>
389 DnuTildaEff() const
390 {
391  return volScalarField::New
392  (
393  "DnuTildaEff",
394  (nuTilda_ + this->nu())/sigmaNut_
395  );
396 }
397 
398 
399 template<class BasicMomentumTransportModel>
401 {
402  const volScalarField chi(this->chi());
403  const volScalarField fv1(this->fv1(chi));
404 
405  volScalarField dTildaExtrapolated
406  (
407  IOobject
408  (
409  "dTildaExtrapolated",
410  this->mesh_.time().timeName(),
411  this->mesh_
412  ),
413  this->mesh_,
414  dimLength,
415  zeroGradientFvPatchScalarField::typeName
416  );
417  dTildaExtrapolated.ref() = dTilda(chi, fv1, fvc::grad(this->U_));
418  dTildaExtrapolated.correctBoundaryConditions();
419 
421  (
423  (
424  modelName("k"),
425  sqr(this->nut()/ck_/dTildaExtrapolated)
426  )
427  );
428 
429  const fvPatchList& patches = this->mesh_.boundary();
430  volScalarField::Boundary& kBf = tk.ref().boundaryFieldRef();
431 
432  forAll(patches, patchi)
433  {
434  if (isA<wallFvPatch>(patches[patchi]))
435  {
436  kBf[patchi] = 0;
437  }
438  }
439 
440  return tk;
441 }
442 
443 
444 template<class BasicMomentumTransportModel>
446 {
447  if (!this->turbulence_)
448  {
449  return;
450  }
451 
452  // Local references
453  const alphaField& alpha = this->alpha_;
454  const rhoField& rho = this->rho_;
455  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
456  const volVectorField& U = this->U_;
457  fv::options& fvOptions(fv::options::New(this->mesh_));
458 
460 
461  const volScalarField chi(this->chi());
462  const volScalarField fv1(this->fv1(chi));
463 
464  tmp<volTensorField> tgradU = fvc::grad(U);
465  const volScalarField::Internal Omega(this->Omega(tgradU()));
466  const volScalarField::Internal dTilda(this->dTilda(chi, fv1, tgradU()));
467  const volScalarField::Internal Stilda
468  (
469  this->Stilda(chi, fv1, Omega, dTilda)
470  );
471  tgradU.clear();
472 
473  tmp<fvScalarMatrix> nuTildaEqn
474  (
475  fvm::ddt(alpha, rho, nuTilda_)
476  + fvm::div(alphaRhoPhi, nuTilda_)
477  - fvm::laplacian(alpha*rho*DnuTildaEff(), nuTilda_)
478  - Cb2_/sigmaNut_*alpha*rho*magSqr(fvc::grad(nuTilda_))
479  ==
480  Cb1_*alpha()*rho()*Stilda*nuTilda_()
481  - fvm::Sp
482  (
483  Cw1_*alpha()*rho()*fw(Stilda, dTilda)*nuTilda_()/sqr(dTilda),
484  nuTilda_
485  )
486  + fvOptions(alpha, rho, nuTilda_)
487  );
488 
489  nuTildaEqn.ref().relax();
490  fvOptions.constrain(nuTildaEqn.ref());
491  solve(nuTildaEqn);
492  fvOptions.correct(nuTilda_);
493  bound(nuTilda_, dimensionedScalar(nuTilda_.dimensions(), 0));
494  nuTilda_.correctBoundaryConditions();
495 
496  correctNut();
497 
498  // Optionally cache the LESRegion field
499  cacheLESRegion(dTilda);
500 }
501 
502 
503 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
504 
505 } // End namespace LESModels
506 } // End namespace Foam
507 
508 // ************************************************************************* //
Eddy viscosity LES SGS model base class.
scalar delta
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
BasicMomentumTransportModel::rhoField rhoField
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:230
fv::options & fvOptions
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &)
Return a temporary field constructed from name, mesh.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedTensor skew(const dimensionedTensor &dt)
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
dimensionedScalar sqrt(const dimensionedScalar &ds)
bool readIfPresent(const dictionary &)
Update the value of dimensioned<Type> if found in the dictionary.
virtual bool read()
Read model coefficients if they have changed.
Generic dimensioned Type class.
patches[0]
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Finite-volume options.
Definition: fvOptions.H:52
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
dimensionedScalar neg(const dimensionedScalar &ds)
phi
Definition: pEqn.H:104
static const wallDist & New(const fvMesh &mesh)
Definition: MeshObject.C:44
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
scalar y
A class for handling words, derived from string.
Definition: word.H:59
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
virtual tmp< volScalarField::Internal > dTilda(const volScalarField::Internal &chi, const volScalarField::Internal &fv1, const volTensorField::Internal &gradU) const
Length scale.
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
tmp< volScalarField::Internal > r(const volScalarField::Internal &nur, const volScalarField::Internal &Omega, const volScalarField::Internal &dTilda) const
tmp< volScalarField > fv1(const volScalarField &chi) const
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Internal & ref()
Return a reference to the dimensioned internal field.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
BasicMomentumTransportModel::alphaField alphaField
label patchi
U
Definition: pEqn.H:72
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
SpalartAllmarasDES(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
tmp< volScalarField::Internal > Omega(const volTensorField::Internal &gradU) const
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
BasicMomentumTransportModel::transportModel transportModel
virtual void cacheLESRegion(const volScalarField::Internal &dTilda) const
Cache the LES region indicator field.
void correctBoundaryConditions()
Correct boundary field.
dimensionedScalar pow6(const dimensionedScalar &ds)
tmp< volScalarField::Internal > Stilda(const volScalarField::Internal &chi, const volScalarField::Internal &fv1, const volScalarField::Internal &Omega, const volScalarField::Internal &dTilda) const
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:101
const dimensionedVector & g
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
tmp< volScalarField > chi() const
virtual void correct()
Correct nuTilda and related properties.
tmp< volScalarField::Internal > fw(const volScalarField::Internal &Omega, const volScalarField::Internal &dTilda) const
tmp< volScalarField::Internal > fv2(const volScalarField::Internal &chi, const volScalarField::Internal &fv1) const
scalar nut
Namespace for OpenFOAM.
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.