kEpsilonPhitF.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) 2019-2021 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::RASModels::kEpsilonPhitF
28 
29 Group
30  grpRASTurbulence
31 
32 Description
33  The k-epsilon-phit-f turbulence closure model for incompressible and
34  compressible flows.
35 
36  The model is a three-transport-equation linear-eddy-viscosity turbulence
37  closure model alongside an elliptic relaxation equation.
38 
39  \heading Input fields
40  \plaintable
41  k | Turbulent kinetic energy [m2/s2]
42  epsilon | Turbulent kinetic energy dissipation rate [m2/s3]
43  phit | Normalised wall-normal fluctuating velocity scale [-]
44  f | Elliptic relaxation factor [1/s]
45  \endplaintable
46 
47  Reference:
48  \verbatim
49  Standard model (Tag:LUU):
50  Laurence, D. R., Uribe, J. C., & Utyuzhnikov, S. V. (2005).
51  A robust formulation of the v2-f model.
52  Flow, Turbulence and Combustion, 73(3-4), 169–185.
53  DOI:10.1007/s10494-005-1974-8
54  \endverbatim
55 
56  The default model coefficients are (LUU:Eqs. 19-20):
57  \verbatim
58  kEpsilonPhitFCoeffs
59  {
60  includeNu true; // include nu in (LUU: Eq. 17), see Notes
61  Cmu 0.22; // Turbulent viscosity constant
62  Ceps1a 1.4; // Model constant for epsilon
63  Ceps1b 1.0; // Model constant for epsilon
64  Ceps1c 0.05; // Model constant for epsilon
65  Ceps2 1.9; // Model constant for epsilon
66  Cf1 1.4; // Model constant for f
67  Cf2 0.3; // Model constant for f
68  CL 0.25; // Model constant for L
69  Ceta 110.0; // Model constant for L
70  CT 6.0; // Model constant for T
71  sigmaK 1.0; // Turbulent Prandtl number for k
72  sigmaEps 1.3; // Turbulent Prandtl number for epsilon
73  sigmaPhit 1.0; // Turbulent Prandtl number for phit = sigmaK
74  }
75  \endverbatim
76 
77 Note
78  The name of the original variable replacing 'v2' is 'phi' (LUU:Eq. 14).
79  However, the name 'phi' preexisted in OpenFOAM; therefore, this name was
80  replaced by 'phit' herein.
81 
82  Including \c nu in \c DphitEff even though it is not present in (LUU:Eq. 17)
83  provided higher level of resemblance to benchmarks for the tests considered,
84  particularly for the peak skin friction (yet, pressure-related predictions
85  were unaffected). Users can switch off \c nu in \c DphitEff by using
86  \c includeNu entry in \c kEpsilonPhitFCoeffs as shown above in order to
87  follow the reference paper thereat. \c includeNu is left \c true by default.
88  See GitLab issue #1560.
89 
90 SourceFiles
91  kEpsilonPhitF.C
92 
93 SeeAlso
94  kEpsilon.C
95 
96 \*---------------------------------------------------------------------------*/
97 
98 #ifndef kEpsilonPhitF_H
99 #define kEpsilonPhitF_H
100 
101 #include "RASModel.H"
102 #include "eddyViscosity.H"
103 
104 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
105 
106 namespace Foam
107 {
108 namespace RASModels
109 {
110 
111 /*---------------------------------------------------------------------------*\
112  Class kEpsilonPhitF Declaration
113 \*---------------------------------------------------------------------------*/
114 
115 template<class BasicTurbulenceModel>
116 class kEpsilonPhitF
117 :
118  public eddyViscosity<RASModel<BasicTurbulenceModel>>
119 {
120  // Private Member Functions
121 
122  //- No copy construct
123  kEpsilonPhitF(const kEpsilonPhitF&) = delete;
124 
125  //- No copy assignment
126  void operator=(const kEpsilonPhitF&) = delete;
128 
129 protected:
130 
131  // Protected Data
132 
134 
135  // Model coefficients
136 
150 
151 
152  // Fields
154  //- Turbulent kinetic energy [m2/s2]
157  //- Turbulent kinetic energy dissipation rate [m2/s3]
160  //- Normalised wall-normal fluctuating velocity scale [-]
163  //- Elliptic relaxation factor [1/s]
165 
166  //- Turbulent time scale [s]
168 
169 
170  // Bounding values
171 
176 
178  // Protected Member Functions
179 
180  //- Update nut with the latest available k, phit, and T
181  virtual void correctNut();
183  //- Return the turbulent time scale, T
184  tmp<volScalarField> Ts() const;
185 
186  //- Return the turbulent length scale, L
188 
189 
190 public:
191 
192  typedef typename BasicTurbulenceModel::alphaField alphaField;
193  typedef typename BasicTurbulenceModel::rhoField rhoField;
194  typedef typename BasicTurbulenceModel::transportModel transportModel;
195 
196 
197  //- Runtime type information
198  TypeName("kEpsilonPhitF");
201  // Constructors
202 
203  //- Construct from components
205  (
206  const alphaField& alpha,
207  const rhoField& rho,
208  const volVectorField& U,
209  const surfaceScalarField& alphaRhoPhi,
210  const surfaceScalarField& phi,
211  const transportModel& transport,
212  const word& propertiesName = turbulenceModel::propertiesName,
213  const word& type = typeName
214  );
215 
216 
217  //- Destructor
218  virtual ~kEpsilonPhitF() = default;
219 
220 
221  // Member Functions
222 
223  //- Re-read model coefficients if they have changed
224  virtual bool read();
226  //- Return the effective diffusivity for k (LUU:Eq. 3)
227  tmp<volScalarField> DkEff() const
228  {
229  return tmp<volScalarField>
230  (
231  new volScalarField
232  (
233  "DkEff",
234  this->nut_/sigmaK_ + this->nu()
235  )
236  );
237  }
238 
239  //- Return the effective diffusivity for epsilon (LUU:Eq. 4)
241  {
242  return tmp<volScalarField>
243  (
244  new volScalarField
245  (
246  "DepsilonEff",
247  this->nut_/sigmaEps_ + this->nu()
248  )
249  );
250  }
251 
252  //- Return the effective diffusivity for phit (LUU:Eq. 17)
254  {
255  auto tfld =
256  tmp<volScalarField>::New("DphitEff", this->nut_/sigmaPhit_);
257 
258  if (includeNu_)
259  {
260  tfld.ref() += this->nu();
261  }
262 
263  return tfld;
264  }
265 
266  //- Return the turbulent kinetic energy field
267  virtual tmp<volScalarField> k() const
268  {
269  return k_;
270  }
271 
272  //- Return the turbulent kinetic energy dissipation rate field
273  virtual tmp<volScalarField> epsilon() const
274  {
275  return epsilon_;
276  }
277 
278  //- Return the normalised wall-normal fluctuating velocity scale field
279  virtual tmp<volScalarField> phit() const
280  {
281  return phit_;
282  }
284  //- Return the elliptic relaxation factor field
285  virtual tmp<volScalarField> f() const
286  {
287  return f_;
288  }
289 
290  //- Solve the transport equations and correct the turbulent viscosity
291  virtual void correct();
292 };
293 
294 
295 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
296 
297 } // End namespace RASModels
298 } // End namespace Foam
299 
300 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
301 
302 #ifdef NoRepository
303  #include "kEpsilonPhitF.C"
304 #endif
305 
306 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
307 
308 #endif
309 
310 // ************************************************************************* //
virtual void correct()
Solve the transport equations and correct the turbulent viscosity.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon (LUU:Eq. 4)
volScalarField phit_
Normalised wall-normal fluctuating velocity scale [-].
volScalarField f_
Elliptic relaxation factor [1/s].
BasicTurbulenceModel::alphaField alphaField
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, any/none. Also accepts 0/1 as a string and shortcuts t/f, y/n.
Definition: Switch.H:77
tmp< volScalarField > Ls() const
Return the turbulent length scale, L.
Definition: kEpsilonPhitF.C:67
volScalarField T_
Turbulent time scale [s].
virtual tmp< volScalarField > phit() const
Return the normalised wall-normal fluctuating velocity scale field.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
BasicTurbulenceModel::rhoField rhoField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual tmp< volScalarField > epsilon() const
Return the turbulent kinetic energy dissipation rate field.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
virtual void correctNut()
Update nut with the latest available k, phit, and T.
Definition: kEpsilonPhitF.C:35
BasicTurbulenceModel::transportModel transportModel
TypeName("kEpsilonPhitF")
Runtime type information.
U
Definition: pEqn.H:72
virtual bool read()
Re-read model coefficients if they have changed.
tmp< volScalarField > DphitEff() const
Return the effective diffusivity for phit (LUU:Eq. 17)
tmp< volScalarField > Ts() const
Return the turbulent time scale, T.
Definition: kEpsilonPhitF.C:47
The k-epsilon-phit-f turbulence closure model for incompressible and compressible flows...
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k (LUU:Eq. 3)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
volScalarField epsilon_
Turbulent kinetic energy dissipation rate [m2/s3].
volScalarField k_
Turbulent kinetic energy [m2/s2].
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual tmp< volScalarField > f() const
Return the elliptic relaxation factor field.
volScalarField & nu
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy field.
Namespace for OpenFOAM.
virtual ~kEpsilonPhitF()=default
Destructor.