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-2023 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().lookupType<momentumTransportModel>(internalField().group());
46 
47  const scalarField& y = turbModel.y()[patchi];
48  const tmp<volScalarField> tk = turbModel.k();
49  const volScalarField& k = tk();
50  const tmp<scalarField> tnuw = turbModel.nu(patchi);
51  const scalarField& nuw = tnuw();
52 
53  const scalar Cmu25 = pow025(Cmu_);
54 
55  tmp<scalarField> tnutw(new scalarField(patch().size(), 0.0));
56  scalarField& nutw = tnutw.ref();
57 
58  forAll(nutw, facei)
59  {
60  label celli = patch().faceCells()[facei];
61 
62  scalar yPlus = Cmu25*y[facei]*sqrt(k[celli])/nuw[facei];
63 
64  if (yPlus > yPlusLam_)
65  {
66  nutw[facei] = nuw[facei]*(yPlus*kappa_/log(E_*yPlus) - 1.0);
67  }
68  }
69 
70  return tnutw;
71 }
72 
73 
74 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
75 
77 (
78  const fvPatch& p,
80  const dictionary& dict
81 )
82 :
84 {}
85 
86 
88 (
90  const fvPatch& p,
92  const fvPatchFieldMapper& mapper
93 )
94 :
95  nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper)
96 {}
97 
98 
100 (
103 )
104 :
106 {}
107 
108 
109 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
110 
112 {
113  const label patchi = patch().index();
114 
115  const momentumTransportModel& turbModel =
116  db().lookupType<momentumTransportModel>(internalField().group());
117 
118  const scalarField& y = turbModel.y()[patchi];
119 
120  const tmp<volScalarField> tk = turbModel.k();
121  const volScalarField& k = tk();
122  tmp<scalarField> kwc = k.boundaryField()[patchi].patchInternalField();
123  const tmp<scalarField> tnuw = turbModel.nu(patchi);
124  const scalarField& nuw = tnuw();
125 
126  return pow025(Cmu_)*y*sqrt(kwc)/nuw;
127 }
128 
129 
130 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
131 
133 (
136 );
137 
138 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139 
140 } // End namespace Foam
141 
142 // ************************************************************************* //
scalar y
label k
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
Abstract base class for turbulence models (RAS, LES and laminar).
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
const volScalarField::Boundary & y() const
Return the near wall distance.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions,...
scalar yPlusLam_
Y+ at the edge of the laminar sublayer.
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
Turbulent viscosity wall-function boundary condition for high Reynolds number flows based on near-wal...
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
virtual tmp< scalarField > nut() const
Calculate the turbulence viscosity.
nutkWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
label patchi
Namespace for OpenFOAM.
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
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar sqrt(const dimensionedScalar &ds)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dimensionedScalar pow025(const dimensionedScalar &ds)
dictionary dict
volScalarField & p