nutkWallFunctionFvPatchScalarField.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) 2011-2016, 2019 OpenFOAM Foundation
9  Copyright (C) 2019-2021 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 
30 #include "turbulenceModel.H"
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
33 #include "wallFvPatch.H"
35 
36 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
37 
39 calcNut() const
40 {
41  const label patchi = patch().index();
42 
43  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
44  (
46  (
48  internalField().group()
49  )
50  );
51 
52  const scalarField& y = turbModel.y()[patchi];
53  const tmp<volScalarField> tk = turbModel.k();
54  const volScalarField& k = tk();
55  const tmp<scalarField> tnuw = turbModel.nu(patchi);
56  const scalarField& nuw = tnuw();
57 
58  const scalar Cmu25 = pow025(Cmu_);
59 
60  tmp<scalarField> tnutw(new scalarField(patch().size(), Zero));
61  scalarField& nutw = tnutw.ref();
62 
63  forAll(nutw, facei)
64  {
65  const label celli = patch().faceCells()[facei];
66 
67  const scalar yPlus = Cmu25*y[facei]*sqrt(k[celli])/nuw[facei];
68 
69  // Viscous sublayer contribution
70  const scalar nutVis = 0.0;
71 
72  // Inertial sublayer contribution
73  const scalar nutLog =
74  nuw[facei]*(yPlus*kappa_/log(max(E_*yPlus, 1 + 1e-4)) - 1.0);
75 
76  nutw[facei] = blend(nutVis, nutLog, yPlus);
77  }
78 
79  return tnutw;
80 }
81 
82 
83 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
84 
86 (
87  const fvPatch& p,
89 )
90 :
92 {}
93 
94 
96 (
98  const fvPatch& p,
100  const fvPatchFieldMapper& mapper
101 )
102 :
103  nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper)
104 {}
105 
106 
108 (
109  const fvPatch& p,
111  const dictionary& dict
112 )
113 :
115 {}
116 
117 
119 (
121 )
122 :
124 {}
125 
126 
128 (
131 )
132 :
134 {}
135 
136 
137 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
138 
140 yPlus() const
141 {
142  const label patchi = patch().index();
143 
144  const auto& turbModel = db().lookupObject<turbulenceModel>
145  (
147  (
149  internalField().group()
150  )
151  );
152 
153  const scalarField& y = turbModel.y()[patchi];
154 
155  tmp<volScalarField> tk = turbModel.k();
156  const volScalarField& k = tk();
157  tmp<scalarField> tkwc = k.boundaryField()[patchi].patchInternalField();
158  const scalarField& kwc = tkwc();
159 
160  tmp<scalarField> tnuw = turbModel.nu(patchi);
161  const scalarField& nuw = tnuw();
162 
163  tmp<scalarField> tnuEff = turbModel.nuEff(patchi);
164  const scalarField& nuEff = tnuEff();
165 
166  const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
167  const scalarField magGradUw(mag(Uw.snGrad()));
168 
169  const scalar Cmu25 = pow025(Cmu_);
170 
171  auto tyPlus = tmp<scalarField>::New(patch().size(), Zero);
172  auto& yPlus = tyPlus.ref();
173 
174  forAll(yPlus, facei)
175  {
176  // inertial sublayer
177  yPlus[facei] = Cmu25*y[facei]*sqrt(kwc[facei])/nuw[facei];
178 
179  if (yPlusLam_ > yPlus[facei])
180  {
181  // viscous sublayer
182  yPlus[facei] =
183  y[facei]*sqrt(nuEff[facei]*magGradUw[facei])/nuw[facei];
184  }
185  }
186 
187  return tyPlus;
188 }
189 
190 
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192 
193 namespace Foam
194 {
196  (
199  );
200 }
201 
202 
203 // ************************************************************************* //
Foam::fvPatchField< vector >
volFields.H
Foam::fvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:225
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
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::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
wallFvPatch.H
fvPatchFieldMapper.H
Foam::turbulenceModel::y
const nearWallDist & y() const
Return the near wall distances.
Definition: turbulenceModel.H:159
Foam::nutkWallFunctionFvPatchScalarField::calcNut
virtual tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
Definition: nutkWallFunctionFvPatchScalarField.C:39
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::nutWallFunctionFvPatchScalarField::kappa_
scalar kappa_
von Kármán constant
Definition: nutWallFunctionFvPatchScalarField.H:290
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:133
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::Field< scalar >
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::nutkWallFunctionFvPatchScalarField::yPlus
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
Definition: nutkWallFunctionFvPatchScalarField.C:140
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
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::nutWallFunctionFvPatchScalarField::E_
scalar E_
Wall roughness parameter.
Definition: nutWallFunctionFvPatchScalarField.H:293
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::nutWallFunctionFvPatchScalarField::Cmu_
scalar Cmu_
Empirical model coefficient.
Definition: nutWallFunctionFvPatchScalarField.H:287
Foam::turbulenceModel::k
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
U
U
Definition: pEqn.H:72
Foam::nutkWallFunctionFvPatchScalarField::nutkWallFunctionFvPatchScalarField
nutkWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: nutkWallFunctionFvPatchScalarField.C:86
Foam::nutWallFunctionFvPatchScalarField::blend
scalar blend(const scalar nutVis, const scalar nutLog, const scalar yPlus) const
Return the blended nut according to the chosen blending treatment.
Definition: nutWallFunctionFvPatchScalarField.C:267
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
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
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
yPlus
scalar yPlus
Definition: evaluateNearWall.H:16
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
nutkWallFunctionFvPatchScalarField.H
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
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
Foam::nutkWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the turbulent viscosity, i.e....
Definition: nutkWallFunctionFvPatchScalarField.H:94