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 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 
37 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38 
40 calcNut() const
41 {
42  const label patchi = patch().index();
43 
44  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
45  (
47  (
49  internalField().group()
50  )
51  );
52 
53  const scalarField& y = turbModel.y()[patchi];
54  const tmp<volScalarField> tk = turbModel.k();
55  const volScalarField& k = tk();
56  const tmp<scalarField> tnuw = turbModel.nu(patchi);
57  const scalarField& nuw = tnuw();
58 
59  const scalar Cmu25 = pow025(Cmu_);
60 
61  tmp<scalarField> tnutw(new scalarField(patch().size(), Zero));
62  scalarField& nutw = tnutw.ref();
63 
64  forAll(nutw, facei)
65  {
66  const label celli = patch().faceCells()[facei];
67 
68  const scalar yPlus = Cmu25*y[facei]*sqrt(k[celli])/nuw[facei];
69 
70  // Viscous sublayer contribution
71  const scalar nutVis = 0.0;
72 
73  // Inertial sublayer contribution
74  const scalar nutLog =
75  nuw[facei]*(yPlus*kappa_/log(max(E_*yPlus, 1 + 1e-4)) - 1.0);
76 
77  nutw[facei] = blend(nutVis, nutLog, yPlus);
78  }
79 
80  return tnutw;
81 }
82 
83 
84 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
85 
87 (
88  const fvPatch& p,
90 )
91 :
93 {}
94 
95 
97 (
99  const fvPatch& p,
101  const fvPatchFieldMapper& mapper
102 )
103 :
104  nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper)
105 {}
106 
107 
109 (
110  const fvPatch& p,
112  const dictionary& dict
113 )
114 :
116 {}
117 
118 
120 (
122 )
123 :
125 {}
126 
127 
129 (
132 )
133 :
135 {}
136 
137 
138 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
139 
141 yPlus() const
142 {
143  const label patchi = patch().index();
144 
145  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
146  (
148  (
150  internalField().group()
151  )
152  );
153 
154  const scalarField& y = turbModel.y()[patchi];
155 
156  const tmp<volScalarField> tk = turbModel.k();
157  const volScalarField& k = tk();
158  tmp<scalarField> kwc = k.boundaryField()[patchi].patchInternalField();
159  tmp<scalarField> nuw = turbModel.nu(patchi);
160 
161  return pow025(Cmu_)*y*sqrt(kwc)/nuw;
162 }
163 
164 
165 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166 
167 namespace Foam
168 {
170  (
173  );
174 }
175 
176 
177 // ************************************************************************* //
Foam::fvPatchField< scalar >
volFields.H
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:40
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:228
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:141
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:121
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
Foam::nutkWallFunctionFvPatchScalarField::nutkWallFunctionFvPatchScalarField
nutkWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: nutkWallFunctionFvPatchScalarField.C:87
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
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
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