v2f.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) 2012-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::RASModels::v2f
26 
27 Group
28  grpRASTurbulence
29 
30 Description
31  Lien and Kalitzin's v2-f turbulence model for incompressible and
32  compressible flows, with a limit imposed on the turbulent viscosity given
33  by Davidson et al.
34 
35  The model solves for turbulence kinetic energy k and turbulence dissipation
36  rate epsilon, with additional equations for the turbulence stress normal to
37  streamlines, v2, and elliptic damping function, f.
38 
39  The variant implemented employs N=6, such that f=0 on walls.
40 
41  Wall boundary conditions are:
42 
43  k = kLowReWallFunction
44  epsilon = epsilonLowReWallFunction
45  v2 = v2WallFunction
46  f = fWallFunction
47 
48  These are applicable to both low- and high-Reynolds number flows.
49 
50  Inlet values can be approximated by:
51 
52  v2 = 2/3 k
53  f = zero-gradient
54 
55  References:
56  \verbatim
57  Lien, F. S., & Kalitzin, G. (2001).
58  Computations of transonic flow with the v2f turbulence model.
59  International Journal of Heat and Fluid Flow, 22(1), 53-61.
60 
61  Davidson, L., Nielsen, P., & Sveningsson, A. (2003).
62  Modifications of the v2-f model for computing the flow in a
63  3D wall jet.
64  Turbulence, Heat and Mass Transfer, 4, 577-584
65  \endverbatim
66 
67  The default model coefficients are
68  \verbatim
69  v2fCoeffs
70  {
71  Cmu 0.22;
72  CmuKEps 0.09;
73  C1 1.4;
74  C2 0.3;
75  CL 0.23;
76  Ceta 70;
77  Ceps2 1.9;
78  Ceps3 -0.33;
79  sigmaEps 1.3;
80  sigmaK 1;
81  }
82  \endverbatim
83 
84 Note
85  If the kLowReWallFunction is employed, a velocity variant of the turbulent
86  viscosity wall function should be used, e.g. nutUWallFunction. Turbulence
87  k variants (nutk...) for this case will not behave correctly.
88 
89 SeeAlso
90  Foam::RASModels::v2fBase
91  Foam::RASModels::kEpsilon
92  Foam::kLowReWallFunctionFvPatchScalarField
93  Foam::epsilonLowReWallFunctionFvPatchScalarField
94  Foam::v2WallFunctionFvPatchScalarField
95  Foam::fWallFunctionFvPatchScalarField
96 
97 SourceFiles
98  v2f.C
99 
100 \*---------------------------------------------------------------------------*/
101 
102 #ifndef v2f_H
103 #define v2f_H
104 
105 #include "v2fBase.H"
106 #include "RASModel.H"
107 #include "eddyViscosity.H"
108 
109 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
110 
111 namespace Foam
112 {
113 namespace RASModels
114 {
115 
116 /*---------------------------------------------------------------------------*\
117  Class v2f Declaration
118 \*---------------------------------------------------------------------------*/
119 
120 template<class BasicTurbulenceModel>
121 class v2f
122 :
123  public eddyViscosity<RASModel<BasicTurbulenceModel> >,
124  public v2fBase
125 {
126 
127 protected:
128 
129  // Protected data
130 
131  // Model coefficients
132 
143 
144 
145  // Fields
146 
147  //- Turbulence kinetic energy
149 
150  //- Turbulence dissipation
152 
153  //- Turbulence stress normal to streamlines
155 
156  //- Damping function
158 
159 
160  // Bounding values
161 
164 
165 
166  // Protected Member Functions
167 
168  virtual void correctNut();
169 
170  //- Return time scale, Ts
171  tmp<volScalarField> Ts() const;
172 
173  //- Return length scale, Ls
174  tmp<volScalarField> Ls() const;
175 
176 
177 public:
178 
179  typedef typename BasicTurbulenceModel::alphaField alphaField;
180  typedef typename BasicTurbulenceModel::rhoField rhoField;
181  typedef typename BasicTurbulenceModel::transportModel transportModel;
182 
183 
184  //- Runtime type information
185  TypeName("v2f");
186 
187 
188  // Constructors
189 
190  //- Construct from components
191  v2f
192  (
193  const alphaField& alpha,
194  const rhoField& rho,
195  const volVectorField& U,
196  const surfaceScalarField& alphaRhoPhi,
197  const surfaceScalarField& phi,
198  const transportModel& transport,
199  const word& propertiesName = turbulenceModel::propertiesName,
200  const word& type = typeName
201  );
202 
203 
204  //- Destructor
205  virtual ~v2f()
206  {}
207 
208 
209  // Member Functions
210 
211  //- Re-read model coefficients if they have changed
212  virtual bool read();
213 
214  //- Return the effective diffusivity for k
215  tmp<volScalarField> DkEff() const
216  {
217  return tmp<volScalarField>
218  (
219  new volScalarField
220  (
221  "DkEff",
222  this->nut_/sigmaK_ + this->nu()
223  )
224  );
225  }
226 
227  //- Return the effective diffusivity for epsilon
229  {
230  return tmp<volScalarField>
231  (
232  new volScalarField
233  (
234  "DepsilonEff",
235  this->nut_/sigmaEps_ + this->nu()
236  )
237  );
238  }
239 
240  //- Return the turbulence kinetic energy
241  virtual tmp<volScalarField> k() const
242  {
243  return k_;
244  }
245 
246  //- Return the turbulence kinetic energy dissipation rate
247  virtual tmp<volScalarField> epsilon() const
248  {
249  return epsilon_;
250  }
251 
252  //- Return turbulence stress normal to streamlines
253  virtual tmp<volScalarField> v2() const
254  {
255  return v2_;
256  }
257 
258  //- Return the damping function
259  virtual tmp<volScalarField> f() const
260  {
261  return f_;
262  }
263 
264  //- Solve the turbulence equations and correct the turbulence viscosity
265  virtual void correct();
266 };
267 
268 
269 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
270 
271 } // End namespace RASModels
272 } // End namespace Foam
273 
274 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
275 
276 #ifdef NoRepository
277 # include "v2f.C"
278 #endif
279 
280 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281 
282 #endif
283 
284 // ************************************************************************* //
Foam::RASModels::v2f::fMin_
dimensionedScalar fMin_
Definition: v2f.H:162
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::RASModels::v2f::k
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: v2f.H:240
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
v2fBase.H
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::RASModels::v2f::epsilon_
volScalarField epsilon_
Turbulence dissipation.
Definition: v2f.H:150
Foam::RASModels::v2f::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: v2f.H:180
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:97
Foam::RASModels::v2f::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: v2f.H:178
Foam::RASModels::v2f::sigmaK_
dimensionedScalar sigmaK_
Definition: v2f.H:140
Foam::RASModels::v2f::Cmu_
dimensionedScalar Cmu_
Definition: v2f.H:132
Foam::RASModels::v2f::Ceps3_
dimensionedScalar Ceps3_
Definition: v2f.H:139
Foam::RASModels::v2f::sigmaEps_
dimensionedScalar sigmaEps_
Definition: v2f.H:141
Foam::RASModels::v2f::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: v2f.C:280
Foam::RASModels::v2f::CL_
dimensionedScalar CL_
Definition: v2f.H:136
U
U
Definition: pEqn.H:46
Foam::RASModels::v2f::Ts
tmp< volScalarField > Ts() const
Return time scale, Ts.
Definition: v2f.C:39
nu
volScalarField & nu
Definition: readMechanicalProperties.H:179
eddyViscosity.H
Foam::RASModels::v2f
Lien and Kalitzin's v2-f turbulence model for incompressible and compressible flows,...
Definition: v2f.H:120
Foam::RASModels::v2f::TypeName
TypeName("v2f")
Runtime type information.
Foam::RASModels::v2f::DepsilonEff
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: v2f.H:227
Foam::RASModels::v2f::v2f
v2f(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)
Construct from components.
Definition: v2f.C:99
Foam::RASModels::v2f::f
virtual tmp< volScalarField > f() const
Return the damping function.
Definition: v2f.H:258
Foam::RASModels::v2f::CmuKEps_
dimensionedScalar CmuKEps_
Definition: v2f.H:133
Foam::RASModels::v2f::C2_
dimensionedScalar C2_
Definition: v2f.H:135
Foam::RASModels::v2f::DkEff
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: v2f.H:214
Foam::RASModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:100
RASModel.H
Foam::RASModels::v2f::Ceta_
dimensionedScalar Ceta_
Definition: v2f.H:137
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::RASModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:99
Foam::RASModels::v2f::f_
volScalarField f_
Damping function.
Definition: v2f.H:156
rho
rho
Definition: pEqn.H:3
Foam::RASModels::v2f::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: v2f.C:305
Foam::RASModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:101
Foam::RASModels::v2f::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: v2f.H:179
Foam::RASModels::v2f::epsilon
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: v2f.H:246
Foam::RASModels::v2f::C1_
dimensionedScalar C1_
Definition: v2f.H:134
Foam::RASModels::v2f::k_
volScalarField k_
Turbulence kinetic energy.
Definition: v2f.H:147
Foam::RASModels::v2f::v2_
volScalarField v2_
Turbulence stress normal to streamlines.
Definition: v2f.H:153
Foam::RASModels::v2f::correctNut
virtual void correctNut()
Definition: v2f.C:86
Foam::RASModels::v2f::~v2f
virtual ~v2f()
Destructor.
Definition: v2f.H:204
Foam::eddyViscosity
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
Foam::RASModels::v2f::v2Min_
dimensionedScalar v2Min_
Definition: v2f.H:161
Foam::eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_
volScalarField nut_
Definition: eddyViscosity.H:63
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
v2f.C
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::RASModels::v2f::v2
virtual tmp< volScalarField > v2() const
Return turbulence stress normal to streamlines.
Definition: v2f.H:252
Foam::RASModels::v2fBase
Abstract base-class for v2-f models to provide BCs access to the v2 and f fields.
Definition: v2fBase.H:57
Foam::RASModels::v2f::Ceps2_
dimensionedScalar Ceps2_
Definition: v2f.H:138
Foam::RASModels::v2f::Ls
tmp< volScalarField > Ls() const
Return length scale, Ls.
Definition: v2f.C:61