alphatPhaseJayatillekeWallFunctionFvPatchScalarField.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-2020 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  alphatPhaseJayatillekeWallFunctionFvPatchScalarField
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 alphatPhaseJayatillekeWallFunction;
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  alphatPhaseJayatillekeWallFunctionFvPatchScalarField.C
59 
60 \*---------------------------------------------------------------------------*/
61 
62 #ifndef compressible_alphatPhaseJayatillekeWallFunctionFvPatchScalarField_H
63 #define compressible_alphatPhaseJayatillekeWallFunctionFvPatchScalarField_H
64 
67 
68 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69 
70 namespace Foam
71 {
72 namespace compressible
73 {
74 
75 /*---------------------------------------------------------------------------*\
76  Class alphatPhaseJayatillekeWallFunctionFvPatchScalarField Declaration
77 \*---------------------------------------------------------------------------*/
78 
79 class alphatPhaseJayatillekeWallFunctionFvPatchScalarField
80 :
81  public fixedValueFvPatchScalarField
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 
112 public:
113 
114  //- Runtime type information
115  TypeName("compressible::alphatPhaseJayatillekeWallFunction");
118  // Constructors
119 
120  //- Construct from patch and internal field
122  (
123  const fvPatch&,
125  );
126 
127  //- Construct from patch, internal field and dictionary
129  (
130  const fvPatch&,
132  const dictionary&
133  );
134 
135  //- Construct by mapping given
136  // alphatPhaseJayatillekeWallFunctionFvPatchScalarField
137  // onto a new patch
139  (
141  const fvPatch&,
143  const fvPatchFieldMapper&
144  );
145 
146  //- Disallow copy without setting internal field reference
148  (
150  ) = delete;
151 
152  //- Copy constructor setting internal field reference
154  (
157  );
158 
159  //- Construct and return a clone setting internal field reference
161  (
163  ) const
164  {
166  (
168  (
169  *this,
170  iF
171  )
172  );
173  }
174 
175 
176  // Member Functions
177 
178  // Evaluation functions
179 
180  //- Evaluate the turbulent thermal diffusivity
181  tmp<scalarField> calcAlphat(const scalarField& prevAlphat) const;
182 
183  //- Update the coefficients associated with the patch field
184  virtual void updateCoeffs();
185 
186 
187  // I-O
188 
189  //- Write
190  virtual void write(Ostream&) const;
191 };
192 
193 
194 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195 
196 } // End namespace compressible
197 } // End namespace Foam
198 
199 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200 
201 #endif
202 
203 // ************************************************************************* //
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual tmp< fvPatchScalarField > clone(const DimensionedField< scalar, volMesh > &iF) const
Construct and return a clone setting internal field reference.
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
tmp< scalarField > calcAlphat(const scalarField &prevAlphat) const
Evaluate the turbulent thermal diffusivity.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
TypeName("compressible::alphatPhaseJayatillekeWallFunction")
Runtime type information.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions...
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Foam::fvPatchFieldMapper.
tmp< scalarField > yPlusTherm(const nutWallFunctionFvPatchScalarField &nutw, const scalarField &P, const scalarField &Prat) const
Calculate y+ at the edge of the thermal laminar sublayer.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
tmp< scalarField > Psmooth(const scalarField &Prat) const
&#39;P&#39; function
alphatPhaseJayatillekeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
This boundary condition provides a thermal wall function for turbulent thermal diffusivity (usuallyal...
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.