Loading...
Searching...
No Matches
alphatJayatillekeWallFunctionFvPatchScalarField.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2017-2022 OpenCFD Ltd
10-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField
29
30Group
31 grpIcoWallFunctions
32
33Description
34 This boundary condition provides a kinematic turbulent thermal conductivity
35 for using wall functions, using the Jayatilleke 'P' function.
36
37Usage
38 Example of the boundary condition specification:
39 \verbatim
40 <patchName>
41 {
42 // Mandatory entries
43 type alphatJayatillekeWallFunction;
44
45 // Optional entries
46 Prt <scalar>;
47 kappa <scalar>;
48 E <scalar>;
49
50 // Inherited entries
51 ...
52 }
53 \endverbatim
54
55 where the entries mean:
56 \table
57 Property | Description | Type | Reqd | Deflt
58 type | Type name: alphatJayatillekeWallFunction | word | yes | -
59 Prt | Turbulent Prandtl number | scalar | no | 0.85
60 kappa | von Karman constant | scalar | no | 0.41
61 E | Wall roughness parameter | scalar | no | 9.8
62 \endtable
63
64 The inherited entries are elaborated in:
65 - \link fixedValueFvPatchFields.H \endlink
66
67Note
68 - The units of kinematic turbulent thermal conductivity are [m2/s].
69
70SourceFiles
71 alphatJayatillekeWallFunctionFvPatchScalarField.C
72
73\*---------------------------------------------------------------------------*/
74
75#ifndef incompressible_alphatJayatillekeWallFunctionFvPatchScalarField_H
76#define incompressible_alphatJayatillekeWallFunctionFvPatchScalarField_H
77
79#include "turbulenceModel.H"
80
81// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
82
83namespace Foam
84{
85namespace incompressible
86{
87
88/*---------------------------------------------------------------------------*\
89 Class alphatJayatillekeWallFunctionFvPatchScalarField Declaration
90\*---------------------------------------------------------------------------*/
91
93:
94 public fixedValueFvPatchScalarField
95{
96protected:
97
98 // Protected Data
99
100 //- Turbulent Prandtl number
101 scalar Prt_;
102
103 //- Von Karman constant
104 scalar kappa_;
105
106 //- E coefficient
107 scalar E_;
108
109
110 // Solution parameters
111
112 static scalar tolerance_;
113 static label maxIters_;
114
115
116 // Protected Member Functions
118 //- Check the type of the patch
119 virtual void checkType();
120
121 //- Return the patch y+
122 tmp<scalarField> yPlus(const turbulenceModel& turbModel) const;
123
124 //- `P' function
125 scalar Psmooth(const scalar Prat) const;
126
127 //- Calculate y+ at the edge of the thermal laminar sublayer
129 (
130 const scalar P,
131 const scalar Prat
132 ) const;
134
135public:
136
137 //- Runtime type information
138 TypeName("alphatJayatillekeWallFunction");
139
140
141 // Constructors
142
143 //- Construct from patch and internal field
145 (
146 const fvPatch&,
148 );
149
150 //- Construct from patch, internal field and dictionary
152 (
153 const fvPatch&,
155 const dictionary&
156 );
157
158 //- Construct by mapping given
159 //- alphatJayatillekeWallFunctionFvPatchScalarField
160 //- onto a new patch
162 (
164 const fvPatch&,
166 const fvPatchFieldMapper&
167 );
168
169 //- Construct as copy
171 (
173 );
174
175 //- Construct as copy setting internal field reference
177 (
180 );
181
182 //- Return a clone
183 virtual tmp<fvPatchField<scalar>> clone() const
184 {
185 return fvPatchField<scalar>::Clone(*this);
186 }
187
188 //- Clone with an internal field reference
190 (
192 ) const
193 {
194 return fvPatchField<scalar>::Clone(*this, iF);
195 }
196
197
198 // Member Functions
199
200 // Evaluation
201
202 //- Update the coefficients associated with the patch field
203 virtual void updateCoeffs();
204
205
206 // I-O
207
208 //- Write
209 virtual void write(Ostream&) const;
210};
211
212
213// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
214
215} // End namespace incompressible
216} // End namespace Foam
217
218// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219
220#endif
221
222// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A FieldMapper for finite-volume patch fields.
static tmp< fvPatchField< Type > > Clone(const DerivedPatchField &pf, Args &&... args)
Clone a patch field, optionally with internal field reference etc.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
This boundary condition provides a kinematic turbulent thermal conductivity for using wall functions,...
TypeName("alphatJayatillekeWallFunction")
Runtime type information.
virtual tmp< fvPatchField< scalar > > clone(const DimensionedField< scalar, volMesh > &iF) const
Clone with an internal field reference.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
scalar yPlusTherm(const scalar P, const scalar Prat) const
`P' function */ scalar Psmooth(const scalar Prat) const;
alphatJayatillekeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A class for managing temporary objects.
Definition tmp.H:75
scalar yPlus
IncompressibleTurbulenceModel< transportModel > turbulenceModel
Namespace for OpenFOAM.
runTime write()
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68