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-2021 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 "momentumTransportModel.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 momentumTransportModel& turbModel =
45  db().lookupObject<momentumTransportModel>
46  (
48  (
49  momentumTransportModel::typeName,
50  internalField().group()
51  )
52  );
53 
54  const scalarField& y = turbModel.y()[patchi];
55  const tmp<volScalarField> tk = turbModel.k();
56  const volScalarField& k = tk();
57  const tmp<scalarField> tnuw = turbModel.nu(patchi);
58  const scalarField& nuw = tnuw();
59 
60  const scalar Cmu25 = pow025(Cmu_);
61 
62  tmp<scalarField> tnutw(new scalarField(patch().size(), 0.0));
63  scalarField& nutw = tnutw.ref();
64 
65  forAll(nutw, facei)
66  {
67  label celli = patch().faceCells()[facei];
68 
69  scalar yPlus = Cmu25*y[facei]*sqrt(k[celli])/nuw[facei];
70 
71  if (yPlus > yPlusLam_)
72  {
73  nutw[facei] = nuw[facei]*(yPlus*kappa_/log(E_*yPlus) - 1.0);
74  }
75  }
76 
77  return tnutw;
78 }
79 
80 
81 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
82 
84 (
85  const fvPatch& p,
87 )
88 :
90 {}
91 
92 
94 (
95  const fvPatch& p,
97  const dictionary& dict
98 )
99 :
101 {}
102 
103 
105 (
107  const fvPatch& p,
109  const fvPatchFieldMapper& mapper
110 )
111 :
112  nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper)
113 {}
114 
115 
117 (
120 )
121 :
123 {}
124 
125 
126 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127 
129 {
130  const label patchi = patch().index();
131 
132  const momentumTransportModel& turbModel =
133  db().lookupObject<momentumTransportModel>
134  (
136  (
137  momentumTransportModel::typeName,
138  internalField().group()
139  )
140  );
141 
142  const scalarField& y = turbModel.y()[patchi];
143 
144  const tmp<volScalarField> tk = turbModel.k();
145  const volScalarField& k = tk();
146  tmp<scalarField> kwc = k.boundaryField()[patchi].patchInternalField();
147  const tmp<scalarField> tnuw = turbModel.nu(patchi);
148  const scalarField& nuw = tnuw();
149 
150  return pow025(Cmu_)*y*sqrt(kwc)/nuw;
151 }
152 
153 
154 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
155 
157 (
160 );
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
164 } // End namespace Foam
165 
166 // ************************************************************************* //
const char *const group
Group name for atomic constants.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
dimensionedScalar log(const dimensionedScalar &ds)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
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:181
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
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.
const volScalarField::Boundary & y() const
Return the near wall distance.
Macros for easy insertion into run-time selection tables.
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
scalar y
Turbulent viscosity wall-function boundary condition for high Reynolds number flows based on near-wal...
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Abstract base class for turbulence models (RAS, LES and laminar).
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)
virtual tmp< scalarField > nut() const
Calculate the turbulence viscosity.
nutkWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Namespace for OpenFOAM.