nutkWallFunctionFvPatchScalarField.C
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) 2011-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 \*---------------------------------------------------------------------------*/
25 
27 #include "turbulenceModel.H"
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "wallFvPatch.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
39 
41 {
42  const label patchi = patch().index();
43 
44  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
45  (
47  (
49  internalField().group()
50  )
51  );
52 
53  const scalarField& y = turbModel.y()[patchi];
54  const tmp<volScalarField> tk = turbModel.k();
55  const volScalarField& k = tk();
56  const tmp<scalarField> tnuw = turbModel.nu(patchi);
57  const scalarField& nuw = tnuw();
58 
59  const scalar Cmu25 = pow025(Cmu_);
60 
61  tmp<scalarField> tnutw(new scalarField(patch().size(), 0.0));
62  scalarField& nutw = tnutw.ref();
63 
64  forAll(nutw, facei)
65  {
66  label celli = patch().faceCells()[facei];
67 
68  scalar yPlus = Cmu25*y[facei]*sqrt(k[celli])/nuw[facei];
69 
70  if (yPlus > yPlusLam_)
71  {
72  nutw[facei] = nuw[facei]*(yPlus*kappa_/log(E_*yPlus) - 1.0);
73  }
74  }
75 
76  return tnutw;
77 }
78 
79 
80 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
81 
83 (
84  const fvPatch& p,
86 )
87 :
89 {}
90 
91 
93 (
95  const fvPatch& p,
97  const fvPatchFieldMapper& mapper
98 )
99 :
100  nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper)
101 {}
102 
103 
105 (
106  const fvPatch& p,
108  const dictionary& dict
109 )
110 :
112 {}
113 
114 
116 (
118 )
119 :
121 {}
122 
123 
125 (
128 )
129 :
131 {}
132 
133 
134 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
135 
137 {
138  const label patchi = patch().index();
139 
140  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
141  (
143  (
145  internalField().group()
146  )
147  );
148 
149  const scalarField& y = turbModel.y()[patchi];
150 
151  const tmp<volScalarField> tk = turbModel.k();
152  const volScalarField& k = tk();
153  tmp<scalarField> kwc = k.boundaryField()[patchi].patchInternalField();
154  const tmp<scalarField> tnuw = turbModel.nu(patchi);
155  const scalarField& nuw = tnuw();
156 
157  return pow025(Cmu_)*y*sqrt(kwc)/nuw;
158 }
159 
160 
161 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
162 
164 (
167 );
168 
169 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
170 
171 } // End namespace Foam
172 
173 // ************************************************************************* //
const char *const group
Group name for atomic constants.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
virtual tmp< scalarField > calcNut() const
Calculate the turbulence viscosity.
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
dimensionedScalar log(const dimensionedScalar &ds)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
dimensionedScalar pow025(const dimensionedScalar &ds)
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
label k
Boltzmann constant.
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
scalar y
Turbulent viscosity wall-function boundary condition for high Reynolds number flows based on near-wal...
static const word propertiesName
Default name of the turbulence properties dictionary.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
const nearWallDist & y() const
Return the near wall distances.
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
label patchi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
nutWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A class for managing temporary objects.
Definition: PtrList.H:53
scalar yPlusLam_
Y+ at the edge of the laminar sublayer.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
nutkWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Namespace for OpenFOAM.