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;
127 
128 
129 protected:
130 
131  // Protected Data
132 
134 
135  // Model coefficients
136 
150 
151 
152  // Fields
153 
154  //- Turbulent kinetic energy [m2/s2]
156 
157  //- Turbulent kinetic energy dissipation rate [m2/s3]
159 
160  //- Normalised wall-normal fluctuating velocity scale [-]
162 
163  //- Elliptic relaxation factor [1/s]
165 
166  //- Turbulent time scale [s]
168 
169 
170  // Bounding values
171 
176 
177 
178  // Protected Member Functions
179 
180  //- Update nut with the latest available k, phit, and T
181  virtual void correctNut();
182 
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");
199 
200 
201  // Constructors
202 
203  //- Construct from components
205  (
206  const alphaField& alpha,
207  const rhoField& rho,
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();
225 
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  {
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  }
283 
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 // ************************************************************************* //
Foam::RASModels::kEpsilonPhitF::sigmaK_
dimensionedScalar sigmaK_
Definition: kEpsilonPhitF.H:162
Foam::RASModels::kEpsilonPhitF::Ceta_
dimensionedScalar Ceta_
Definition: kEpsilonPhitF.H:160
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:77
Foam::RASModels::kEpsilonPhitF::epsilon_
volScalarField epsilon_
Turbulent kinetic energy dissipation rate [m2/s3].
Definition: kEpsilonPhitF.H:173
Foam::RASModels::kEpsilonPhitF::correct
virtual void correct()
Solve the transport equations and correct the turbulent viscosity.
Definition: kEpsilonPhitF.C:375
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::RASModels::kEpsilonPhitF::DepsilonEff
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon (LUU:Eq. 4)
Definition: kEpsilonPhitF.H:255
Foam::RASModels::kEpsilonPhitF::fMin_
dimensionedScalar fMin_
Definition: kEpsilonPhitF.H:188
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::RASModels::kEpsilonPhitF::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: kEpsilonPhitF.H:207
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::RASModels::kEpsilonPhitF::f_
volScalarField f_
Elliptic relaxation factor [1/s].
Definition: kEpsilonPhitF.H:179
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::RASModels::kEpsilonPhitF::sigmaEps_
dimensionedScalar sigmaEps_
Definition: kEpsilonPhitF.H:163
Foam::RASModels::kEpsilonPhitF::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: kEpsilonPhitF.H:208
Foam::RASModels::kEpsilonPhitF::CL_
dimensionedScalar CL_
Definition: kEpsilonPhitF.H:159
rho
rho
Definition: readInitialConditions.H:88
Foam::RASModels::kEpsilonPhitF::T_
volScalarField T_
Turbulent time scale [s].
Definition: kEpsilonPhitF.H:182
Foam::RASModels::kEpsilonPhitF::Cf2_
dimensionedScalar Cf2_
Definition: kEpsilonPhitF.H:158
Foam::RASModels::kEpsilonPhitF::Ls
tmp< volScalarField > Ls() const
Return the turbulent length scale, L.
Definition: kEpsilonPhitF.C:74
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
eddyViscosity.H
Foam::RASModels::kEpsilonPhitF::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: kEpsilonPhitF.H:209
Foam::RASModels::kEpsilonPhitF::correctNut
virtual void correctNut()
Update nut with the latest available k, phit, and T.
Definition: kEpsilonPhitF.C:42
Foam::RASModels::kEpsilonPhitF::Cf1_
dimensionedScalar Cf1_
Definition: kEpsilonPhitF.H:157
Foam::RASModels::kEpsilonPhitF::phit
virtual tmp< volScalarField > phit() const
Return the normalised wall-normal fluctuating velocity scale field.
Definition: kEpsilonPhitF.H:294
Foam::RASModels::kEpsilonPhitF::CT_
dimensionedScalar CT_
Definition: kEpsilonPhitF.H:161
Foam::RASModels::kEpsilonPhitF::TypeName
TypeName("kEpsilonPhitF")
Runtime type information.
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::RASModels::kEpsilonPhitF::Ceps1b_
dimensionedScalar Ceps1b_
Definition: kEpsilonPhitF.H:154
Foam::RASModels::kEpsilonPhitF::phitMin_
dimensionedScalar phitMin_
Definition: kEpsilonPhitF.H:187
Foam::RASModels::kEpsilonPhitF::Ceps2_
dimensionedScalar Ceps2_
Definition: kEpsilonPhitF.H:156
RASModel.H
Foam::dimensioned< scalar >
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::RASModels::kEpsilonPhitF::epsilon
virtual tmp< volScalarField > epsilon() const
Return the turbulent kinetic energy dissipation rate field.
Definition: kEpsilonPhitF.H:288
Foam::RASModels::kEpsilonPhitF::sigmaPhit_
dimensionedScalar sigmaPhit_
Definition: kEpsilonPhitF.H:164
Foam::RASModels::kEpsilonPhitF::Ceps1a_
dimensionedScalar Ceps1a_
Definition: kEpsilonPhitF.H:153
U
U
Definition: pEqn.H:72
Foam::RASModels::kEpsilonPhitF::DphitEff
tmp< volScalarField > DphitEff() const
Return the effective diffusivity for phit (LUU:Eq. 17)
Definition: kEpsilonPhitF.H:268
Foam::RASModels::kEpsilonPhitF::DkEff
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k (LUU:Eq. 3)
Definition: kEpsilonPhitF.H:242
Foam::RASModels::kEpsilonPhitF::includeNu_
Switch includeNu_
Definition: kEpsilonPhitF.H:148
Foam::RASModels::kEpsilonPhitF::k
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy field.
Definition: kEpsilonPhitF.H:282
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::RASModels::kEpsilonPhitF::Cmu_
dimensionedScalar Cmu_
Definition: kEpsilonPhitF.H:152
Foam::eddyViscosity
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:55
Foam::RASModels::kEpsilonPhitF::~kEpsilonPhitF
virtual ~kEpsilonPhitF()=default
Destructor.
Foam::RASModels::kEpsilonPhitF::k_
volScalarField k_
Turbulent kinetic energy [m2/s2].
Definition: kEpsilonPhitF.H:170
Foam::eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_
volScalarField nut_
Definition: eddyViscosity.H:66
Foam::RASModels::kEpsilonPhitF::phit_
volScalarField phit_
Normalised wall-normal fluctuating velocity scale [-].
Definition: kEpsilonPhitF.H:176
kEpsilonPhitF.C
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::RASModels::kEpsilonPhitF::Ceps1c_
dimensionedScalar Ceps1c_
Definition: kEpsilonPhitF.H:155
Foam::RASModels::kEpsilonPhitF
The k-epsilon-phit-f turbulence closure model for incompressible and compressible flows.
Definition: kEpsilonPhitF.H:131
Foam::RASModels::kEpsilonPhitF::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kEpsilonPhitF.C:348
Foam::RASModels::kEpsilonPhitF::f
virtual tmp< volScalarField > f() const
Return the elliptic relaxation factor field.
Definition: kEpsilonPhitF.H:300
Foam::RASModels::kEpsilonPhitF::L2Min_
dimensionedScalar L2Min_
Definition: kEpsilonPhitF.H:190
Foam::RASModels::kEpsilonPhitF::TMin_
dimensionedScalar TMin_
Definition: kEpsilonPhitF.H:189
Foam::RASModels::kEpsilonPhitF::Ts
tmp< volScalarField > Ts() const
Return the turbulent time scale, T.
Definition: kEpsilonPhitF.C:54