nutUBlendedWallFunctionFvPatchScalarField.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) 2016-2022 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::nutUBlendedWallFunctionFvPatchScalarField
28 
29 Group
30  grpWallFunctions
31 
32 Description
33  This boundary condition provides a wall function for the turbulent
34  viscosity (i.e. \c nut) based on velocity (i.e. \c U) using a
35  binomial-function wall-function blending method between the viscous
36  and inertial sublayer predictions of \c nut for low- and high-Reynolds
37  number applications.
38 
39  \f[
40  u_\tau = (u_{\tau,vis}^n + u_{\tau,log}^n)^{1/n}
41  \f]
42 
43  where
44  \vartable
45  u_\tau | Friction velocity
46  u_{\tau,vis} | Friction velocity in the viscous sublayer
47  u_{\tau,log} | Friction velocity in the inertial sublayer
48  \endvartable
49 
50  Reference:
51  \verbatim
52  See the section that describes 'automatic wall treatment':
53  Menter, F., Ferreira, J. C., Esch, T., Konno, B. (2003).
54  The SST turbulence model with improved wall treatment
55  for heat transfer predictions in gas turbines.
56  In Proceedings of the International Gas Turbine Congress.
57  November, 2003. Tokyo, Japan. pp. 2-7.
58  \endverbatim
59 
60 Usage
61  Example of the boundary condition specification:
62  \verbatim
63  <patchName>
64  {
65  // Mandatory entries
66  type nutUBlendedWallFunction;
67 
68  // Optional entries
69  n 4.0;
70 
71  // Inherited entries
72  ...
73  }
74  \endverbatim
75 
76  where the entries mean:
77  \table
78  Property | Description | Type | Reqd | Deflt
79  type | Type name: nutUBlendedWallFunction | word | yes | -
80  n | Blending factor | scalar | no | 4.0
81  \endtable
82 
83  The inherited entries are elaborated in:
84  - \link nutWallFunctionFvPatchScalarField.H \endlink
85 
86 Note
87  - The full 'automatic wall treatment' description also requires use of the
88  \link omegaWallFunctionFvPatchScalarField.H \endlink with the \c blending
89  option \c binomial or with the deprecated \c blended flag set to \c on.
90  - Suffers from non-exact restart since \c correctNut() (called through
91  \c turbulence->validate) returns a slightly different value every time
92  it is called.
93  See \link nutUSpaldingWallFunctionFvPatchScalarField.C \endlink.
94 
95 SourceFiles
96  nutUBlendedWallFunctionFvPatchScalarField.C
97 
98 \*---------------------------------------------------------------------------*/
99 
100 #ifndef nutUBlendedWallFunctionFvPatchScalarField_H
101 #define nutUBlendedWallFunctionFvPatchScalarField_H
102 
104 
105 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
106 
107 namespace Foam
108 {
109 
110 /*---------------------------------------------------------------------------*\
111  Class nutUBlendedWallFunctionFvPatchScalarField Declaration
112 \*---------------------------------------------------------------------------*/
113 
114 class nutUBlendedWallFunctionFvPatchScalarField
115 :
116  public nutWallFunctionFvPatchScalarField
117 {
118 protected:
119 
120  // Protected Data
121 
122  //- Model coefficient; default = 4
123  scalar n_;
124 
125 
126  // Protected Member Functions
127 
128  //- Calculate the turbulent viscosity
129  virtual tmp<scalarField> calcNut() const;
130 
131  //- Calculate the friction velocity
132  tmp<scalarField> calcUTau(const scalarField& magGradU) const;
133 
134  //- Write local wall function variables
135  void writeLocalEntries(Ostream&) const;
136 
137 
138 public:
140  //- Runtime type information
141  TypeName("nutUBlendedWallFunction");
142 
143 
144  // Constructors
145 
146  //- Construct from patch and internal field
148  (
149  const fvPatch&,
151  );
152 
153  //- Construct from patch, internal field and dictionary
155  (
156  const fvPatch&,
158  const dictionary&
159  );
160 
161  //- Construct by mapping given
162  //- nutUBlendedWallFunctionFvPatchScalarField
163  //- onto a new patch
165  (
167  const fvPatch&,
169  const fvPatchFieldMapper&
170  );
171 
172  //- Construct as copy
174  (
176  );
177 
178  //- Construct and return a clone
179  virtual tmp<fvPatchScalarField> clone() const
180  {
182  (
184  );
185  }
186 
187  //- Construct as copy setting internal field reference
189  (
192  );
193 
194  //- Construct and return a clone setting internal field reference
196  (
198  ) const
199  {
201  (
203  );
204  }
205 
206 
207  // Member Functions
208 
209  // Evaluation
210 
211  //- Calculate and return the yPlus at the boundary
212  virtual tmp<scalarField> yPlus() const;
213 
214 
215  // I-O
216 
217  //- Write
218  virtual void write(Ostream& os) const;
219 };
220 
221 
222 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223 
224 } // End namespace Foam
225 
226 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
227 
228 #endif
229 
230 // ************************************************************************* //
virtual tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
nutUBlendedWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
This boundary condition provides a wall function for the turbulent viscosity (i.e. nut) based on velocity (i.e. U) using a binomial-function wall-function blending method between the viscous and inertial sublayer predictions of nut for low- and high-Reynolds number applications.
tmp< scalarField > calcUTau(const scalarField &magGradU) const
Calculate the friction velocity.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
virtual tmp< fvPatchScalarField > clone() const
Construct and return a clone.
TypeName("nutUBlendedWallFunction")
Runtime type information.
OBJstream os(runTime.globalPath()/outputName)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Namespace for OpenFOAM.
void writeLocalEntries(Ostream &) const
Write local wall function variables.