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-2019 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 "turbulenceModel.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 turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
74  (
76  (
78  internalField().group()
79  )
80  );
81 
82  const scalarField& y = turbModel.y()[patchi];
83 
84  const tmp<volScalarField> tk = turbModel.k();
85  const volScalarField& k = tk();
86 
87  const tmp<scalarField> tnuw = turbModel.nu(patchi);
88  const scalarField& nuw = tnuw();
89 
90  const scalar Cmu25 = pow025(Cmu_);
91 
92  tmp<scalarField> tnutw(new scalarField(*this));
93  scalarField& nutw = tnutw.ref();
94 
95  forAll(nutw, facei)
96  {
97  const label celli = patch().faceCells()[facei];
98 
99  const scalar uStar = Cmu25*sqrt(k[celli]);
100  const scalar KsPlus = uStar*Ks_[facei]/nuw[facei];
101  const scalar E = this->E(KsPlus, Cs_[facei]);
102  const scalar yPlusMin = constant::mathematical::e/E;
103  const scalar yPlus = max(uStar*y[facei]/nuw[facei], yPlusMin);
104 
105  // To avoid oscillations limit the change in the wall viscosity
106  nutw[facei] =
107  max
108  (
109  min
110  (
111  nuw[facei]*max(yPlus*kappa_/log(E*yPlus) - 1, 0),
112  max(2*nutw[facei], nuw[facei])
113  ),
114  0.5*nutw[facei]
115  );
116 
117  if (debug)
118  {
119  Info<< "yPlus = " << yPlus
120  << ", KsPlus = " << KsPlus
121  << ", E = " << E
122  << ", yPlusMin " << yPlusMin
123  << ", yPlusLam " << yPlusLam(kappa_, E)
124  << ", nutw = " << nutw[facei]
125  << endl;
126  }
127  }
128 
129  return tnutw;
130 }
131 
132 
133 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
134 
136 (
137  const fvPatch& p,
139 )
140 :
142  Ks_(p.size(), 0.0),
143  Cs_(p.size(), 0.0)
144 {}
145 
146 
148 (
149  const fvPatch& p,
151  const dictionary& dict
152 )
153 :
155  Ks_("Ks", dict, p.size()),
156  Cs_("Cs", dict, p.size())
157 {}
158 
159 
161 (
163  const fvPatch& p,
165  const fvPatchFieldMapper& mapper
166 )
167 :
168  nutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
169  Ks_(mapper(ptf.Ks_)),
170  Cs_(mapper(ptf.Cs_))
171 {}
172 
173 
175 (
177 )
178 :
180  Ks_(rwfpsf.Ks_),
181  Cs_(rwfpsf.Cs_)
182 {}
183 
184 
186 (
189 )
190 :
192  Ks_(rwfpsf.Ks_),
193  Cs_(rwfpsf.Cs_)
194 {}
195 
196 
197 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
198 
200 (
201  const fvPatchFieldMapper& m
202 )
203 {
204  nutkWallFunctionFvPatchScalarField::autoMap(m);
205  m(Ks_, Ks_);
206  m(Cs_, Cs_);
207 }
208 
209 
211 (
212  const fvPatchScalarField& ptf,
213  const labelList& addr
214 )
215 {
216  nutkWallFunctionFvPatchScalarField::rmap(ptf, addr);
217 
219  refCast<const nutkRoughWallFunctionFvPatchScalarField>(ptf);
220 
221  Ks_.rmap(nrwfpsf.Ks_, addr);
222  Cs_.rmap(nrwfpsf.Cs_, addr);
223 }
224 
225 
227 {
229  writeLocalEntries(os);
230  writeEntry(os, "Cs", Cs_);
231  writeEntry(os, "Ks", Ks_);
232  writeEntry(os, "value", *this);
233 }
234 
235 
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
237 
239 (
242 );
243 
244 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
245 
246 } // End namespace Foam
247 
248 // ************************************************************************* //
const char *const group
Group name for atomic constants.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual tmp< scalarField > nut() const
Calculate the turbulence viscosity.
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
dimensionedScalar log(const dimensionedScalar &ds)
scalar yPlusLam() const
Return the Y+ at the edge of the laminar sublayer.
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
dimensionedScalar pow025(const dimensionedScalar &ds)
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions ...
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
label k
Boltzmann constant.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:280
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
scalar y
static const word propertiesName
Default name of the turbulence properties dictionary.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
const nearWallDist & y() const
Return the near wall distances.
virtual label size() const
Return size.
Definition: fvPatch.H:155
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
dimensionedScalar sin(const dimensionedScalar &ds)
nutkRoughWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:381
label patchi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
messageStream Info
A class for managing temporary objects.
Definition: PtrList.H:53
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
nutkWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Namespace for OpenFOAM.