kLowReWallFunctionFvPatchScalarField.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) 2012-2016, 2019 OpenFOAM Foundation
9  Copyright (C) 2019-2022 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 
31 #include "turbulenceModel.H"
33 
34 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
35 
37 (
38  Ostream& os
39 ) const
40 {
41  os.writeEntryIfDifferent<scalar>("Ceps2", 1.9, Ceps2_);
42  os.writeEntryIfDifferent<scalar>("Ck", -0.416, Ck_);
43  os.writeEntryIfDifferent<scalar>("Bk", 8.366, Bk_);
44  os.writeEntryIfDifferent<scalar>("C", 11.0, C_);
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
52 (
53  const fvPatch& p,
55 )
56 :
57  fixedValueFvPatchField<scalar>(p, iF),
58  Ceps2_(1.9),
59  Ck_(-0.416),
60  Bk_(8.366),
61  C_(11.0),
62  wallCoeffs_()
63 {}
64 
65 
67 (
69  const fvPatch& p,
71  const fvPatchFieldMapper& mapper
72 )
73 :
74  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
75  Ceps2_(ptf.Ceps2_),
76  Ck_(ptf.Ck_),
77  Bk_(ptf.Bk_),
78  C_(ptf.C_),
79  wallCoeffs_(ptf.wallCoeffs_)
80 {}
81 
82 
84 (
85  const fvPatch& p,
87  const dictionary& dict
88 )
89 :
90  fixedValueFvPatchField<scalar>(p, iF, dict),
91  Ceps2_
92  (
93  dict.getCheckOrDefault<scalar>
94  (
95  "Ceps2",
96  1.9,
97  scalarMinMax::ge(SMALL)
98  )
99  ),
100  Ck_(dict.getOrDefault<scalar>("Ck", -0.416)),
101  Bk_(dict.getOrDefault<scalar>("Bk", 8.366)),
102  C_(dict.getOrDefault<scalar>("C", 11.0)),
103  wallCoeffs_(dict)
104 {}
105 
106 
108 (
110 )
111 :
112  fixedValueFvPatchField<scalar>(kwfpsf),
113  Ceps2_(kwfpsf.Ceps2_),
114  Ck_(kwfpsf.Ck_),
115  Bk_(kwfpsf.Bk_),
116  C_(kwfpsf.C_),
117  wallCoeffs_(kwfpsf.wallCoeffs_)
118 {}
119 
120 
122 (
125 )
126 :
127  fixedValueFvPatchField<scalar>(kwfpsf, iF),
128  Ceps2_(kwfpsf.Ceps2_),
129  Ck_(kwfpsf.Ck_),
130  Bk_(kwfpsf.Bk_),
131  C_(kwfpsf.C_),
132  wallCoeffs_(kwfpsf.wallCoeffs_)
133 {}
134 
135 
136 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
137 
139 {
140  if (updated())
141  {
142  return;
143  }
144 
145  const label patchi = patch().index();
146 
147  const auto& turbModel = db().lookupObject<turbulenceModel>
148  (
150  (
152  internalField().group()
153  )
154  );
155 
156  const scalarField& y = turbModel.y()[patchi];
157 
158  const tmp<scalarField> tnuw = turbModel.nu(patchi);
159  const scalarField& nuw = tnuw();
160 
161  const tmp<volScalarField> tk = turbModel.k();
162  const volScalarField& k = tk();
163 
164  const scalar Cmu25 = pow025(wallCoeffs_.Cmu());
165  const scalar kappa = wallCoeffs_.kappa();
166  const scalar yPlusLam = wallCoeffs_.yPlusLam();
167 
168  scalarField& kw = *this;
169 
170  // Set k wall values
171  forAll(kw, facei)
172  {
173  const label celli = patch().faceCells()[facei];
174  const scalar uTau = Cmu25*sqrt(k[celli]);
175  const scalar yPlus = uTau*y[facei]/nuw[facei];
176 
177  if (yPlus > yPlusLam)
178  {
179  kw[facei] = Ck_/kappa*log(yPlus) + Bk_;
180  }
181  else
182  {
183  const scalar Cf =
184  1.0/sqr(yPlus + C_) + 2.0*yPlus/pow3(C_) - 1.0/sqr(C_);
185  kw[facei] = 2400.0/sqr(Ceps2_)*Cf;
186  }
187 
188  kw[facei] *= sqr(uTau);
189  }
190 
191  // Limit kw to avoid failure of the turbulence model due to division by kw
192  kw = max(kw, SMALL);
193 
195 
196  // TODO: perform averaging for cells sharing more than one boundary face
197 }
198 
199 
201 (
202  Ostream& os
203 ) const
204 {
206  writeLocalEntries(os);
208 }
209 
210 
211 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
212 
213 namespace Foam
214 {
216  (
218  kLowReWallFunctionFvPatchScalarField
219  );
220 }
221 
222 
223 // ************************************************************************* //
wallFunctionCoefficients wallCoeffs_
Wall-function coefficients.
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
dictionary dict
dimensionedScalar log(const dimensionedScalar &ds)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
This boundary condition provides a wall function for the turbulent kinetic energy (i...
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
dimensionedScalar pow025(const dimensionedScalar &ds)
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
label k
Boltzmann constant.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
constexpr const char *const group
Group name for atomic constants.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
scalar uTau
scalar y
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
fvPatchField< scalar > fvPatchScalarField
static const word propertiesName
Default name of the turbulence properties dictionary.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
void writeEntries(Ostream &) const
Write wall-function coefficients as dictionary entries.
OBJstream os(runTime.globalPath()/outputName)
dimensionedScalar pow3(const dimensionedScalar &ds)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
scalar yPlus
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:309
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
volScalarField & p
void writeLocalEntries(Ostream &) const
Write local wall function variables.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
kLowReWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Namespace for OpenFOAM.