Loading...
Searching...
No Matches
nutUBlendedWallFunctionFvPatchScalarField.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) 2016-2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
29#include "turbulenceModel.H"
31#include "volFields.H"
33
34// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
35
38{
39 const label patchi = patch().index();
40
41 const auto& turbModel = db().lookupObject<turbulenceModel>
42 (
44 (
46 internalField().group()
47 )
48 );
49
50 const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
51 const scalarField magGradU(mag(Uw.snGrad()));
52
53 const tmp<scalarField> tnuw = turbModel.nu(patchi);
54 const scalarField& nuw = tnuw();
55
56 return max
57 (
58 scalar(0),
59 sqr(calcUTau(magGradU))/(magGradU + ROOTVSMALL) - nuw
60 );
61}
62
63
66(
67 const scalarField& magGradU
68) const
69{
70 const label patchi = patch().index();
71
72 const auto& turbModel = db().lookupObject<turbulenceModel>
73 (
75 (
77 internalField().group()
78 )
79 );
80
81 const scalarField& y = turbModel.y()[patchi];
82
83 const tmp<scalarField> tnuw = turbModel.nu(patchi);
84 const scalarField& nuw = tnuw();
85
86 const vectorField n(patch().nf());
87 const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
88 vectorField Up(Uw.patchInternalField() - Uw);
89 Up -= n*(n & Up);
90 const scalarField magUp(mag(Up));
91
92 auto tuTaup = tmp<scalarField>::New(patch().size(), Zero);
93 auto& uTaup = tuTaup.ref();
94
95 const scalarField& nutw = *this;
96
97 const scalar kappa = wallCoeffs_.kappa();
98 const scalar E = wallCoeffs_.E();
99
100 forAll(uTaup, facei)
101 {
102 scalar ut = sqrt((nutw[facei] + nuw[facei])*magGradU[facei]);
103 if (mag(ut) > ROOTVSMALL)
104 {
105 scalar error = GREAT;
106 label iter = 0;
107 while (iter++ < 10 && error > 0.001)
108 {
109 const scalar yPlus = y[facei]*ut/nuw[facei];
110 const scalar uTauVis = magUp[facei]/yPlus;
111 const scalar uTauLog =
112 kappa*magUp[facei]/log(max(E*yPlus, 1 + 1e-4));
113
114 const scalar utNew =
115 pow(pow(uTauVis, n_) + pow(uTauLog, n_), 1.0/n_);
116 error = mag(ut - utNew)/(ut + ROOTVSMALL);
117 ut = 0.5*(ut + utNew);
118 }
119 }
120 uTaup[facei] = ut;
121 }
122
123 return tuTaup;
124}
125
126
128(
129 Ostream& os
130) const
132 os.writeEntryIfDifferent<scalar>("n", 4.0, n_);
133}
134
135
136// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
137
140(
141 const fvPatch& p,
156 const fvPatchFieldMapper& mapper
169 const dictionary& dict
179(
193)
194:
196 n_(wfpsf.n_)
197{}
198
199
200// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
201
204{
205 const label patchi = patch().index();
206
207 const auto& turbModel = db().lookupObject<turbulenceModel>
208 (
210 (
212 internalField().group()
213 )
214 );
215
216 const scalarField& y = turbModel.y()[patchi];
217
218 const tmp<scalarField> tnuw = turbModel.nu(patchi);
219 const scalarField& nuw = tnuw();
220
221 const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
222 const scalarField magGradU(mag(Uw.snGrad()));
223
224 return y*calcUTau(magGradU)/nuw;
225}
226
227
229(
230 Ostream& os
231) const
232{
234 writeLocalEntries(os);
237
238
239// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240
241namespace Foam
242{
244 (
246 nutUBlendedWallFunctionFvPatchScalarField
247 );
248}
249
250
251// ************************************************************************* //
scalar y
label n
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:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Definition error.H:74
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 function for 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.
tmp< scalarField > calcUTau(const scalarField &magGradU) const
Calculate the friction velocity.
nutUBlendedWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
The class nutWallFunction is an abstract base class that hosts calculation methods and common functi...
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
wallFunctionCoefficients wallCoeffs_
Wall-function coefficients.
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.
U
Definition pEqn.H:72
volScalarField & p
scalar yPlus
scalar magUp
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
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.
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
fvPatchField< vector > fvPatchVectorField
fvPatchField< scalar > fvPatchScalarField
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299