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-2018 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 
66 
67 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68 
69 namespace Foam
70 {
71 namespace compressible
72 {
73 
74 /*---------------------------------------------------------------------------*\
75  Class alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField Declaration
76 \*---------------------------------------------------------------------------*/
77 
78 class alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
79 :
80  public alphatPhaseChangeWallFunctionFvPatchScalarField
81 {
82 
83 protected:
84 
85  // Protected data
86 
87  //- Turbulent Prandtl number
88  scalar Prt_;
89 
90  //- Cmu coefficient
91  scalar Cmu_;
92 
93  //- Von Karman constant
94  scalar kappa_;
95 
96  //- E coefficient
97  scalar E_;
98 
99  // Solution parameters
100 
101  static scalar maxExp_;
102  static scalar tolerance_;
103  static label maxIters_;
104 
105 
106  // Protected Member Functions
107 
108  //- Check the type of the patch
109  void checkType();
110 
111  //- 'P' function
112  tmp<scalarField> Psmooth(const scalarField& Prat) const;
113 
114  //- Calculate y+ at the edge of the thermal laminar sublayer
116  (
117  const scalarField& P,
118  const scalarField& Prat
119  ) const;
121  //- Update turbulent thermal diffusivity
123  (
124  const scalarField& prevAlphat
125  ) const;
127 
128 public:
129 
130  //- Runtime type information
131  TypeName("compressible::alphatPhaseChangeJayatillekeWallFunction");
132 
133 
134  // Constructors
135 
136  //- Construct from patch and internal field
138  (
139  const fvPatch&,
141  );
142 
143  //- Construct from patch, internal field and dictionary
145  (
146  const fvPatch&,
148  const dictionary&
149  );
150 
151  //- Construct by mapping given
152  // alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
153  // onto a new patch
155  (
157  const fvPatch&,
159  const fvPatchFieldMapper&
160  );
161 
162  //- Construct as copy
164  (
166  );
167 
168  //- Construct and return a clone
169  virtual tmp<fvPatchScalarField> clone() const
170  {
172  (
174  (
175  *this
176  )
177  );
178  }
179 
180  //- Construct as copy setting internal field reference
182  (
185  );
186 
187  //- Construct and return a clone setting internal field reference
189  (
191  ) const
192  {
194  (
196  (
197  *this,
198  iF
199  )
200  );
201  }
202 
203 
204  // Member functions
205 
206  // Evaluation functions
207 
208  //- Update the coefficients associated with the patch field
209  virtual void updateCoeffs();
210 
211 
212  // I-O
213 
214  //- Write
215  virtual void write(Ostream&) const;
216 };
217 
218 
219 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
220 
221 } // End namespace compressible
222 } // End namespace Foam
223 
224 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
225 
226 #endif
227 
228 // ************************************************************************* //
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:137
tmp< scalarField > calcAlphat(const scalarField &prevAlphat) const
Update turbulent thermal diffusivity.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
tmp< scalarField > yPlusTherm(const scalarField &P, const scalarField &Prat) const
Calculate y+ at the edge of the thermal laminar sublayer.
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.