kOmegaSSTIDDES.H
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) 2015 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 Class
25  Foam::LESModels::kOmegaSSTIDDES
26 
27 Group
28  grpDESTurbulence
29 
30 Description
31  k-omega-SST IDDES turbulence model for incompressible and compressible
32  flows
33 
34  Reference:
35  \verbatim
36  Gritskevich, M.S., Garbaruk, A.V., Schuetze, J., Menter, F.R. (2011)
37  Development of DDES and IDDES Formulations for the k-omega
38  Shear Stress Transport Model, Flow, Turbulence and Combustion,
39  pp. 1-19
40  \endverbatim
41 
42 SourceFiles
43  kOmegaSSTIDDES.C
44 
45 \*---------------------------------------------------------------------------*/
46 
47 #ifndef kOmegaSSTIDDES_H
48 #define kOmegaSSTIDDES_H
49 
50 #include "kOmegaSSTDES.H"
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 namespace Foam
55 {
56 namespace LESModels
57 {
58 
59 /*---------------------------------------------------------------------------*\
60  class kOmegaSSTIDDES Declaration
61 \*---------------------------------------------------------------------------*/
62 
63 template<class BasicTurbulenceModel>
64 class kOmegaSSTIDDES
65 :
66  public kOmegaSSTDES<BasicTurbulenceModel>
67 {
68  // Private Member Functions
69 
70  //- Check that the supplied delta is an IDDESDelta
71  const IDDESDelta& setDelta() const;
72 
73  tmp<volScalarField> alpha() const;
74  tmp<volScalarField> ft(const volScalarField& magGradU) const;
75  tmp<volScalarField> fl(const volScalarField& magGradU) const;
76 
78  (
79  const volScalarField& nur,
80  const volScalarField& magGradU
81  ) const;
82 
83  //- Delay function
84  tmp<volScalarField> fdt(const volScalarField& magGradU) const;
85 
86  // Disallow default bitwise copy construct and assignment
89 
90 
91 protected:
92 
93  // Protected data
94 
95  // Model coefficients
96 
101 
102  // Fields
103 
104  const IDDESDelta& IDDESDelta_;
105 
106  // Protected Member Functions
107 
108  //- Length scale
110  (
111  const volScalarField& magGradU,
112  const volScalarField& CDES
113  ) const;
114 
115 
116 public:
117 
118  typedef typename BasicTurbulenceModel::alphaField alphaField;
119  typedef typename BasicTurbulenceModel::rhoField rhoField;
120  typedef typename BasicTurbulenceModel::transportModel transportModel;
121 
122 
123  //- Runtime type information
124  TypeName("kOmegaSSTIDDES");
125 
126 
127  // Constructors
128 
129  //- Construct from components
131  (
132  const alphaField& alpha,
133  const rhoField& rho,
134  const volVectorField& U,
135  const surfaceScalarField& alphaRhoPhi,
136  const surfaceScalarField& phi,
137  const transportModel& transport,
138  const word& propertiesName = turbulenceModel::propertiesName,
139  const word& type = typeName
140  );
141 
142 
143  //- Destructor
144  virtual ~kOmegaSSTIDDES()
145  {}
146 
147 
148  // Member Functions
149 
150  //- Re-read model coefficients if they have changed
151  virtual bool read();
152 };
153 
154 
155 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
156 
157 } // End namespace LESModels
158 } // End namespace Foam
159 
160 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
161 
162 #ifdef NoRepository
163 # include "kOmegaSSTIDDES.C"
164 #endif
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 
168 #endif
169 
170 // ************************************************************************* //
Foam::LESModels::kOmegaSSTIDDES::fl
tmp< volScalarField > fl(const volScalarField &magGradU) const
Definition: kOmegaSSTIDDES.C:70
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::LESModels::kOmegaSSTDES
k-omega-SST DES turbulence model for incompressible and compressible flows
Definition: kOmegaSSTDES.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::LESModels::kOmegaSSTIDDES::fdt
tmp< volScalarField > fdt(const volScalarField &magGradU) const
Delay function.
Definition: kOmegaSSTIDDES.C:111
Foam::LESModels::kOmegaSSTIDDES::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: kOmegaSSTIDDES.H:117
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:97
kOmegaSSTIDDES.C
Foam::LESModels::IDDESDelta
Definition: IDDESDelta.H:52
Foam::LESModels::kOmegaSSTIDDES::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: kOmegaSSTIDDES.H:119
Foam::LESModels::kOmegaSSTIDDES
k-omega-SST IDDES turbulence model for incompressible and compressible flows
Definition: kOmegaSSTIDDES.H:63
Foam::LESModels::kOmegaSSTIDDES::rd
tmp< volScalarField > rd(const volScalarField &nur, const volScalarField &magGradU) const
Definition: kOmegaSSTIDDES.C:80
Foam::LESModels::kOmegaSSTIDDES::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kOmegaSSTIDDES.C:234
U
U
Definition: pEqn.H:46
Foam::LESModels::kOmegaSSTIDDES::Cdt2_
dimensionedScalar Cdt2_
Definition: kOmegaSSTIDDES.H:97
Foam::LESModels::kOmegaSSTIDDES::Ct_
dimensionedScalar Ct_
Definition: kOmegaSSTIDDES.H:99
Foam::LESModels::kOmegaSSTIDDES::ft
tmp< volScalarField > ft(const volScalarField &magGradU) const
Definition: kOmegaSSTIDDES.C:60
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::LESModels::kOmegaSSTIDDES::IDDESDelta_
const IDDESDelta & IDDESDelta_
Definition: kOmegaSSTIDDES.H:103
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::LESModels::kOmegaSSTDES::CDES
virtual tmp< volScalarField > CDES(const volScalarField &F1) const
Blending for CDES parameter.
Definition: kOmegaSSTDES.H:90
Foam::LESModels::kOmegaSSTDES::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: kOmegaSSTDES.H:110
Foam::LESModels::kOmegaSSTIDDES::Cl_
dimensionedScalar Cl_
Definition: kOmegaSSTIDDES.H:98
rho
rho
Definition: pEqn.H:3
Foam::LESModels::kOmegaSSTIDDES::TypeName
TypeName("kOmegaSSTIDDES")
Runtime type information.
Foam::LESModels::kOmegaSSTIDDES::dTilda
virtual tmp< volScalarField > dTilda(const volScalarField &magGradU, const volScalarField &CDES) const
Length scale.
Definition: kOmegaSSTIDDES.C:121
Foam::LESModels::kOmegaSSTIDDES::~kOmegaSSTIDDES
virtual ~kOmegaSSTIDDES()
Destructor.
Definition: kOmegaSSTIDDES.H:143
Foam::LESModels::kOmegaSSTIDDES::operator=
kOmegaSSTIDDES & operator=(const kOmegaSSTIDDES &)
Foam::LESModels::kOmegaSSTDES::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: kOmegaSSTDES.H:108
Foam::LESModels::kOmegaSSTIDDES::alpha
tmp< volScalarField > alpha() const
Definition: kOmegaSSTIDDES.C:52
Foam::LESModels::kOmegaSSTIDDES::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: kOmegaSSTIDDES.H:118
Foam::LESModels::kOmegaSSTIDDES::setDelta
const IDDESDelta & setDelta() const
Check that the supplied delta is an IDDESDelta.
Definition: kOmegaSSTIDDES.C:38
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::LESModels::kOmegaSSTIDDES::kOmegaSSTIDDES
kOmegaSSTIDDES(const kOmegaSSTIDDES &)
Foam::LESModels::kOmegaSSTIDDES::Cdt1_
dimensionedScalar Cdt1_
Definition: kOmegaSSTIDDES.H:96
Foam::LESModels::kOmegaSSTDES::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: kOmegaSSTDES.H:109
kOmegaSSTDES.H