SpalartAllmarasDDES.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-2018 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 "SpalartAllmarasDDES.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 namespace LESModels
33 {
34 
35 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
36 
37 template<class BasicTurbulenceModel>
38 tmp<volScalarField> SpalartAllmarasDDES<BasicTurbulenceModel>::rd
39 (
40  const volScalarField& magGradU
41 ) const
42 {
43  tmp<volScalarField> tr
44  (
45  min
46  (
47  this->nuEff()
48  /(
49  max
50  (
51  magGradU,
52  dimensionedScalar(magGradU.dimensions(), small)
53  )
54  *sqr(this->kappa_*this->y_)
55  ),
56  scalar(10)
57  )
58  );
59  tr.ref().boundaryFieldRef() == 0.0;
60 
61  return tr;
62 }
63 
64 
65 template<class BasicTurbulenceModel>
66 tmp<volScalarField> SpalartAllmarasDDES<BasicTurbulenceModel>::fd
67 (
68  const volScalarField& magGradU
69 ) const
70 {
71  return 1 - tanh(pow3(8*rd(magGradU)));
72 }
73 
74 
75 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
76 
77 template<class BasicTurbulenceModel>
79 (
80  const volScalarField& chi,
81  const volScalarField& fv1,
82  const volTensorField& gradU
83 ) const
84 {
85  return max
86  (
87  this->y_
88  - fd(mag(gradU))
89  *max
90  (
91  this->y_ - this->CDES_*this->delta(),
93  ),
95  );
96 }
97 
98 
99 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
100 
101 template<class BasicTurbulenceModel>
103 (
104  const alphaField& alpha,
105  const rhoField& rho,
106  const volVectorField& U,
107  const surfaceScalarField& alphaRhoPhi,
108  const surfaceScalarField& phi,
109  const transportModel& transport,
110  const word& propertiesName,
111  const word& type
112 )
113 :
115  (
116  alpha,
117  rho,
118  U,
119  alphaRhoPhi,
120  phi,
121  transport,
122  propertiesName
123  )
124 {}
125 
126 
127 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
128 
129 } // End namespace LESModels
130 } // End namespace Foam
131 
132 // ************************************************************************* //
scalar delta
dimensionedScalar tanh(const dimensionedScalar &ds)
surfaceScalarField & phi
BasicTurbulenceModel::transportModel transportModel
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
SpalartAllmarasDDES(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
BasicTurbulenceModel::alphaField alphaField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
A class for handling words, derived from string.
Definition: word.H:59
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
Length scale.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar pow3(const dimensionedScalar &ds)
U
Definition: pEqn.H:72
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
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
BasicTurbulenceModel::rhoField rhoField
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))
Namespace for OpenFOAM.