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::compressible::alphatJayatillekeWallFunctionFvPatchScalarField
29
30Group
31 grpCmpWallFunctions
32
33Description
34 This boundary condition provides a thermal wall function for turbulent
35 thermal diffusivity (usually\c alphat) based on the Jayatilleke model.
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
67SourceFiles
68 alphatJayatillekeWallFunctionFvPatchScalarField.C
69
70\*---------------------------------------------------------------------------*/
71
72#ifndef compressible_alphatJayatillekeWallFunctionFvPatchScalarField_H
73#define compressible_alphatJayatillekeWallFunctionFvPatchScalarField_H
74
77
78// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
79
80namespace Foam
81{
82namespace compressible
83{
84
85/*---------------------------------------------------------------------------*\
86 Class alphatJayatillekeWallFunctionFvPatchScalarField Declaration
87\*---------------------------------------------------------------------------*/
88
90:
91 public fixedValueFvPatchScalarField
92{
93 // Private Data
94
95 //- Turbulent Prandtl number
96 scalar Prt_;
97
98 //- Von Karman constant
99 scalar kappa_;
100
101 //- E coefficient
102 scalar E_;
103
104
105 // Solution parameters
106
107 static scalar maxExp_;
108 static scalar tolerance_;
109 static label maxIters_;
110
111
112 // Private Member Functions
113
114 //- Check the type of the patch
115 void checkType();
116
117 //- Return the patch y+
119 (
120 const compressible::turbulenceModel& turbModel
121 ) const;
122
123 //- `P' function
124 scalar Psmooth(const scalar Prat) const;
125
126 //- Calculate y+ at the edge of the thermal laminar sublayer
127 scalar yPlusTherm
128 (
129 const scalar P,
130 const scalar Prat
131 ) const;
132
133
134public:
135
136 //- Runtime type information
137 TypeName("compressible::alphatJayatillekeWallFunction");
138
139
140 // Constructors
141
142 //- Construct from patch and internal field
144 (
145 const fvPatch&,
147 );
148
149 //- Construct from patch, internal field and dictionary
151 (
152 const fvPatch&,
154 const dictionary&
155 );
156
157 //- Construct by mapping given an
158 //- alphatJayatillekeWallFunctionFvPatchScalarField
159 //- onto a new patch
161 (
163 const fvPatch&,
165 const fvPatchFieldMapper&
166 );
167
168 //- Construct as copy
170 (
172 );
173
174 //- Construct as copy setting internal field reference
176 (
179 );
180
181 //- Return a clone
182 virtual tmp<fvPatchField<scalar>> clone() const
183 {
184 return fvPatchField<scalar>::Clone(*this);
185 }
186
187 //- Clone with an internal field reference
189 (
191 ) const
192 {
193 return fvPatchField<scalar>::Clone(*this, iF);
194 }
195
196
197 // Member Functions
198
199 // Evaluation
200
201 //- Update the coefficients associated with the patch field
202 virtual void updateCoeffs();
203
204
205 // I-O
206
207 //- Write
208 void write(Ostream&) const;
209};
210
211
212// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213
214} // End namespace compressible
215} // End namespace Foam
216
217// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
218
219#endif
220
221// ************************************************************************* //
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
This boundary condition provides a thermal wall function for turbulent thermal diffusivity (usuallyal...
TypeName("compressible::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.
alphatJayatillekeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
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
A class for managing temporary objects.
Definition tmp.H:75
scalar yPlus
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Namespace for OpenFOAM.
runTime write()
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68