alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.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) 2015-2019 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::
26  alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
27 
28 Description
29  This boundary condition provides a thermal wall function for turbulent
30  thermal diffusivity (usually\c alphat) based on the Jayatilleke model for
31  the Eulerian multiphase solvers.
32 
33 Usage
34  \table
35  Property | Description | Required | Default value
36  Prt | Turbulent Prandtl number | no | 0.85
37  Cmu | Model coefficient | no | 0.09
38  kappa | von Karman constant | no | 0.41
39  E | Model coefficient | no | 9.8
40  \endtable
41 
42  Example of the boundary condition specification:
43  \verbatim
44  <patchName>
45  {
46  type alphatPhaseChangeJayatillekeWallFunction;
47  Prt 0.85;
48  kappa 0.41;
49  E 9.8;
50  value uniform 0; // optional value entry
51  }
52  \endverbatim
53 
54 See also
55  Foam::compressible::alphatPhaseChangeWallFunctionFvPatchScalarField
56 
57 SourceFiles
58  alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.C
59 
60 \*---------------------------------------------------------------------------*/
61 
62 #ifndef compressible_alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField_H
63 #define compressible_alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField_H
64 
67 
68 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69 
70 namespace Foam
71 {
72 namespace compressible
73 {
74 
75 /*---------------------------------------------------------------------------*\
76  Class alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField Declaration
77 \*---------------------------------------------------------------------------*/
78 
79 class alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
80 :
81  public alphatPhaseChangeWallFunctionFvPatchScalarField
82 {
83 
84 protected:
85 
86  // Protected data
87 
88  //- Turbulent Prandtl number
89  scalar Prt_;
90 
91  // Solution parameters
92 
93  static scalar maxExp_;
94  static scalar tolerance_;
95  static label maxIters_;
96 
97 
98  // Protected Member Functions
99 
100  //- 'P' function
101  tmp<scalarField> Psmooth(const scalarField& Prat) const;
103  //- Calculate y+ at the edge of the thermal laminar sublayer
105  (
107  const scalarField& P,
108  const scalarField& Prat
109  ) const;
110 
111  //- Update turbulent thermal diffusivity
113  (
114  const scalarField& prevAlphat
115  ) const;
118 public:
119 
120  //- Runtime type information
121  TypeName("compressible::alphatPhaseChangeJayatillekeWallFunction");
122 
123 
124  // Constructors
125 
126  //- Construct from patch and internal field
128  (
129  const fvPatch&,
131  );
132 
133  //- Construct from patch, internal field and dictionary
135  (
136  const fvPatch&,
138  const dictionary&
139  );
140 
141  //- Construct by mapping given
142  // alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
143  // onto a new patch
145  (
147  const fvPatch&,
149  const fvPatchFieldMapper&
150  );
151 
152  //- Copy constructor
154  (
156  );
157 
158  //- Construct and return a clone
159  virtual tmp<fvPatchScalarField> clone() const
160  {
162  (
164  (
165  *this
166  )
167  );
168  }
169 
170  //- Copy constructor setting internal field reference
172  (
175  );
176 
177  //- Construct and return a clone setting internal field reference
179  (
181  ) const
182  {
184  (
186  (
187  *this,
188  iF
189  )
190  );
191  }
192 
193 
194  // Member Functions
195 
196  // Evaluation functions
197 
198  //- Update the coefficients associated with the patch field
199  virtual void updateCoeffs();
200 
201 
202  // I-O
203 
204  //- Write
205  virtual void write(Ostream&) const;
206 };
207 
208 
209 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
210 
211 } // End namespace compressible
212 } // End namespace Foam
213 
214 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215 
216 #endif
217 
218 // ************************************************************************* //
tmp< scalarField > Psmooth(const scalarField &Prat) const
&#39;P&#39; function
TypeName("compressible::alphatPhaseChangeJayatillekeWallFunction")
Runtime type information.
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
tmp< scalarField > calcAlphat(const scalarField &prevAlphat) const
Update turbulent thermal diffusivity.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
tmp< scalarField > yPlusTherm(const nutWallFunctionFvPatchScalarField &nutw, const scalarField &P, const scalarField &Prat) const
Calculate y+ at the edge of the thermal laminar sublayer.
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions...
alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
This boundary condition provides a thermal wall function for turbulent thermal diffusivity (usuallyal...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.