kLowReWallFunctionFvPatchScalarField.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) 2012-2020 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 
28 #include "momentumTransportModel.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
39 (
40  const fvPatch& p,
42 )
43 :
45  Ceps2_(1.9)
46 {}
47 
48 
50 (
52  const fvPatch& p,
54  const fvPatchFieldMapper& mapper
55 )
56 :
57  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
58  Ceps2_(ptf.Ceps2_)
59 {}
60 
61 
63 (
64  const fvPatch& p,
66  const dictionary& dict
67 )
68 :
70  Ceps2_(dict.lookupOrDefault<scalar>("Ceps2", 1.9))
71 {}
72 
73 
75 (
77 )
78 :
80  Ceps2_(kwfpsf.Ceps2_)
81 {}
82 
83 
85 (
88 )
89 :
91  Ceps2_(kwfpsf.Ceps2_)
92 {}
93 
94 
95 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96 
98 {
99  if (updated())
100  {
101  return;
102  }
103 
104  const label patchi = patch().index();
105 
106  const momentumTransportModel& turbModel =
108  (
110  (
111  momentumTransportModel::typeName,
112  internalField().group()
113  )
114  );
115 
117  nutWallFunctionFvPatchScalarField::nutw(turbModel, patchi);
118 
119  const scalarField& y = turbModel.y()[patchi];
120 
121  const tmp<volScalarField> tk = turbModel.k();
122  const volScalarField& k = tk();
123 
124  const tmp<scalarField> tnuw = turbModel.nu(patchi);
125  const scalarField& nuw = tnuw();
126 
127  const scalar Cmu25 = pow025(nutw.Cmu());
128 
129  scalarField& kw = *this;
130 
131  // Set k wall values
132  forAll(kw, facei)
133  {
134  label celli = patch().faceCells()[facei];
135 
136  scalar uTau = Cmu25*sqrt(k[celli]);
137 
138  scalar yPlus = uTau*y[facei]/nuw[facei];
139 
140  if (yPlus > nutw.yPlusLam())
141  {
142  scalar Ck = -0.416;
143  scalar Bk = 8.366;
144  kw[facei] = Ck/nutw.kappa()*log(yPlus) + Bk;
145  }
146  else
147  {
148  scalar C = 11.0;
149  scalar Cf = (1.0/sqr(yPlus + C) + 2.0*yPlus/pow3(C) - 1.0/sqr(C));
150  kw[facei] = 2400.0/sqr(Ceps2_)*Cf;
151  }
152 
153  kw[facei] *= sqr(uTau);
154  }
155 
156  // Limit kw to avoid failure of the turbulence model due to division by kw
157  kw = max(kw, small);
158 
160 
161  // TODO: perform averaging for cells sharing more than one boundary face
162 }
163 
164 
166 {
167  writeEntry(os, "Ceps2", Ceps2_);
169 }
170 
171 
172 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
173 
175 (
178 );
179 
180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181 
182 } // End namespace Foam
183 
184 // ************************************************************************* //
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
const char *const group
Group name for atomic constants.
dictionary dict
Graphite solid properties.
Definition: C.H:48
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
virtual void write(Ostream &) const
Write.
dimensionedScalar log(const dimensionedScalar &ds)
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:173
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:356
dimensionedSymmTensor sqr(const dimensionedVector &dv)
This boundary condition provides a turbulence kinetic energy wall function condition for low- and hig...
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
dimensionedScalar pow025(const dimensionedScalar &ds)
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
const nearWallDist & y() const
Return the near wall distances.
label k
Boltzmann constant.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Macros for easy insertion into run-time selection tables.
scalar uTau
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
scalar y
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions...
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:93
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:325
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
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 pow3(const dimensionedScalar &ds)
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
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:229
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:167
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
kLowReWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Namespace for OpenFOAM.
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:332