nutkFilmWallFunctionFvPatchScalarField.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField
28 
29 Group
30  grpSurfaceFilmBoundaryConditions grpCmpWallFunctions
31 
32 Description
33  This boundary condition provides a turbulent viscosity condition when
34  using wall functions, based on turbulence kinetic energy, for use with
35  surface film models.
36 
37 Usage
38  Example of the boundary condition specification:
39  \verbatim
40  <patchName>
41  {
42  type nutkFilmWallFunction;
43  value uniform 0;
44  }
45  \endverbatim
46 
47 See also
48  Foam::nutkWallFunctionFvPatchScalarField
49 
50 SourceFiles
51  nutkFilmWallFunctionFvPatchScalarField.C
52 
53 \*---------------------------------------------------------------------------*/
54 
55 #ifndef nutkFilmWallFunctionFvPatchScalarField_H
56 #define nutkFilmWallFunctionFvPatchScalarField_H
57 
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62 namespace Foam
63 {
64 namespace compressible
65 {
66 namespace RASModels
67 {
68 
69 /*---------------------------------------------------------------------------*\
70  Class nutkFilmWallFunctionFvPatchScalarField Declaration
71 \*---------------------------------------------------------------------------*/
72 
74 :
76 {
77 protected:
78 
79  // Protected data
80 
81  //- Name of film region
83 
84  //- B Coefficient (default = 5.5)
85  scalar B_;
86 
87  //- y+ value for laminar -> turbulent transition (default = 11.05)
88  scalar yPlusCrit_;
89 
90 
91  // Protected member functions
92 
93  //- Calculate the turbulence viscosity
94  virtual tmp<scalarField> calcNut() const;
95 
96  //- Calculate the friction velocity
97  virtual tmp<scalarField> calcUTau(const scalarField& magGradU) const;
98 
99 
100 public:
101 
102  //- Runtime type information
103  TypeName("nutkFilmWallFunction");
104 
105 
106  // Constructors
107 
108  //- Construct from patch and internal field
110  (
111  const fvPatch&,
113  );
114 
115  //- Construct from patch, internal field and dictionary
117  (
118  const fvPatch&,
120  const dictionary&
121  );
122 
123  //- Construct by mapping given
124  // nutkFilmWallFunctionFvPatchScalarField
125  // onto a new patch
127  (
129  const fvPatch&,
131  const fvPatchFieldMapper&
132  );
133 
134  //- Construct as copy
136  (
138  );
139 
140  //- Construct and return a clone
141  virtual tmp<fvPatchScalarField> clone() const
142  {
144  (
146  );
147  }
148 
149  //- Construct as copy setting internal field reference
151  (
154  );
155 
156  //- Construct and return a clone setting internal field reference
158  (
160  ) const
161  {
163  (
165  );
166  }
167 
168 
169  // Member functions
170 
171  // Evaluation functions
172 
173  //- Calculate and return the yPlus at the boundary
174  virtual tmp<scalarField> yPlus() const;
175 
176 
177  // I-O
178 
179  //- Write
180  virtual void write(Ostream& os) const;
181 };
182 
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
186 } // End namespace RASModels
187 } // End namespace compressible
188 } // End namespace Foam
189 
190 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
191 
192 #endif
193 
194 // ************************************************************************* //
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::write
virtual void write(Ostream &os) const
Write.
Definition: nutkFilmWallFunctionFvPatchScalarField.C:249
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::yPlusCrit_
scalar yPlusCrit_
y+ value for laminar -> turbulent transition (default = 11.05)
Definition: nutkFilmWallFunctionFvPatchScalarField.H:87
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::yPlus
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
Definition: nutkFilmWallFunctionFvPatchScalarField.C:227
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::calcNut
virtual tmp< scalarField > calcNut() const
Calculate the turbulence viscosity.
Definition: nutkFilmWallFunctionFvPatchScalarField.C:127
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::B_
scalar B_
B Coefficient (default = 5.5)
Definition: nutkFilmWallFunctionFvPatchScalarField.H:84
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::nutkFilmWallFunctionFvPatchScalarField
nutkFilmWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: nutkFilmWallFunctionFvPatchScalarField.C:156
Foam::Field< scalar >
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
compressible
bool compressible
Definition: pEqn.H:2
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::calcUTau
virtual tmp< scalarField > calcUTau(const scalarField &magGradU) const
Calculate the friction velocity.
Definition: nutkFilmWallFunctionFvPatchScalarField.C:50
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::filmRegionName_
word filmRegionName_
Name of film region.
Definition: nutkFilmWallFunctionFvPatchScalarField.H:81
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::TypeName
TypeName("nutkFilmWallFunction")
Runtime type information.
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField
This boundary condition provides a turbulent viscosity condition when using wall functions,...
Definition: nutkFilmWallFunctionFvPatchScalarField.H:72
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
nutkWallFunctionFvPatchScalarField.H
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::clone
virtual tmp< fvPatchScalarField > clone() const
Construct and return a clone.
Definition: nutkFilmWallFunctionFvPatchScalarField.H:140
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::nutkWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the turbulent viscosity, i.e....
Definition: nutkWallFunctionFvPatchScalarField.H:94