nutWallFunctionFvPatchScalarField.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-2024 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 "fieldMapper.H"
29 #include "wallFvPatch.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
37 }
38 
39 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
40 
42 {
43  if (!isA<wallFvPatch>(patch()))
44  {
46  << "Invalid wall function specification" << nl
47  << " Patch type for patch " << patch().name()
48  << " must be wall" << nl
49  << " Current patch type is " << patch().type() << nl << endl
50  << abort(FatalError);
51  }
52 }
53 
54 
56 (
57  Ostream& os
58 ) const
59 {
60  writeEntry(os, "Cmu", Cmu_);
61  writeEntry(os, "kappa", kappa_);
62  writeEntry(os, "E", E_);
63 }
64 
65 
66 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
67 
69 (
70  const fvPatch& p,
72  const dictionary& dict
73 )
74 :
75  fixedValueFvPatchScalarField(p, iF, dict),
76  Cmu_(dict.lookupOrDefault<scalar>("Cmu", dimless, 0.09)),
77  kappa_(dict.lookupOrDefault<scalar>("kappa", dimless, 0.41)),
78  E_(dict.lookupOrDefault<scalar>("E", dimless, 9.8)),
79  yPlusLam_(yPlusLam(kappa_, E_))
80 {
81  checkType();
82 }
83 
84 
86 (
88  const fvPatch& p,
90  const fieldMapper& mapper
91 )
92 :
93  fixedValueFvPatchScalarField(ptf, p, iF, mapper, false),
94  Cmu_(ptf.Cmu_),
95  kappa_(ptf.kappa_),
96  E_(ptf.E_),
97  yPlusLam_(ptf.yPlusLam_)
98 {
99  checkType();
100  mapper(*this, ptf, [&](){ return this->patchInternalField(); });
101 }
102 
103 
105 (
108 )
109 :
110  fixedValueFvPatchScalarField(wfpsf, iF),
111  Cmu_(wfpsf.Cmu_),
112  kappa_(wfpsf.kappa_),
113  E_(wfpsf.E_),
114  yPlusLam_(wfpsf.yPlusLam_)
115 {
116  checkType();
117 }
118 
119 
120 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
121 
124 (
125  const momentumTransportModel& turbModel,
126  const label patchi
127 )
128 {
129  return
130  refCast<const nutWallFunctionFvPatchScalarField>
131  (
132  turbModel.nut()().boundaryField()[patchi]
133  );
134 }
135 
136 
138 (
139  const scalar kappa,
140  const scalar E
141 )
142 {
143  const scalar ypl0 = 12.0;
144 
145  scalar ypl = ypl0;
146 
147  for (int i=0; i<10; i++)
148  {
149  ypl = log(max(E*ypl, 1))/kappa;
150  }
151 
152  // If ypl is greater than the starting value recalculate from a larger value
153  // to ensure the result is slightly larger rather than slightly smaller
154  // than exact value
155  if (ypl > ypl0)
156  {
157  ypl += 1;
158 
159  for (int i=0; i<10; i++)
160  {
161  ypl = log(max(E*ypl, 1))/kappa;
162  }
163  }
164 
165  return ypl;
166 }
167 
168 
170 {
171  return yPlusLam_;
172 }
173 
174 
176 (
177  const fvPatchScalarField& ptf,
178  const fieldMapper& mapper
179 )
180 {
181  mapper(*this, ptf, [&](){ return this->patchInternalField(); });
182 }
183 
184 
186 {
187  if (updated())
188  {
189  return;
190  }
191 
192  operator==(nut());
193 
194  fixedValueFvPatchScalarField::updateCoeffs();
195 }
196 
197 
199 {
201  writeLocalEntries(os);
202  writeEntry(os, "value", *this);
203 }
204 
205 
206 // ************************************************************************* //
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...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
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 > nut() const =0
Return the turbulence viscosity.
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions,...
nutWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
scalar yPlusLam() const
Return the Y+ at the edge of the laminar sublayer.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
virtual void checkType()
Check the type of the patch.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const scalar nut
label patchi
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
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
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
dimensionedScalar log(const dimensionedScalar &ds)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
error FatalError
static const char nl
Definition: Ostream.H:266
dictionary dict
volScalarField & p