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 | Copyright (C) 2011-2016 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 // * * * * * * * * * * * * Protected 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 
47  if (KsPlus < 90.0)
48  {
49  return pow
50  (
51  (KsPlus - 2.25)/87.75 + Cs*KsPlus,
52  sin(0.4258*(log(KsPlus) - 0.811))
53  );
54  }
55  else
56  {
57  return (1.0 + Cs*KsPlus);
58  }
59 }
60 
61 
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  label faceCelli = patch().faceCells()[facei];
88 
89  scalar uStar = Cmu25*sqrt(k[faceCelli]);
90  scalar yPlus = uStar*y[facei]/nuw[facei];
91  scalar KsPlus = uStar*Ks_[facei]/nuw[facei];
92 
93  scalar Edash = E_;
94  if (KsPlus > 2.25)
95  {
96  Edash /= fnRough(KsPlus, Cs_[facei]);
97  }
98 
99  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 
131 (
132  const fvPatch& p,
134 )
135 :
137  Ks_(p.size(), 0.0),
138  Cs_(p.size(), 0.0)
139 {}
140 
141 
143 (
145  const fvPatch& p,
147  const fvPatchFieldMapper& mapper
148 )
149 :
150  nutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
151  Ks_(ptf.Ks_, mapper),
152  Cs_(ptf.Cs_, mapper)
153 {}
154 
155 
157 (
158  const fvPatch& p,
160  const dictionary& dict
161 )
162 :
164  Ks_("Ks", dict, p.size()),
165  Cs_("Cs", dict, p.size())
166 {}
167 
168 
170 (
172 )
173 :
175  Ks_(rwfpsf.Ks_),
176  Cs_(rwfpsf.Cs_)
177 {}
178 
179 
181 (
184 )
185 :
187  Ks_(rwfpsf.Ks_),
188  Cs_(rwfpsf.Cs_)
189 {}
190 
191 
192 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
193 
195 (
196  const fvPatchFieldMapper& m
197 )
198 {
199  nutkWallFunctionFvPatchScalarField::autoMap(m);
200  Ks_.autoMap(m);
201  Cs_.autoMap(m);
202 }
203 
204 
206 (
207  const fvPatchScalarField& ptf,
208  const labelList& addr
209 )
210 {
211  nutkWallFunctionFvPatchScalarField::rmap(ptf, addr);
212 
214  refCast<const nutkRoughWallFunctionFvPatchScalarField>(ptf);
215 
216  Ks_.rmap(nrwfpsf.Ks_, addr);
217  Cs_.rmap(nrwfpsf.Cs_, addr);
218 }
219 
220 
222 {
224  writeLocalEntries(os);
225  Cs_.writeEntry("Cs", os);
226  Ks_.writeEntry("Ks", os);
227  writeEntry("value", os);
228 }
229 
230 
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232 
234 (
237 );
238 
239 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240 
241 } // End namespace Foam
242 
243 // ************************************************************************* //
virtual scalar fnRough(const scalar KsPlus, const scalar Cs) const
Compute the roughness function.
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:428
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)
const double e
Elementary charge.
Definition: doubleFloat.H:78
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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:253
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 ...
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:719
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
label k
Boltzmann constant.
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.
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
scalar y
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
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.
virtual label size() const
Return size.
Definition: fvPatch.H:161
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.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
void autoMap(const FieldMapper &map, const bool applyFlip=true)
Map from self.
Definition: Field.C:514
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:577
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:54
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
virtual tmp< scalarField > calcNut() const
Calculate the turbulence viscosity.
const nearWallDist & y() const
Return the near wall distances.
nutkWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Namespace for OpenFOAM.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:363