nutUSpaldingWallFunctionFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "momentumTransportModel.H"
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38 
40 {
41  const label patchi = patch().index();
42 
43  const momentumTransportModel& turbModel =
44  db().lookupType<momentumTransportModel>(internalField().group());
45 
46  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
47  const scalarField magGradU(mag(Uw.snGrad()));
48  const tmp<scalarField> tnuw = turbModel.nu(patchi);
49  const scalarField& nuw = tnuw();
50 
51  return max
52  (
53  scalar(0),
54  sqr(calcUTau(magGradU))/(magGradU + rootVSmall) - nuw
55  );
56 }
57 
58 
60 (
61  const scalarField& magGradU
62 ) const
63 {
64  const label patchi = patch().index();
65 
66  const momentumTransportModel& turbModel =
67  db().lookupType<momentumTransportModel>(internalField().group());
68 
69  const scalarField& y = turbModel.y()[patchi];
70 
71  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
72  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
73 
74  const tmp<scalarField> tnuw = turbModel.nu(patchi);
75  const scalarField& nuw = tnuw();
76 
77  const scalarField& nutw = *this;
78 
79  tmp<scalarField> tuTau(new scalarField(patch().size(), 0.0));
80  scalarField& uTau = tuTau.ref();
81 
82  forAll(uTau, facei)
83  {
84  scalar ut = sqrt((nutw[facei] + nuw[facei])*magGradU[facei]);
85 
86  if (ut > rootVSmall)
87  {
88  int iter = 0;
89  scalar err = great;
90 
91  do
92  {
93  scalar kUu = min(kappa_*magUp[facei]/ut, 50);
94  scalar fkUu = exp(kUu) - 1 - kUu*(1 + 0.5*kUu);
95 
96  scalar f =
97  - ut*y[facei]/nuw[facei]
98  + magUp[facei]/ut
99  + 1/E_*(fkUu - 1.0/6.0*kUu*sqr(kUu));
100 
101  scalar df =
102  y[facei]/nuw[facei]
103  + magUp[facei]/sqr(ut)
104  + 1/E_*kUu*fkUu/ut;
105 
106  scalar uTauNew = ut + f/df;
107  err = mag((ut - uTauNew)/ut);
108  ut = uTauNew;
109 
110  } while (ut > rootVSmall && err > 0.01 && ++iter < 10);
111 
112  uTau[facei] = max(0.0, ut);
113  }
114  }
115 
116  return tuTau;
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
121 
124 (
125  const fvPatch& p,
127  const dictionary& dict
128 )
129 :
131 {}
132 
133 
136 (
138  const fvPatch& p,
140  const fvPatchFieldMapper& mapper
141 )
142 :
143  nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper)
144 {}
145 
146 
149 (
152 )
153 :
155 {}
156 
157 
158 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
159 
161 {
162  const label patchi = patch().index();
163 
164  const momentumTransportModel& turbModel =
165  db().lookupType<momentumTransportModel>(internalField().group());
166 
167  const scalarField& y = turbModel.y()[patchi];
168  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
169  const tmp<scalarField> tnuw = turbModel.nu(patchi);
170  const scalarField& nuw = tnuw();
171 
172  return y*calcUTau(mag(Uw.snGrad()))/nuw;
173 }
174 
175 
177 {
179  writeLocalEntries(os);
180  writeEntry(os, "value", *this);
181 }
182 
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
187 (
190 );
191 
192 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
193 
194 } // End namespace Foam
195 
196 // ************************************************************************* //
scalar y
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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...
const Boundary & boundaryField() const
Return const-reference to the boundary field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:172
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:164
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
Abstract base class for turbulence models (RAS, LES and laminar).
const volVectorField & U() const
Access function to velocity field.
const volScalarField::Boundary & y() const
Return the near wall distance.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions ...
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
virtual tmp< scalarField > nut() const
Calculate the turbulence viscosity.
virtual tmp< scalarField > calcUTau(const scalarField &magGradU) const
Calculate the friction velocity.
nutUSpaldingWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions,...
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
const scalar uTau
const scalar magUp
label patchi
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
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:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
labelList f(nPoints)
dictionary dict
volScalarField & p