nutkRoughWallFunctionFvPatchScalarField.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 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
38 
40 (
41  const scalar KsPlus,
42  const scalar Cs
43 ) const
44 {
45  // Return fn based on non-dimensional roughness height
46  if (KsPlus < 2.25)
47  {
48  return E_;
49  }
50  else if (KsPlus < 90)
51  {
52  return
53  E_
54  /pow
55  (
56  (KsPlus - 2.25)/87.75 + Cs*KsPlus,
57  sin(0.4258*(log(KsPlus) - 0.811))
58  );
59  }
60  else
61  {
62  return E_/(1 + Cs*KsPlus);
63  }
64 }
65 
66 
67 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
68 
70 {
71  const label patchi = patch().index();
72 
73  const momentumTransportModel& turbModel =
74  db().lookupType<momentumTransportModel>(internalField().group());
75 
76  const scalarField& y = turbModel.y()[patchi];
77 
78  const tmp<volScalarField> tk = turbModel.k();
79  const volScalarField& k = tk();
80 
81  const tmp<scalarField> tnuw = turbModel.nu(patchi);
82  const scalarField& nuw = tnuw();
83 
84  const scalar Cmu25 = pow025(Cmu_);
85 
86  tmp<scalarField> tnutw(new scalarField(*this));
87  scalarField& nutw = tnutw.ref();
88 
89  forAll(nutw, facei)
90  {
91  const label celli = patch().faceCells()[facei];
92 
93  const scalar uStar = Cmu25*sqrt(k[celli]);
94  const scalar KsPlus = uStar*Ks_[facei]/nuw[facei];
95  const scalar E = this->E(KsPlus, Cs_[facei]);
96  const scalar yPlusMin = constant::mathematical::e/E;
97  const scalar yPlus = max(uStar*y[facei]/nuw[facei], yPlusMin);
98 
99  // To avoid oscillations limit the change in the wall viscosity
100  nutw[facei] =
101  max
102  (
103  min
104  (
105  nuw[facei]*max(yPlus*kappa_/log(E*yPlus) - 1, 0),
106  max(2*nutw[facei], nuw[facei])
107  ),
108  0.5*nutw[facei]
109  );
110 
111  if (debug)
112  {
113  Info<< "yPlus = " << yPlus
114  << ", KsPlus = " << KsPlus
115  << ", E = " << E
116  << ", yPlusMin " << yPlusMin
117  << ", yPlusLam " << yPlusLam(kappa_, E)
118  << ", nutw = " << nutw[facei]
119  << endl;
120  }
121  }
122 
123  return tnutw;
124 }
125 
126 
127 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
128 
130 (
131  const fvPatch& p,
133  const dictionary& dict
134 )
135 :
137  Ks_("Ks", dict, p.size()),
138  Cs_("Cs", dict, p.size())
139 {}
140 
141 
143 (
145  const fvPatch& p,
147  const fvPatchFieldMapper& mapper
148 )
149 :
150  nutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
151  Ks_(mapper(ptf.Ks_)),
152  Cs_(mapper(ptf.Cs_))
153 {}
154 
155 
157 (
160 )
161 :
163  Ks_(rwfpsf.Ks_),
164  Cs_(rwfpsf.Cs_)
165 {}
166 
167 
168 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
169 
171 (
172  const fvPatchScalarField& ptf,
173  const fvPatchFieldMapper& mapper
174 )
175 {
176  nutkWallFunctionFvPatchScalarField::map(ptf, mapper);
177 
179  refCast<const nutkRoughWallFunctionFvPatchScalarField>(ptf);
180 
181  mapper(Ks_, nrwfpsf.Ks_);
182  mapper(Cs_, nrwfpsf.Cs_);
183 }
184 
185 
187 (
188  const fvPatchScalarField& ptf
189 )
190 {
191  nutkWallFunctionFvPatchScalarField::reset(ptf);
192 
194  refCast<const nutkRoughWallFunctionFvPatchScalarField>(ptf);
195 
196  Ks_.reset(nrwfpsf.Ks_);
197  Cs_.reset(nrwfpsf.Cs_);
198 }
199 
200 
202 {
204  writeLocalEntries(os);
205  writeEntry(os, "Cs", Cs_);
206  writeEntry(os, "Ks", Ks_);
207  writeEntry(os, "value", *this);
208 }
209 
210 
211 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
212 
214 (
217 );
218 
219 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
220 
221 } // End namespace Foam
222 
223 // ************************************************************************* //
scalar y
label k
#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...
void reset(const Field< Type > &)
Reset the field values to the given field.
Definition: Field.C:432
Generic GeometricField class.
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
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
Abstract base class for turbulence models (RAS, LES and laminar).
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
const volScalarField::Boundary & y() const
Return the near wall distance.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
scalar yPlusLam() const
Return the Y+ at the edge of the laminar sublayer.
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions ...
virtual tmp< scalarField > nut() const
Calculate the turbulence viscosity.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
const scalarField & Cs() const
Return the roughness constant scale.
virtual void map(const fvPatchScalarField &, const fvPatchFieldMapper &)
Map the given fvPatchField onto this fvPatchField.
nutkRoughWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
Turbulent viscosity wall-function boundary condition for high Reynolds number flows based on near-wal...
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
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
label patchi
Namespace for OpenFOAM.
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
dimensionedScalar sin(const dimensionedScalar &ds)
messageStream Info
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar sqrt(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dimensionedScalar pow025(const dimensionedScalar &ds)
dictionary dict
volScalarField & p