alphatJayatillekeWallFunctionFvPatchScalarField.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::compressible::alphatJayatillekeWallFunctionFvPatchScalarField
26 
27 Description
28  This boundary condition provides a thermal wall function for turbulent
29  thermal diffusivity (usually\c alphat) based on the Jayatilleke model.
30 
31 Usage
32  \table
33  Property | Description | Required | Default value
34  Prt | turbulent Prandtl number | no | 0.85
35  \endtable
36 
37  Example of the boundary condition specification:
38  \verbatim
39  <patchName>
40  {
41  type alphatJayatillekeWallFunction;
42  Prt 0.85;
43  value uniform 0;
44  }
45  \endverbatim
46 
47  Note that other model constants (i.e., Cmu, kappa and E) are obtained from
48  the corresponding turbulent viscosity boundary condition.
49 
50 See also
51  Foam::fixedValueFvPatchField
52 
53 SourceFiles
54  alphatJayatillekeWallFunctionFvPatchScalarField.C
55 
56 \*---------------------------------------------------------------------------*/
57 
58 #ifndef alphatJayatillekeWallFunctionFvPatchScalarField_H
59 #define alphatJayatillekeWallFunctionFvPatchScalarField_H
60 
63 
64 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
66 namespace Foam
67 {
68 
69 class fluidThermophysicalTransportModel;
70 
71 namespace compressible
72 {
73 
74 /*---------------------------------------------------------------------------*\
75  Class alphatJayatillekeWallFunctionFvPatchScalarField Declaration
76 \*---------------------------------------------------------------------------*/
77 
78 class alphatJayatillekeWallFunctionFvPatchScalarField
79 :
80  public fixedValueFvPatchScalarField
81 {
82  // Private Static Data
83 
84  //- Solution tolerance
85  static const scalar tolerance_;
86 
87  //- Maximum number of solution iterations
88  static const label maxIters_;
89 
90 
91  // Private Data
92 
93  //- Turbulent Prandtl number
94  scalar Prt_;
95 
96 
97 public:
98 
99  //- Runtime type information
100  TypeName("compressible::alphatJayatillekeWallFunction");
101 
102 
103  // Constructors
104 
105  //- Construct from patch, internal field and dictionary
107  (
108  const fvPatch&,
110  const dictionary&
111  );
112 
113  //- Construct by mapping given an
114  // alphatJayatillekeWallFunctionFvPatchScalarField
115  // onto a new patch
117  (
119  const fvPatch&,
121  const fieldMapper&
122  );
123 
124  //- Disallow copy without setting internal field reference
126  (
128  ) = delete;
129 
130  //- Copy constructor setting internal field reference
132  (
135  );
136 
137  //- Construct and return a clone setting internal field reference
139  (
141  ) const
142  {
144  (
146  (
147  *this,
148  iF
149  )
150  );
151  }
152 
153 
154  // Member Functions
155 
156  // Mapping functions
157 
158  //- Map the given fvPatchField onto this fvPatchField
159  virtual void map(const fvPatchScalarField&, const fieldMapper&);
160 
161 
162  // Evaluation functions
163 
164  //- Calculate the smoothing function
165  static tmp<scalarField> P(const scalarField& Prat);
166 
167  //- Calculate y+ at the edge of the thermal laminar sublayer
169  (
171  const scalarField& P,
172  const scalarField& Prat
173  );
174 
175  //- Calculate the turbulent thermal diffusivity
176  static tmp<scalarField> alphat
177  (
179  const scalar Prt,
180  const label patchi
181  );
182 
183  //- Update the coefficients associated with the patch field
184  virtual void updateCoeffs();
185 
186 
187  // I-O
188 
189  //- Write
190  void write(Ostream&) const;
191 };
192 
193 
194 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195 
196 } // End namespace compressible
197 } // End namespace Foam
198 
199 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200 
201 #endif
202 
203 // ************************************************************************* //
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:57
This boundary condition provides a thermal wall function for turbulent thermal diffusivity (usuallyal...
TypeName("compressible::alphatJayatillekeWallFunction")
Runtime type information.
static tmp< scalarField > alphat(const fluidThermophysicalTransportModel &ttm, const scalar Prt, const label patchi)
Calculate the turbulent thermal diffusivity.
virtual tmp< fvPatchScalarField > clone(const DimensionedField< scalar, volMesh > &iF) const
Construct and return a clone setting internal field reference.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
static tmp< scalarField > yPlusTherm(const nutWallFunctionFvPatchScalarField &nutw, const scalarField &P, const scalarField &Prat)
Calculate y+ at the edge of the thermal laminar sublayer.
alphatJayatillekeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
static tmp< scalarField > P(const scalarField &Prat)
Calculate the smoothing function.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class for fluid thermophysical transport models RAS, LES and laminar.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions,...
A class for managing temporary objects.
Definition: tmp.H:55
label patchi
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59