nutWallFunctionFvPatchScalarField.C
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, 2019 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "fvPatchFieldMapper.H"
31 #include "volFields.H"
32 #include "wallFvPatch.H"
33 #include "turbulenceModel.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(nutWallFunctionFvPatchScalarField, 0);
41 }
42 
43 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
44 
46 {
47  if (!isA<wallFvPatch>(patch()))
48  {
50  << "Invalid wall function specification" << nl
51  << " Patch type for patch " << patch().name()
52  << " must be wall" << nl
53  << " Current patch type is " << patch().type() << nl << endl
54  << abort(FatalError);
55  }
56 }
57 
58 
60 (
61  const turbulenceModel& turb
62 ) const
63 {
64  if (UName_ == word::null)
65  {
66  return turb.U();
67  }
68  else
69  {
70  return db().lookupObject<volVectorField>(UName_);
71  }
72 }
73 
74 
76 (
77  Ostream& os
78 ) const
79 {
80  if (UName_ != word::null)
81  {
82  os.writeEntry("U", UName_);
83  }
84 
85  os.writeEntry("Cmu", Cmu_);
86  os.writeEntry("kappa", kappa_);
87  os.writeEntry("E", E_);
88 }
89 
90 
91 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
92 
94 (
95  const fvPatch& p,
97 )
98 :
99  fixedValueFvPatchScalarField(p, iF),
100  UName_(word::null),
101  Cmu_(0.09),
102  kappa_(0.41),
103  E_(9.8),
104  yPlusLam_(yPlusLam(kappa_, E_))
105 {
106  checkType();
107 }
108 
109 
111 (
113  const fvPatch& p,
115  const fvPatchFieldMapper& mapper
116 )
117 :
118  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
119  UName_(ptf.UName_),
120  Cmu_(ptf.Cmu_),
121  kappa_(ptf.kappa_),
122  E_(ptf.E_),
123  yPlusLam_(ptf.yPlusLam_)
124 {
125  checkType();
126 }
127 
128 
130 (
131  const fvPatch& p,
133  const dictionary& dict
134 )
135 :
136  fixedValueFvPatchScalarField(p, iF, dict),
137  UName_(dict.getOrDefault<word>("U", word::null)),
138  Cmu_(dict.getOrDefault<scalar>("Cmu", 0.09)),
139  kappa_(dict.getOrDefault<scalar>("kappa", 0.41)),
140  E_(dict.getOrDefault<scalar>("E", 9.8)),
141  yPlusLam_(yPlusLam(kappa_, E_))
142 {
143  checkType();
144 }
145 
146 
148 (
150 )
151 :
152  fixedValueFvPatchScalarField(wfpsf),
153  UName_(wfpsf.UName_),
154  Cmu_(wfpsf.Cmu_),
155  kappa_(wfpsf.kappa_),
156  E_(wfpsf.E_),
157  yPlusLam_(wfpsf.yPlusLam_)
158 {
159  checkType();
160 }
161 
162 
164 (
167 )
168 :
169  fixedValueFvPatchScalarField(wfpsf, iF),
170  UName_(wfpsf.UName_),
171  Cmu_(wfpsf.Cmu_),
172  kappa_(wfpsf.kappa_),
173  E_(wfpsf.E_),
174  yPlusLam_(wfpsf.yPlusLam_)
175 {
176  checkType();
177 }
178 
179 
180 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
181 
184 (
185  const turbulenceModel& turbModel,
186  const label patchi
187 )
188 {
189  return
190  refCast<const nutWallFunctionFvPatchScalarField>
191  (
192  turbModel.nut()().boundaryField()[patchi],
193  patchi
194  );
195 }
196 
197 
199 (
200  const scalar kappa,
201  const scalar E
202 )
203 {
204  scalar ypl = 11.0;
205 
206  for (int i = 0; i < 10; ++i)
207  {
208  ypl = log(max(E*ypl, 1))/kappa;
209  }
210 
211  return ypl;
212 }
213 
214 
216 {
217  return yPlusLam_;
218 }
219 
220 
222 {
223  if (updated())
224  {
225  return;
226  }
227 
228  operator==(calcNut());
229 
230  fixedValueFvPatchScalarField::updateCoeffs();
231 }
232 
233 
235 (
236  Ostream& os
237 ) const
238 {
240  writeLocalEntries(os);
241  writeEntry("value", os);
242 }
243 
244 
245 // ************************************************************************* //
Foam::nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField
nutWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: nutWallFunctionFvPatchScalarField.C:94
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:364
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::nutWallFunctionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: nutWallFunctionFvPatchScalarField.C:221
Foam::turbulenceModel::nut
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
Foam::nutWallFunctionFvPatchScalarField::yPlusLam_
scalar yPlusLam_
Estimated y+ value at the edge of the viscous sublayer.
Definition: nutWallFunctionFvPatchScalarField.H:133
nutWallFunctionFvPatchScalarField.H
wallFvPatch.H
Foam::nutWallFunctionFvPatchScalarField::U
virtual const volVectorField & U(const turbulenceModel &turb) const
Definition: nutWallFunctionFvPatchScalarField.C:60
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
fvPatchFieldMapper.H
Foam::constant::electromagnetic::kappa
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Foam::nutWallFunctionFvPatchScalarField::kappa_
scalar kappa_
Von Karman constant.
Definition: nutWallFunctionFvPatchScalarField.H:127
Foam::label
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:62
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::nutWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: nutWallFunctionFvPatchScalarField.C:235
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:63
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::nutWallFunctionFvPatchScalarField::E_
scalar E_
E coefficient.
Definition: nutWallFunctionFvPatchScalarField.H:130
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::nutWallFunctionFvPatchScalarField::Cmu_
scalar Cmu_
Cmu coefficient.
Definition: nutWallFunctionFvPatchScalarField.H:124
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::nutWallFunctionFvPatchScalarField::UName_
word UName_
Name of velocity field.
Definition: nutWallFunctionFvPatchScalarField.H:121
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::nutWallFunctionFvPatchScalarField::nutw
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
Definition: nutWallFunctionFvPatchScalarField.C:184
Foam::word::null
static const word null
An empty word.
Definition: word.H:77
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:219
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:122
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::nutWallFunctionFvPatchScalarField
The class nutWallFunction is a base class that parents the derived boundary conditions which provide ...
Definition: nutWallFunctionFvPatchScalarField.H:110
turb
compressible::turbulenceModel & turb
Definition: setRegionFluidFields.H:10
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
turbulenceModel.H
Foam::nutWallFunctionFvPatchScalarField::writeLocalEntries
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
Definition: nutWallFunctionFvPatchScalarField.C:76
Foam::nutWallFunctionFvPatchScalarField::yPlusLam
scalar yPlusLam() const
Return the y+ at the edge of the viscous sublayer.
Definition: nutWallFunctionFvPatchScalarField.C:215
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::nutWallFunctionFvPatchScalarField::checkType
virtual void checkType()
Check the type of the patch.
Definition: nutWallFunctionFvPatchScalarField.C:45