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 | Copyright (C) 2015-2016 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 Group
29  grpCmpWallFunctions
30 
31 Description
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 
36 Usage
37  \table
38  Property | Description | Required | Default value
39  Prt | Turbulent Prandtl number | no | 0.85
40  Cmu | Model coefficient | no | 0.09
41  kappa | von Karman constant | no | 0.41
42  E | Model coefficient | no | 9.8
43  \endtable
44 
45  Example of the boundary condition specification:
46  \verbatim
47  <patchName>
48  {
49  type alphatPhaseChangeJayatillekeWallFunction;
50  Prt 0.85;
51  kappa 0.41;
52  E 9.8;
53  value uniform 0; // optional value entry
54  }
55  \endverbatim
56 
57 See also
58  Foam::compressible::alphatPhaseChangeWallFunctionFvPatchScalarField
59 
60 SourceFiles
61  alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.C
62 
63 \*---------------------------------------------------------------------------*/
64 
65 #ifndef compressiblealphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField_H
66 #define compressiblealphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField_H
67 
69 
70 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71 
72 namespace Foam
73 {
74 namespace compressible
75 {
76 
77 /*---------------------------------------------------------------------------*\
78  Class alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField Declaration
79 \*---------------------------------------------------------------------------*/
80 
81 class alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
82 :
83  public alphatPhaseChangeWallFunctionFvPatchScalarField
84 {
85 
86 protected:
87 
88  // Protected data
89 
90  //- Turbulent Prandtl number
91  scalar Prt_;
92 
93  //- Cmu coefficient
94  scalar Cmu_;
95 
96  //- Von Karman constant
97  scalar kappa_;
98 
99  //- E coefficient
100  scalar E_;
101 
102  // Solution parameters
103 
104  static scalar maxExp_;
105  static scalar tolerance_;
106  static label maxIters_;
107 
108 
109  // Protected Member Functions
110 
111  //- Check the type of the patch
112  void checkType();
113 
114  //- 'P' function
115  tmp<scalarField> Psmooth(const scalarField& Prat) const;
116 
117  //- Calculate y+ at the edge of the thermal laminar sublayer
119  (
120  const scalarField& P,
121  const scalarField& Prat
122  ) const;
124  //- Update turbulent thermal diffusivity
126  (
127  const scalarField& prevAlphat
128  ) const;
130 
131 public:
132 
133  //- Runtime type information
134  TypeName("compressible::alphatPhaseChangeJayatillekeWallFunction");
135 
136 
137  // Constructors
138 
139  //- Construct from patch and internal field
141  (
142  const fvPatch&,
144  );
145 
146  //- Construct from patch, internal field and dictionary
148  (
149  const fvPatch&,
151  const dictionary&
152  );
153 
154  //- Construct by mapping given
155  // alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
156  // onto a new patch
158  (
160  const fvPatch&,
162  const fvPatchFieldMapper&
163  );
164 
165  //- Construct as copy
167  (
169  );
170 
171  //- Construct and return a clone
172  virtual tmp<fvPatchScalarField> clone() const
173  {
175  (
177  (
178  *this
179  )
180  );
181  }
182 
183  //- Construct as copy setting internal field reference
185  (
188  );
189 
190  //- Construct and return a clone setting internal field reference
192  (
194  ) const
195  {
197  (
199  (
200  *this,
201  iF
202  )
203  );
204  }
205 
206 
207  // Member functions
208 
209  // Evaluation functions
210 
211  //- Update the coefficients associated with the patch field
212  virtual void updateCoeffs();
213 
214 
215  // I-O
216 
217  //- Write
218  virtual void write(Ostream&) const;
219 };
220 
221 
222 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223 
224 } // End namespace compressible
225 } // End namespace Foam
226 
227 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
228 
229 #endif
230 
231 // ************************************************************************* //
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
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
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
tmp< scalarField > calcAlphat(const scalarField &prevAlphat) const
Update turbulent thermal diffusivity.
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...
tmp< scalarField > yPlusTherm(const scalarField &P, const scalarField &Prat) const
Calculate y+ at the edge of the thermal laminar sublayer.
A class for managing temporary objects.
Definition: PtrList.H:54
tmp< scalarField > Psmooth(const scalarField &Prat) const
&#39;P&#39; function
Namespace for OpenFOAM.
bool compressible
Definition: pEqn.H:30