Loading...
Searching...
No Matches
nutUTabulatedWallFunctionFvPatchScalarField.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 OpenFOAM Foundation
9 Copyright (C) 2019-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"
32#include "volFields.H"
34
35// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
36
39{
40 const label patchi = patch().index();
41
42 const auto& turbModel = db().lookupObject<turbulenceModel>
43 (
45 (
47 internalField().group()
48 )
49 );
50
51 const scalarField& y = turbModel.y()[patchi];
52
53 const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
54 const scalarField magUp(mag(Uw.patchInternalField() - Uw));
55 const scalarField magGradU(mag(Uw.snGrad()));
56
57 const tmp<scalarField> tnuw = turbModel.nu(patchi);
58 const scalarField& nuw = tnuw();
59
60 return
61 max
62 (
63 scalar(0),
64 sqr(magUp/(calcUPlus(magUp*y/nuw) + ROOTVSMALL))
65 /(magGradU + ROOTVSMALL)
66 - nuw
67 );
68}
69
70
73(
74 const scalarField& Rey
75) const
76{
77 auto tuPlus = tmp<scalarField>::New(patch().size(), Zero);
78 auto& uPlus = tuPlus.ref();
79
80 forAll(uPlus, facei)
81 {
82 uPlus[facei] = uPlusTable_.interpolateLog10(Rey[facei]);
83 }
84
85 return tuPlus;
86}
87
88
90(
91 Ostream& os
92) const
94 os.writeEntry("uPlusTable", uPlusTableName_);
95}
96
97
98// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
99
102(
103 const fvPatch& p,
105)
106:
108 uPlusTableName_("undefined-uPlusTableName"),
109 uPlusTable_
110 (
112 (
113 uPlusTableName_,
114 patch().boundaryMesh().mesh().time().constant(),
115 patch().boundaryMesh().mesh(),
116 IOobject::NO_READ,
117 IOobject::NO_WRITE,
131 const fvPatchFieldMapper& mapper
132)
134 nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
137{}
138
139
142(
143 const fvPatch& p,
145 const dictionary& dict
146)
147:
149 uPlusTableName_(dict.get<word>("uPlusTable")),
150 uPlusTable_
151 (
153 (
154 uPlusTableName_,
155 patch().boundaryMesh().mesh().time().constant(),
156 patch().boundaryMesh().mesh(),
157 IOobject::MUST_READ_IF_MODIFIED,
158 IOobject::NO_WRITE,
170)
183)
184:
186 uPlusTableName_(wfpsf.uPlusTableName_),
188{}
189
190
191// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
192
195{
196 const label patchi = patch().index();
197
198 const auto& turbModel = db().lookupObject<turbulenceModel>
199 (
201 (
203 internalField().group()
204 )
205 );
206
207 const scalarField& y = turbModel.y()[patchi];
208
209 const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
210 const scalarField magUp(mag(Uw.patchInternalField() - Uw));
211
212 const tmp<scalarField> tnuw = turbModel.nu(patchi);
213 const scalarField& nuw = tnuw();
215 const scalarField Rey(magUp*y/nuw);
216
217 return Rey/(calcUPlus(Rey) + ROOTVSMALL);
218}
219
220
222(
223 Ostream& os
224) const
225{
227 writeLocalEntries(os);
230
231
232// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233
234namespace Foam
235{
237 (
239 nutUTabulatedWallFunctionFvPatchScalarField
240 );
241}
242
243
244// ************************************************************************* //
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...
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
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:59
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A FieldMapper for finite-volume patch fields.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
This boundary condition provides a wall constraint on the turbulent viscosity (i.e....
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
void writeLocalEntries(Ostream &) const
Write local wall function variables.
nutUTabulatedWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
tmp< scalarField > calcUPlus(const scalarField &Rey) const
Calculate wall u+ from table.
virtual tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
The class nutWallFunction is an abstract base class that hosts calculation methods and common functi...
nutWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
Abstract base class for turbulence models (RAS, LES and laminar).
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
scalar uPlus
scalar magUp
scalar Rey
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
Different types of constants.
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:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
fvPatchField< vector > fvPatchVectorField
fvPatchField< scalar > fvPatchScalarField
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299