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-2022 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().lookupObject<momentumTransportModel>
75  (
77  (
78  momentumTransportModel::typeName,
79  internalField().group()
80  )
81  );
82 
83  const scalarField& y = turbModel.y()[patchi];
84 
85  const tmp<volScalarField> tk = turbModel.k();
86  const volScalarField& k = tk();
87 
88  const tmp<scalarField> tnuw = turbModel.nu(patchi);
89  const scalarField& nuw = tnuw();
90 
91  const scalar Cmu25 = pow025(Cmu_);
92 
93  tmp<scalarField> tnutw(new scalarField(*this));
94  scalarField& nutw = tnutw.ref();
95 
96  forAll(nutw, facei)
97  {
98  const label celli = patch().faceCells()[facei];
99 
100  const scalar uStar = Cmu25*sqrt(k[celli]);
101  const scalar KsPlus = uStar*Ks_[facei]/nuw[facei];
102  const scalar E = this->E(KsPlus, Cs_[facei]);
103  const scalar yPlusMin = constant::mathematical::e/E;
104  const scalar yPlus = max(uStar*y[facei]/nuw[facei], yPlusMin);
105 
106  // To avoid oscillations limit the change in the wall viscosity
107  nutw[facei] =
108  max
109  (
110  min
111  (
112  nuw[facei]*max(yPlus*kappa_/log(E*yPlus) - 1, 0),
113  max(2*nutw[facei], nuw[facei])
114  ),
115  0.5*nutw[facei]
116  );
117 
118  if (debug)
119  {
120  Info<< "yPlus = " << yPlus
121  << ", KsPlus = " << KsPlus
122  << ", E = " << E
123  << ", yPlusMin " << yPlusMin
124  << ", yPlusLam " << yPlusLam(kappa_, E)
125  << ", nutw = " << nutw[facei]
126  << endl;
127  }
128  }
129 
130  return tnutw;
131 }
132 
133 
134 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
135 
137 (
138  const fvPatch& p,
140 )
141 :
143  Ks_(p.size(), 0.0),
144  Cs_(p.size(), 0.0)
145 {}
146 
147 
149 (
150  const fvPatch& p,
152  const dictionary& dict
153 )
154 :
156  Ks_("Ks", dict, p.size()),
157  Cs_("Cs", dict, p.size())
158 {}
159 
160 
162 (
164  const fvPatch& p,
166  const fvPatchFieldMapper& mapper
167 )
168 :
169  nutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
170  Ks_(mapper(ptf.Ks_)),
171  Cs_(mapper(ptf.Cs_))
172 {}
173 
174 
176 (
179 )
180 :
182  Ks_(rwfpsf.Ks_),
183  Cs_(rwfpsf.Cs_)
184 {}
185 
186 
187 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
188 
190 (
191  const fvPatchFieldMapper& m
192 )
193 {
194  nutkWallFunctionFvPatchScalarField::autoMap(m);
195  m(Ks_, Ks_);
196  m(Cs_, Cs_);
197 }
198 
199 
201 (
202  const fvPatchScalarField& ptf,
203  const labelList& addr
204 )
205 {
206  nutkWallFunctionFvPatchScalarField::rmap(ptf, addr);
207 
209  refCast<const nutkRoughWallFunctionFvPatchScalarField>(ptf);
210 
211  Ks_.rmap(nrwfpsf.Ks_, addr);
212  Cs_.rmap(nrwfpsf.Cs_, addr);
213 }
214 
215 
217 (
218  const fvPatchScalarField& ptf
219 )
220 {
221  nutkWallFunctionFvPatchScalarField::reset(ptf);
222 
224  refCast<const nutkRoughWallFunctionFvPatchScalarField>(ptf);
225 
226  Ks_.reset(nrwfpsf.Ks_);
227  Cs_.reset(nrwfpsf.Cs_);
228 }
229 
230 
232 {
234  writeLocalEntries(os);
235  writeEntry(os, "Cs", Cs_);
236  writeEntry(os, "Ks", Ks_);
237  writeEntry(os, "value", *this);
238 }
239 
240 
241 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
242 
244 (
247 );
248 
249 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
250 
251 } // End namespace Foam
252 
253 // ************************************************************************* //
const char *const group
Group name for atomic constants.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual tmp< scalarField > nut() const
Calculate the turbulence viscosity.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
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:156
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
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:243
const volScalarField::Boundary & y() const
Return the near wall distance.
Macros for easy insertion into run-time selection tables.
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
scalar y
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
virtual label size() const
Return size.
Definition: fvPatch.H:157
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
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
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Abstract base class for turbulence models (RAS, LES and laminar).
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:362
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
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
void reset(const Field< Type > &)
Reset the field values to the given field.
Definition: Field.C:432
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
nutkWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Namespace for OpenFOAM.