kLowReWallFunctionFvPatchScalarField.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) 2012-2016, 2019 OpenFOAM Foundation
9  Copyright (C) 2019-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
27 Class
28  Foam::kLowReWallFunctionFvPatchScalarField
29 
30 Group
31  grpWallFunctions
32 
33 Description
34  This boundary condition provides a wall function for the turbulent kinetic
35  energy (i.e. \c k) for low- and high-Reynolds number simulations.
36 
37 Usage
38  Example of the boundary condition specification:
39  \verbatim
40  <patchName>
41  {
42  // Mandatory entries
43  type kLowReWallFunction;
44 
45  // Optional entries
46  Ceps2 1.9;
47  Ck -0.416;
48  Bk 8.366;
49  C 11.0;
50 
51  // Inherited entries
52  ...
53  }
54  \endverbatim
55 
56  where the entries mean:
57  \table
58  Property | Description | Type | Reqd | Deflt
59  type | Type name: kLowReWallFunction | word | yes | -
60  Ceps2 | Model coefficient | scalar | no | 1.9
61  Ck | Model coefficient | scalar | no | -0.416
62  Bk | Model coefficient | scalar | no | 8.366
63  C | Model coefficient | scalar | no | 11.0
64  \endtable
65 
66  The inherited entries are elaborated in:
67  - \link fixedValueFvPatchField.H \endlink
68  - \link wallFunctionCoefficients.H \endlink
69 
70  Viscous and inertial sublayer predictions for \c k are blended in
71  a stepwise manner:
72 
73  \f[
74  k = k_{log} \qquad if \quad y^+ > y^+_{intersection}
75  \f]
76  \f[
77  k = k_{vis} \qquad if \quad y^+ <= y^+_{intersection}
78  \f]
79  where
80  \vartable
81  k_{vis} | k prediction in the viscous sublayer
82  k_{log} | k prediction in the inertial sublayer
83  y^+ | estimated wall-normal height of the cell centre in wall units
84  y^+_{intersection} | estimated \f$y^+\f$ where sublayers intersect
85  \endvartable
86 
87 SourceFiles
88  kLowReWallFunctionFvPatchScalarField.C
89 
90 \*---------------------------------------------------------------------------*/
91 
92 #ifndef kLowReWallFunctionFvPatchScalarField_H
93 #define kLowReWallFunctionFvPatchScalarField_H
94 
95 #include "fixedValueFvPatchField.H"
97 
98 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
99 
100 namespace Foam
101 {
102 
103 /*---------------------------------------------------------------------------*\
104  Class kLowReWallFunctionFvPatchScalarField Declaration
105 \*---------------------------------------------------------------------------*/
106 
107 class kLowReWallFunctionFvPatchScalarField
108 :
109  public fixedValueFvPatchField<scalar>
110 {
111 protected:
112 
113  // Protected Data
114 
115  //- Ceps2 coefficient
116  scalar Ceps2_;
117 
118  //- Ck coefficient
119  scalar Ck_;
120 
121  //- Bk coefficient
122  scalar Bk_;
123 
124  //- C coefficient
125  scalar C_;
126 
127  //- Wall-function coefficients
128  wallFunctionCoefficients wallCoeffs_;
129 
130 
131  // Protected Member Functions
132 
133  //- Write local wall function variables
134  void writeLocalEntries(Ostream&) const;
135 
136 
137 public:
138 
139  //- Runtime type information
140  TypeName("kLowReWallFunction");
141 
142 
143  // Constructors
144 
145  //- Construct from patch and internal field
147  (
148  const fvPatch&,
149  const DimensionedField<scalar, volMesh>&
150  );
151 
152  //- Construct from patch, internal field and dictionary
154  (
155  const fvPatch&,
157  const dictionary&
158  );
159 
160  //- Construct by mapping given kLowReWallFunctionFvPatchScalarField
161  //- onto a new patch
163  (
165  const fvPatch&,
167  const fvPatchFieldMapper&
168  );
169 
170  //- Construct as copy
172  (
174  );
176  //- Construct as copy setting internal field reference
178  (
181  );
182 
183  //- Return a clone
184  virtual tmp<fvPatchField<scalar>> clone() const
185  {
186  return fvPatchField<scalar>::Clone(*this);
187  }
188 
189  //- Clone with an internal field reference
191  (
193  ) const
194  {
195  return fvPatchField<scalar>::Clone(*this, iF);
196  }
197 
198 
199  // Member Functions
200 
201  // Evaluation
202 
203  //- Update the coefficients associated with the patch field
204  virtual void updateCoeffs();
205 
206 
207  // I-O
208 
209  //- Write
210  virtual void write(Ostream&) const;
211 };
212 
213 
214 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215 
216 } // End namespace Foam
217 
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 
220 #endif
221 
222 // ************************************************************************* //
wallFunctionCoefficients wallCoeffs_
Wall-function coefficients.
virtual tmp< fvPatchField< scalar > > clone() const
Return a clone.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
This boundary condition provides a wall function for the turbulent kinetic energy (i...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
TypeName("kLowReWallFunction")
Runtime type information.
static tmp< fvPatchField< Type > > Clone(const DerivedPatchField &pf, Args &&... args)
Clone a patch field, optionally with internal field reference etc.
Definition: fvPatchField.H:597
A FieldMapper for finite-volume patch fields.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A class for managing temporary objects.
Definition: HashPtrTable.H:50
void writeLocalEntries(Ostream &) const
Write local wall function variables.
kLowReWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Namespace for OpenFOAM.