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 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 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
38 (
39  const fvPatch& p,
41 )
42 :
44  Ceps2_(1.9),
45  Ck_(-0.416),
46  Bk_(8.366),
47  C_(11.0)
48 {}
49 
50 
52 (
54  const fvPatch& p,
56  const fvPatchFieldMapper& mapper
57 )
58 :
59  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
60  Ceps2_(ptf.Ceps2_),
61  Ck_(ptf.Ck_),
62  Bk_(ptf.Bk_),
63  C_(ptf.C_)
64 {}
65 
66 
68 (
69  const fvPatch& p,
71  const dictionary& dict
72 )
73 :
75  Ceps2_
76  (
77  dict.getCheckOrDefault<scalar>
78  (
79  "Ceps2",
80  1.9,
81  scalarMinMax::ge(SMALL)
82  )
83  ),
84  Ck_(dict.getOrDefault<scalar>("Ck", -0.416)),
85  Bk_(dict.getOrDefault<scalar>("Bk", 8.366)),
86  C_(dict.getOrDefault<scalar>("C", 11.0))
87 {}
88 
89 
91 (
93 )
94 :
96  Ceps2_(kwfpsf.Ceps2_),
97  Ck_(kwfpsf.Ck_),
98  Bk_(kwfpsf.Bk_),
99  C_(kwfpsf.C_)
100 {}
101 
102 
104 (
107 )
108 :
109  fixedValueFvPatchField<scalar>(kwfpsf, iF),
110  Ceps2_(kwfpsf.Ceps2_),
111  Ck_(kwfpsf.Ck_),
112  Bk_(kwfpsf.Bk_),
113  C_(kwfpsf.C_)
114 {}
115 
116 
117 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118 
120 {
121  if (updated())
122  {
123  return;
124  }
125 
126  const label patchi = patch().index();
127 
128  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
129  (
131  (
133  internalField().group()
134  )
135  );
136 
138  nutWallFunctionFvPatchScalarField::nutw(turbModel, patchi);
139 
140  const scalarField& y = turbModel.y()[patchi];
141 
142  const tmp<scalarField> tnuw = turbModel.nu(patchi);
143  const scalarField& nuw = tnuw();
144 
145  const tmp<volScalarField> tk = turbModel.k();
146  const volScalarField& k = tk();
147 
148  const scalar Cmu25 = pow025(nutw.Cmu());
149 
150  scalarField& kw = *this;
151 
152  // Set k wall values
153  forAll(kw, facei)
154  {
155  const label celli = patch().faceCells()[facei];
156  const scalar uTau = Cmu25*sqrt(k[celli]);
157  const scalar yPlus = uTau*y[facei]/nuw[facei];
158 
159  if (yPlus > nutw.yPlusLam())
160  {
161  kw[facei] = Ck_/nutw.kappa()*log(yPlus) + Bk_;
162  }
163  else
164  {
165  const scalar Cf =
166  1.0/sqr(yPlus + C_) + 2.0*yPlus/pow3(C_) - 1.0/sqr(C_);
167  kw[facei] = 2400.0/sqr(Ceps2_)*Cf;
168  }
169 
170  kw[facei] *= sqr(uTau);
171  }
172 
173  // Limit kw to avoid failure of the turbulence model due to division by kw
174  kw = max(kw, SMALL);
175 
177 
178  // TODO: perform averaging for cells sharing more than one boundary face
179 }
180 
181 
183 (
184  Ostream& os
185 ) const
186 {
187  os.writeEntry("Ceps2", Ceps2_);
188  os.writeEntry("Ck", Ck_);
189  os.writeEntry("Bk", Bk_);
190  os.writeEntry("C", C_);
192 }
193 
194 
195 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196 
197 namespace Foam
198 {
200  (
203  );
204 }
205 
206 
207 // ************************************************************************* //
Foam::fvPatchField< scalar >
Foam::fvPatchField< scalar >::internalField
const DimensionedField< scalar, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:363
Foam::turbulenceModel::nu
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
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::kLowReWallFunctionFvPatchScalarField::Bk_
scalar Bk_
Bk coefficient.
Definition: kLowReWallFunctionFvPatchScalarField.H:180
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
nutWallFunctionFvPatchScalarField.H
Foam::kLowReWallFunctionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: kLowReWallFunctionFvPatchScalarField.C:119
Foam::kLowReWallFunctionFvPatchScalarField::kLowReWallFunctionFvPatchScalarField
kLowReWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: kLowReWallFunctionFvPatchScalarField.C:38
Foam::turbulenceModel::y
const nearWallDist & y() const
Return the near wall distances.
Definition: turbulenceModel.H:159
Foam::fvPatchField< scalar >::updated
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:387
Foam::fixedValueFvPatchField< scalar >
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::fixedValueFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fixedValueFvPatchField.C:156
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:133
Foam::Field< scalar >
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::kLowReWallFunctionFvPatchScalarField::C_
scalar C_
C coefficient.
Definition: kLowReWallFunctionFvPatchScalarField.H:183
Foam::kLowReWallFunctionFvPatchScalarField::Ck_
scalar Ck_
Ck coefficient.
Definition: kLowReWallFunctionFvPatchScalarField.H:177
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::kLowReWallFunctionFvPatchScalarField::Ceps2_
scalar Ceps2_
Ceps2 coefficient.
Definition: kLowReWallFunctionFvPatchScalarField.H:174
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
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
uTau
scalar uTau
Definition: evaluateNearWall.H:14
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::fvPatch::index
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:197
Foam::kLowReWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: kLowReWallFunctionFvPatchScalarField.C:183
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::turbulenceModel::k
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::fvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:321
Foam::kLowReWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the turbulent kinetic energy,...
Definition: kLowReWallFunctionFvPatchScalarField.H:165
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::fvPatchField< scalar >::db
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:206
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
kLowReWallFunctionFvPatchScalarField.H
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
Foam::fvPatch::faceCells
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:113
yPlus
scalar yPlus
Definition: evaluateNearWall.H:16
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::fvPatchField< scalar >::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:357
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 >
Foam::nutWallFunctionFvPatchScalarField
The class nutWallFunction is a base class that parents the derived boundary conditions which provide ...
Definition: nutWallFunctionFvPatchScalarField.H:251
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