Loading...
Searching...
No Matches
kEqn.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2019-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::LESModels::kEqn
29
30Group
31 grpLESTurbulence
32
33Description
34 One equation eddy-viscosity model
35
36 Eddy viscosity SGS model using a modeled balance equation to simulate the
37 behaviour of k.
38
39 Reference:
40 \verbatim
41 Yoshizawa, A. (1986).
42 Statistical theory for compressible turbulent shear flows,
43 with the application to subgrid modeling.
44 Physics of Fluids (1958-1988), 29(7), 2152-2164.
45 \endverbatim
46
47 The default model coefficients are
48 \verbatim
49 kEqnCoeffs
50 {
51 Ck 0.094;
52 Ce 1.048;
53 }
54 \endverbatim
55
56SourceFiles
57 kEqn.C
58
59\*---------------------------------------------------------------------------*/
60
61#ifndef kEqn_H
62#define kEqn_H
63
64#include "LESeddyViscosity.H"
65
66// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67
68namespace Foam
69{
70namespace LESModels
71{
72
73/*---------------------------------------------------------------------------*\
74 Class kEqn Declaration
75\*---------------------------------------------------------------------------*/
76
77template<class BasicTurbulenceModel>
78class kEqn
79:
80 public LESeddyViscosity<BasicTurbulenceModel>
81{
82protected:
84 // Protected data
85
86 // Fields
87
89
90
91 // Model constants
92
94
95
96 // Protected Member Functions
97
98 //- No copy construct
99 kEqn(const kEqn&) = delete;
100
101 //- No copy assignment
102 void operator=(const kEqn&) = delete;
103
105 virtual void correctNut();
106 virtual tmp<fvScalarMatrix> kSource() const;
107
108
109public:
111 typedef typename BasicTurbulenceModel::alphaField alphaField;
112 typedef typename BasicTurbulenceModel::rhoField rhoField;
113 typedef typename BasicTurbulenceModel::transportModel transportModel;
114
115
116 //- Runtime type information
117 TypeName("kEqn");
119
120 // Constructors
121
122 //- Constructor from components
123 kEqn
124 (
125 const alphaField& alpha,
126 const rhoField& rho,
127 const volVectorField& U,
128 const surfaceScalarField& alphaRhoPhi,
129 const surfaceScalarField& phi,
130 const transportModel& transport,
131 const word& propertiesName = turbulenceModel::propertiesName,
132 const word& type = typeName
133 );
134
135
136 //- Destructor
137 virtual ~kEqn() = default;
138
139
140 // Member Functions
141
142 //- Read model coefficients if they have changed
143 virtual bool read();
144
145 //- Return SGS kinetic energy
146 virtual tmp<volScalarField> k() const
147 {
148 return k_;
149 }
151 //- Return the effective diffusivity for k
153 {
156 new volScalarField("DkEff", this->nut_ + this->nu())
157 );
158 }
159
160 //- Correct eddy-Viscosity and related properties
161 virtual void correct();
162};
164
165// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166
167} // End namespace LESModels
168} // End namespace Foam
169
170// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171
172#ifdef NoRepository
173 #include "kEqn.C"
174#endif
175
176// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177
178#endif
179
180// ************************************************************************* //
Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel >::alphaField alphaField
Definition kEqn.H:110
kEqn(const kEqn &)=delete
No copy construct.
Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel >::rhoField rhoField
Definition kEqn.H:111
virtual ~kEqn()=default
Destructor.
virtual void correct()
Correct eddy-Viscosity and related properties.
Definition kEqn.C:132
kEqn(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)
Constructor from components.
Definition kEqn.C:60
void operator=(const kEqn &)=delete
No copy assignment.
virtual void correctNut()
Definition kEqn.C:35
virtual tmp< fvScalarMatrix > kSource() const
Definition kEqn.C:46
TypeName("kEqn")
Runtime type information.
Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel >::transportModel transportModel
Definition kEqn.H:112
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition kEqn.H:155
virtual bool read()
Read model coefficients if they have changed.
Definition kEqn.C:118
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition kEqn.H:163
A class for managing temporary objects.
Definition tmp.H:75
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
Namespace for LES SGS models.
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volScalarField & nu
volScalarField & alpha
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68