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-2021 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 (
51  const fvPatch& p,
53  const dictionary& dict
54 )
55 :
57  Ceps2_(dict.lookupOrDefault<scalar>("Ceps2", 1.9))
58 {}
59 
60 
62 (
64  const fvPatch& p,
66  const fvPatchFieldMapper& mapper
67 )
68 :
69  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
70  Ceps2_(ptf.Ceps2_)
71 {}
72 
73 
75 (
78 )
79 :
81  Ceps2_(kwfpsf.Ceps2_)
82 {}
83 
84 
85 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
86 
88 {
89  if (updated())
90  {
91  return;
92  }
93 
94  const label patchi = patch().index();
95 
96  const momentumTransportModel& turbModel =
98  (
100  (
101  momentumTransportModel::typeName,
102  internalField().group()
103  )
104  );
105 
107  nutWallFunctionFvPatchScalarField::nutw(turbModel, patchi);
108 
109  const scalarField& y = turbModel.y()[patchi];
110 
111  const tmp<volScalarField> tk = turbModel.k();
112  const volScalarField& k = tk();
113 
114  const tmp<scalarField> tnuw = turbModel.nu(patchi);
115  const scalarField& nuw = tnuw();
116 
117  const scalar Cmu25 = pow025(nutw.Cmu());
118 
119  scalarField& kw = *this;
120 
121  // Set k wall values
122  forAll(kw, facei)
123  {
124  label celli = patch().faceCells()[facei];
125 
126  scalar uTau = Cmu25*sqrt(k[celli]);
127 
128  scalar yPlus = uTau*y[facei]/nuw[facei];
129 
130  if (yPlus > nutw.yPlusLam())
131  {
132  scalar Ck = -0.416;
133  scalar Bk = 8.366;
134  kw[facei] = Ck/nutw.kappa()*log(yPlus) + Bk;
135  }
136  else
137  {
138  scalar C = 11.0;
139  scalar Cf = (1.0/sqr(yPlus + C) + 2.0*yPlus/pow3(C) - 1.0/sqr(C));
140  kw[facei] = 2400.0/sqr(Ceps2_)*Cf;
141  }
142 
143  kw[facei] *= sqr(uTau);
144  }
145 
146  // Limit kw to avoid failure of the turbulence model due to division by kw
147  kw = max(kw, small);
148 
150 
151  // TODO: perform averaging for cells sharing more than one boundary face
152 }
153 
154 
156 {
157  writeEntry(os, "Ceps2", Ceps2_);
159 }
160 
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
165 (
168 );
169 
170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171 
172 } // End namespace Foam
173 
174 // ************************************************************************* //
const char *const group
Group name for atomic constants.
dictionary dict
Graphite solid properties.
Definition: C.H:48
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
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:175
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:369
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:63
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
label k
Boltzmann constant.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
const volScalarField::Boundary & y() const
Return the near wall distance.
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:99
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:351
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:216
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:147
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:357