atmNutUWallFunctionFvPatchScalarField.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) 2020 ENERCON GmbH
9 Copyright (C) 2020-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "turbulenceModel.H"
31#include "fvPatchFieldMapper.H"
32#include "volFields.H"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37namespace Foam
38{
39
40// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41
43{
44 const label patchi = patch().index();
45
46 const auto& turbModel = db().lookupObject<turbulenceModel>
47 (
49 (
51 internalField().group()
52 )
53 );
54
55 const scalarField& y = turbModel.y()[patchi];
56
57 const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
58 const vectorField Up(Uw.patchInternalField() - Uw);
59 const scalarField magUpn(mag(Up - (Up & patch().nf())*patch().nf()));
60
61 const tmp<scalarField> tnuw = turbModel.nu(patchi);
62 const scalarField& nuw = tnuw();
63
64 const scalar kappa = wallCoeffs_.kappa();
65
66 const scalar t = db().time().timeOutputValue();
67 const scalarField z0(z0_->value(t));
68
69 #ifdef FULLDEBUG
70 for (const scalar z : z0)
71 {
72 if (z < VSMALL)
73 {
75 << "z0 field can only contain positive values. "
76 << "Please check input field z0."
77 << exit(FatalError);
78 }
79 }
80 #endif
81
82 auto tnutw = tmp<scalarField>::New(*this);
83 auto& nutw = tnutw.ref();
84
85 forAll(nutw, facei)
86 {
87 const scalar Edash = (y[facei] + z0[facei])/(z0[facei] + 1e-4);
88 const scalar uStar = magUpn[facei]*kappa/log(max(Edash, 1.0 + 1e-4));
89
90 nutw[facei] = sqr(uStar)/max(magUpn[facei], 1e-6)*y[facei] - nuw[facei];
91 }
92
93 if (boundNut_)
94 {
95 nutw = max(nutw, scalar(0));
96 }
97
98 return tnutw;
99}
100
101
103(
104 Ostream& os
105) const
106{
107 os.writeEntryIfDifferent<bool>("boundNut", true, boundNut_);
108
109 if (z0_)
110 {
111 z0_->writeData(os);
112 }
113}
114
115
116// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
117
119(
120 const fvPatch& p,
122)
123:
125 boundNut_(true),
126 z0_(nullptr)
127{}
128
129
131(
133 const fvPatch& p,
135 const fvPatchFieldMapper& mapper
136)
137:
138 nutUWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
139 boundNut_(ptf.boundNut_),
140 z0_(ptf.z0_.clone(p.patch()))
141{}
142
143
145(
146 const fvPatch& p,
148 const dictionary& dict
149)
150:
152 boundNut_(dict.getOrDefault<bool>("boundNut", true)),
153 z0_(PatchFunction1<scalar>::New(p.patch(), "z0", dict))
154{}
155
156
158(
160)
161:
163 boundNut_(rwfpsf.boundNut_),
164 z0_(rwfpsf.z0_.clone(this->patch().patch()))
165{}
166
167
169(
172)
173:
175 boundNut_(rwfpsf.boundNut_),
176 z0_(rwfpsf.z0_.clone(this->patch().patch()))
177{}
178
179
180// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
181
183(
184 const fvPatchFieldMapper& m
185)
186{
188
189 if (z0_)
190 {
191 z0_->autoMap(m);
192 }
193}
194
195
197(
198 const fvPatchScalarField& ptf,
199 const labelList& addr
200)
201{
203
205 refCast<const atmNutUWallFunctionFvPatchScalarField>(ptf);
206
207 if (z0_)
208 {
209 z0_->rmap(nrwfpsf.z0_(), addr);
210 }
211}
212
213
215{
218 writeEntry("value", os);
219}
220
221
222// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223
225(
228);
229
230// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
231
232} // End namespace Foam
233
234// ************************************************************************* //
scalar y
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...
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
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:237
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
This boundary condition provides a wall function for the turbulent viscosity (i.e....
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
wallFunctionCoefficients wallCoeffs_
Wall-function coefficients.
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.
U
Definition: pEqn.H:72
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar log(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333