kEqn.H
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 Class
25  Foam::LESModels::kEqn
26 
27 Description
28  One equation eddy-viscosity model
29 
30  Eddy viscosity SGS model using a modeled balance equation to simulate the
31  behaviour of k.
32 
33  Reference:
34  \verbatim
35  Yoshizawa, A. (1986).
36  Statistical theory for compressible turbulent shear flows,
37  with the application to subgrid modeling.
38  Physics of Fluids (1958-1988), 29(7), 2152-2164.
39  \endverbatim
40 
41  The default model coefficients are
42  \verbatim
43  kEqnCoeffs
44  {
45  Ck 0.094;
46  Ce 1.048;
47  }
48  \endverbatim
49 
50 SourceFiles
51  kEqn.C
52 
53 \*---------------------------------------------------------------------------*/
54 
55 #ifndef kEqn_H
56 #define kEqn_H
57 
58 #include "LESeddyViscosity.H"
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62 namespace Foam
63 {
64 namespace LESModels
65 {
66 
67 /*---------------------------------------------------------------------------*\
68  Class kEqn Declaration
69 \*---------------------------------------------------------------------------*/
70 
71 template<class BasicTurbulenceModel>
72 class kEqn
73 :
74  public LESeddyViscosity<BasicTurbulenceModel>
75 {
76  // Private Member Functions
77 
78  // Disallow default bitwise copy construct and assignment
79  kEqn(const kEqn&);
80  void operator=(const kEqn&);
81 
82 
83 protected:
84 
85  // Protected data
86 
87  // Fields
88 
90 
91 
92  // Model constants
93 
95 
96 
97  // Protected Member Functions
98 
99  virtual void correctNut();
100  virtual tmp<fvScalarMatrix> kSource() const;
101 
102 
103 public:
105  typedef typename BasicTurbulenceModel::alphaField alphaField;
106  typedef typename BasicTurbulenceModel::rhoField rhoField;
107  typedef typename BasicTurbulenceModel::transportModel transportModel;
108 
109 
110  //- Runtime type information
111  TypeName("kEqn");
112 
113 
114  // Constructors
115 
116  //- Constructor from components
117  kEqn
118  (
119  const alphaField& alpha,
120  const rhoField& rho,
121  const volVectorField& U,
122  const surfaceScalarField& alphaRhoPhi,
123  const surfaceScalarField& phi,
124  const transportModel& transport,
125  const word& propertiesName = turbulenceModel::propertiesName,
126  const word& type = typeName
127  );
128 
129 
130  //- Destructor
131  virtual ~kEqn()
132  {}
133 
134 
135  // Member Functions
136 
137  //- Read model coefficients if they have changed
138  virtual bool read();
139 
140  //- Return SGS kinetic energy
141  virtual tmp<volScalarField> k() const
142  {
143  return k_;
144  }
145 
146  //- Return sub-grid disipation rate
147  virtual tmp<volScalarField> epsilon() const;
148 
149  //- Return the effective diffusivity for k
150  tmp<volScalarField> DkEff() const
151  {
152  return tmp<volScalarField>
153  (
154  new volScalarField("DkEff", this->nut_ + this->nu())
155  );
156  }
157 
158  //- Correct eddy-Viscosity and related properties
159  virtual void correct();
160 };
161 
162 
163 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164 
165 } // End namespace LESModels
166 } // End namespace Foam
167 
168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 
170 #ifdef NoRepository
171  #include "kEqn.C"
172 #endif
173 
174 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
175 
176 #endif
177 
178 // ************************************************************************* //
Eddy viscosity LES SGS model base class.
virtual tmp< volScalarField > epsilon() const
Return sub-grid disipation rate.
Definition: kEqn.C:142
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition: kEqn.H:140
surfaceScalarField & phi
TypeName("kEqn")
Runtime type information.
One equation eddy-viscosity model.
Definition: kEqn.H:71
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual ~kEqn()
Destructor.
Definition: kEqn.H:130
static const word propertiesName
Default name of the turbulence properties dictionary.
volScalarField k_
Definition: kEqn.H:88
A class for handling words, derived from string.
Definition: word.H:59
BasicTurbulenceModel::transportModel transportModel
Definition: kEqn.H:106
virtual void correct()
Correct eddy-Viscosity and related properties.
Definition: kEqn.C:163
virtual tmp< fvScalarMatrix > kSource() const
Definition: kEqn.C:50
virtual bool read()
Read model coefficients if they have changed.
Definition: kEqn.C:126
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: kEqn.H:149
BasicTurbulenceModel::alphaField alphaField
Definition: kEqn.H:104
virtual void correctNut()
Definition: kEqn.C:39
U
Definition: pEqn.H:72
dimensionedScalar Ck_
Definition: kEqn.H:93
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
volScalarField & nu
BasicTurbulenceModel::rhoField rhoField
Definition: kEqn.H:105
Namespace for OpenFOAM.