atmNutkWallFunctionFvPatchScalarField.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-2018 OpenFOAM Foundation
9 Copyright (C) 2020 ENERCON GmbH
10 Copyright (C) 2020-2022 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
31#include "turbulenceModel.H"
32#include "fvPatchFieldMapper.H"
33#include "volFields.H"
34#include "bound.H"
36
37// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38
39namespace Foam
40{
41
42// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
43
45{
46 const label patchi = patch().index();
47
48 const auto& turbModel = db().lookupObject<turbulenceModel>
49 (
51 (
53 internalField().group()
54 )
55 );
56
57 const scalarField& y = turbModel.y()[patchi];
58
59 const tmp<volScalarField> tk = turbModel.k();
60 const volScalarField& k = tk();
61
62 const tmp<scalarField> tnuw = turbModel.nu(patchi);
63 const scalarField& nuw = tnuw();
64
65 auto tnutw = tmp<scalarField>::New(*this);
66 auto& nutw = tnutw.ref();
67
68 const scalar Cmu25 = pow025(wallCoeffs_.Cmu());
69 const scalar kappa = wallCoeffs_.kappa();
70
71 const scalar t = db().time().timeOutputValue();
72 const scalarField z0(z0_->value(t));
73
74 #ifdef FULLDEBUG
75 for (const scalar z : z0)
76 {
77 if (z < VSMALL)
78 {
80 << "z0 field can only contain positive values. "
81 << "Please check input field z0."
82 << exit(FatalError);
83 }
84 }
85 #endif
86
87 const labelList& faceCells = patch().faceCells();
88
89 // (HW:Eq. 5)
90 forAll(nutw, facei)
91 {
92 const label celli = faceCells[facei];
93
94 const scalar uStar = Cmu25*sqrt(k[celli]);
95 const scalar yPlus = uStar*y[facei]/nuw[facei];
96 const scalar Edash = (y[facei] + z0[facei])/z0[facei];
97
98 nutw[facei] = nuw[facei]*(yPlus*kappa/log(max(Edash, 1 + 1e-4)) - 1);
99 }
100
101 if (boundNut_)
102 {
103 nutw = max(nutw, scalar(0));
104 }
105
106 return tnutw;
107}
108
109
111(
112 Ostream& os
113) const
114{
115 os.writeEntryIfDifferent<bool>("boundNut", false, boundNut_);
116
117 if (z0_)
118 {
119 z0_->writeData(os);
120 }
121}
122
123
124// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
125
127(
128 const fvPatch& p,
130)
131:
133 boundNut_(true),
134 z0_(nullptr)
135{}
136
137
139(
141 const fvPatch& p,
143 const fvPatchFieldMapper& mapper
144)
145:
146 nutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
147 boundNut_(ptf.boundNut_),
148 z0_(ptf.z0_.clone(p.patch()))
149{}
150
151
153(
154 const fvPatch& p,
156 const dictionary& dict
157)
158:
160 boundNut_(dict.getOrDefault<bool>("boundNut", false)),
161 z0_(PatchFunction1<scalar>::New(p.patch(), "z0", dict))
162{}
163
164
166(
168)
169:
171 boundNut_(rwfpsf.boundNut_),
172 z0_(rwfpsf.z0_.clone(this->patch().patch()))
173{}
174
175
177(
180)
181:
183 boundNut_(rwfpsf.boundNut_),
184 z0_(rwfpsf.z0_.clone(this->patch().patch()))
185{}
186
187
188// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
189
191(
192 const fvPatchFieldMapper& m
193)
194{
196
197 if (z0_)
198 {
199 z0_->autoMap(m);
200 }
201}
202
203
205(
206 const fvPatchScalarField& ptf,
207 const labelList& addr
208)
209{
211
212 const auto& nrwfpsf =
213 refCast<const atmNutkWallFunctionFvPatchScalarField>(ptf);
214
215 if (z0_)
216 {
217 z0_->rmap(nrwfpsf.z0_(), addr);
218 }
219}
220
221
223{
226 writeEntry("value", os);
227}
228
229
230// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
231
233(
236);
237
238// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239
240} // End namespace Foam
241
242// ************************************************************************* //
scalar y
label k
Macros for easy insertion into run-time selection tables.
Bound the given scalar field if it has gone unbounded.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:251
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
void rmap(const atmBoundaryLayer &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
This boundary condition provides a wall constraint on the turbulent viscosity (i.e....
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
void writeLocalEntries(Ostream &) const
Write local wall function variables.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:59
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
wallFunctionCoefficients wallCoeffs_
Wall-function coefficients.
This boundary condition provides a wall function for the turbulent viscosity (i.e....
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
A class for managing temporary objects.
Definition: tmp.H:65
Abstract base class for turbulence models (RAS, LES and laminar).
static const word propertiesName
Default name of the turbulence properties dictionary.
scalar kappa() const noexcept
Return the object: kappa.
scalar Cmu() const noexcept
Return the object: Cmu.
volScalarField & p
bool
Definition: EEqn.H:20
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar sqrt(const dimensionedScalar &ds)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dimensionedScalar pow025(const dimensionedScalar &ds)
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333