nutkRoughWallFunctionFvPatchScalarField.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) 2011-2016, 2019 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "turbulenceModel.H"
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
34 
35 
36 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
37 
39 (
40  const scalar KsPlus,
41  const scalar Cs
42 ) const
43 {
44  // Return fn based on non-dimensional roughness height
45 
46  if (KsPlus < 90.0)
47  {
48  return pow
49  (
50  (KsPlus - 2.25)/87.75 + Cs*KsPlus,
51  sin(0.4258*(log(KsPlus) - 0.811))
52  );
53  }
54  else
55  {
56  return (1.0 + Cs*KsPlus);
57  }
58 }
59 
60 
62 calcNut() const
63 {
64  const label patchi = patch().index();
65 
66  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
67  (
69  (
71  internalField().group()
72  )
73  );
74  const scalarField& y = turbModel.y()[patchi];
75  const tmp<volScalarField> tk = turbModel.k();
76  const volScalarField& k = tk();
77  const tmp<scalarField> tnuw = turbModel.nu(patchi);
78  const scalarField& nuw = tnuw();
79 
80  const scalar Cmu25 = pow025(Cmu_);
81 
82  tmp<scalarField> tnutw(new scalarField(*this));
83  scalarField& nutw = tnutw.ref();
84 
85  forAll(nutw, facei)
86  {
87  const label celli = patch().faceCells()[facei];
88 
89  const scalar uStar = Cmu25*sqrt(k[celli]);
90  const scalar yPlus = uStar*y[facei]/nuw[facei];
91  const scalar KsPlus = uStar*Ks_[facei]/nuw[facei];
92 
93  scalar Edash = E_;
94  if (2.25 < KsPlus)
95  {
96  Edash /= fnRough(KsPlus, Cs_[facei]);
97  }
98 
99  const scalar limitingNutw = max(nutw[facei], nuw[facei]);
100 
101  // To avoid oscillations limit the change in the wall viscosity
102  // which is particularly important if it temporarily becomes zero
103  nutw[facei] =
104  max
105  (
106  min
107  (
108  nuw[facei]
109  *(yPlus*kappa_/log(max(Edash*yPlus, 1+1e-4)) - 1),
110  2*limitingNutw
111  ), 0.5*limitingNutw
112  );
113 
114  if (debug)
115  {
116  Info<< "yPlus = " << yPlus
117  << ", KsPlus = " << KsPlus
118  << ", Edash = " << Edash
119  << ", nutw = " << nutw[facei]
120  << endl;
121  }
122  }
123 
124  return tnutw;
125 }
126 
127 
128 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
129 
132 (
133  const fvPatch& p,
135 )
136 :
138  Ks_(p.size(), Zero),
139  Cs_(p.size(), Zero)
140 {}
141 
142 
145 (
147  const fvPatch& p,
149  const fvPatchFieldMapper& mapper
150 )
151 :
152  nutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
153  Ks_(ptf.Ks_, mapper),
154  Cs_(ptf.Cs_, mapper)
155 {}
156 
157 
160 (
161  const fvPatch& p,
163  const dictionary& dict
164 )
165 :
167  Ks_("Ks", dict, p.size()),
168  Cs_("Cs", dict, p.size())
169 {}
170 
171 
174 (
176 )
177 :
179  Ks_(rwfpsf.Ks_),
180  Cs_(rwfpsf.Cs_)
181 {}
182 
183 
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  Ks_.autoMap(m);
206  Cs_.autoMap(m);
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 (
228  Ostream& os
229 ) const
230 {
232  writeLocalEntries(os);
233  Cs_.writeEntry("Cs", os);
234  Ks_.writeEntry("Ks", os);
235  writeEntry("value", os);
236 }
237 
238 
239 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240 
241 namespace Foam
242 {
244  (
247  );
248 }
249 
250 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvPatchField< scalar >
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::turbulenceModel::nu
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
Foam::nutkRoughWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: nutkRoughWallFunctionFvPatchScalarField.C:227
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::constant::atomic::group
constexpr const char *const group
Group name for atomic constants.
Definition: atomicConstants.H:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Foam::nutkRoughWallFunctionFvPatchScalarField::fnRough
virtual scalar fnRough(const scalar KsPlus, const scalar Cs) const
Compute the roughness function.
Definition: nutkRoughWallFunctionFvPatchScalarField.C:39
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
fvPatchFieldMapper.H
Foam::turbulenceModel::y
const nearWallDist & y() const
Return the near wall distances.
Definition: turbulenceModel.H:159
Foam::nutkRoughWallFunctionFvPatchScalarField::rmap
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: nutkRoughWallFunctionFvPatchScalarField.C:211
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::nutWallFunctionFvPatchScalarField::kappa_
scalar kappa_
von Kármán constant
Definition: nutWallFunctionFvPatchScalarField.H:290
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:133
nutkRoughWallFunctionFvPatchScalarField.H
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::nutkRoughWallFunctionFvPatchScalarField::nutkRoughWallFunctionFvPatchScalarField
nutkRoughWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: nutkRoughWallFunctionFvPatchScalarField.C:132
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::nutkRoughWallFunctionFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: nutkRoughWallFunctionFvPatchScalarField.C:200
Foam::nutkWallFunctionFvPatchScalarField::yPlus
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
Definition: nutkWallFunctionFvPatchScalarField.C:140
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:63
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::nutWallFunctionFvPatchScalarField::E_
scalar E_
Wall roughness parameter.
Definition: nutWallFunctionFvPatchScalarField.H:293
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::nutkRoughWallFunctionFvPatchScalarField::Ks_
scalarField Ks_
Roughness height.
Definition: nutkRoughWallFunctionFvPatchScalarField.H:121
Foam::nutWallFunctionFvPatchScalarField::Cmu_
scalar Cmu_
Empirical model coefficient.
Definition: nutWallFunctionFvPatchScalarField.H:287
Foam::turbulenceModel::k
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::nutkRoughWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the turbulent viscosity, i.e....
Definition: nutkRoughWallFunctionFvPatchScalarField.H:112
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::nutWallFunctionFvPatchScalarField::nutw
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
Definition: nutWallFunctionFvPatchScalarField.C:235
Foam::List< label >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Cs
const scalarField & Cs
Definition: solveBulkSurfactant.H:10
Foam::nutkRoughWallFunctionFvPatchScalarField::calcNut
virtual tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
Definition: nutkRoughWallFunctionFvPatchScalarField.C:62
Foam::nutkRoughWallFunctionFvPatchScalarField::Cs_
scalarField Cs_
Roughness constant.
Definition: nutkRoughWallFunctionFvPatchScalarField.H:124
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::GeometricField< scalar, fvPatchField, volMesh >
turbulenceModel.H
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::nutkWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the turbulent viscosity, i.e....
Definition: nutkWallFunctionFvPatchScalarField.H:94