Loading...
Searching...
No Matches
alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.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) 2015-2018 OpenFOAM Foundation
9 Copyright (C) 2021 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::
29 alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
30
31Description
32 This boundary condition provides a thermal wall function for turbulent
33 thermal diffusivity (usually\c alphat) based on the Jayatilleke model for
34 the Eulerian multiphase solvers.
35
36Usage
37 Example of the boundary condition specification:
38 \verbatim
39 <patchName>
40 {
41 // Mandatory entries
42 type alphatPhaseChangeJayatillekeWallFunction;
43
44 // Optional entries
45 Prt <scalar>;
46 Cmu <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: <!--
59 --> compressible::alphatPhaseChangeJayatillekeWallFunction <!--
60 --> | word | yes | -
61 Prt | Turbulent Prandtl number | scalar | no | 0.85
62 Cmu | Empirical model coefficient | scalar | no | 0.09
63 kappa | Von Karman constant | scalar | no | 0.41
64 E | Wall roughness parameter | scalar | no | 9.8
65 \endtable
66
67 The inherited entries are elaborated in:
68 - \link alphatPhaseChangeWallFunctionFvPatchScalarField.H \endlink
69
70See also
71 Foam::compressible::alphatPhaseChangeWallFunctionFvPatchScalarField
72
73SourceFiles
74 alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.C
75
76\*---------------------------------------------------------------------------*/
77
78#ifndef compressible_alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField_H
79#define compressible_alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField_H
80
82
83// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
84
85namespace Foam
86{
87namespace compressible
88{
89
90/*---------------------------------------------------------------------------*\
91 Class alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField Declaration
92\*---------------------------------------------------------------------------*/
93
95:
97{
98protected:
99
100 // Protected Data
101
102 //- Turbulent Prandtl number
103 scalar Prt_;
104
105 //- Empirical model coefficient
106 scalar Cmu_;
107
108 //- Von Karman constant
109 scalar kappa_;
110
111 //- Wall roughness parameter
112 scalar E_;
113
114 // Solution parameters
115
116 //- Absolute tolerance
117 static scalar tolerance_;
118
119 //- Maximum number of iterations
120 static label maxIters_;
121
122
123 // Protected Member Functions
124
125 //- Check the type of the patch
126 void checkType();
127
128 //- 'P' function
129 tmp<scalarField> Psmooth(const scalarField& Prat) const;
130
131 //- Calculate y+ at the edge of the thermal laminar sublayer
133 (
134 const scalarField& P,
135 const scalarField& Prat
136 ) const;
137
138 //- Update turbulent thermal diffusivity
140 (
141 const scalarField& prevAlphat
142 ) const;
143
145public:
146
147 //- Runtime type information
148 TypeName("compressible::alphatPhaseChangeJayatillekeWallFunction");
150
151 // Constructors
152
153 //- Construct from patch and internal field
155 (
156 const fvPatch&,
158 );
159
160 //- Construct from patch, internal field and dictionary
162 (
163 const fvPatch&,
165 const dictionary&
166 );
167
168 //- Construct by mapping given
169 //- alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
170 //- onto a new patch
172 (
174 const fvPatch&,
176 const fvPatchFieldMapper&
177 );
178
179 //- Construct as copy
181 (
183 );
184
185 //- Construct as copy setting internal field reference
187 (
190 );
191
192 //- Return a clone
193 virtual tmp<fvPatchField<scalar>> clone() const
194 {
195 return fvPatchField<scalar>::Clone(*this);
196 }
197
198 //- Clone with an internal field reference
200 (
202 ) const
203 {
205 }
206
207
208 // Member Functions
209
210 // Evaluation
211
212 //- Update the coefficients associated with the patch field
213 virtual void updateCoeffs();
214
215
216 // I-O
217
218 //- Write
219 virtual void write(Ostream&) const;
220};
221
222
223// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
224
225} // End namespace compressible
226} // End namespace Foam
227
228// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229
230#endif
231
232// ************************************************************************* //
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...
virtual tmp< fvPatchField< scalar > > clone(const DimensionedField< scalar, volMesh > &iF) const
Clone with an internal field reference.
tmp< scalarField > calcAlphat(const scalarField &prevAlphat) const
Update turbulent thermal diffusivity.
tmp< scalarField > yPlusTherm(const scalarField &P, const scalarField &Prat) const
Calculate y+ at the edge of the thermal laminar sublayer.
TypeName("compressible::alphatPhaseChangeJayatillekeWallFunction")
Runtime type information.
alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Abstract base-class for all alphatWallFunctions supporting phase-change.
alphatPhaseChangeWallFunctionFvPatchScalarField(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
Namespace for OpenFOAM.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
runTime write()
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68