nutUTabulatedWallFunctionFvPatchScalarField.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 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"
34 
35 
36 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
37 
40 {
41  const label patchi = patch().index();
42 
43  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
44  (
46  (
48  internalField().group()
49  )
50  );
51  const scalarField& y = turbModel.y()[patchi];
52  const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
53  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
54  const scalarField magGradU(mag(Uw.snGrad()));
55  const tmp<scalarField> tnuw = turbModel.nu(patchi);
56  const scalarField& nuw = tnuw();
57 
58  return
59  max
60  (
61  scalar(0),
62  sqr(magUp/(calcUPlus(magUp*y/nuw) + ROOTVSMALL))
63  /(magGradU + ROOTVSMALL)
64  - nuw
65  );
66 }
67 
68 
71 (
72  const scalarField& Rey
73 ) const
74 {
75  tmp<scalarField> tuPlus(new scalarField(patch().size(), Zero));
76  scalarField& uPlus = tuPlus.ref();
77 
78  forAll(uPlus, facei)
79  {
80  uPlus[facei] = uPlusTable_.interpolateLog10(Rey[facei]);
81  }
82 
83  return tuPlus;
84 }
85 
86 
87 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
88 
91 (
92  const fvPatch& p,
94 )
95 :
97  uPlusTableName_("undefined-uPlusTableName"),
98  uPlusTable_
99  (
100  IOobject
101  (
102  uPlusTableName_,
103  patch().boundaryMesh().mesh().time().constant(),
104  patch().boundaryMesh().mesh(),
107  false
108  ),
109  false
110  )
111 {}
112 
113 
116 (
118  const fvPatch& p,
120  const fvPatchFieldMapper& mapper
121 )
122 :
123  nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
124  uPlusTableName_(ptf.uPlusTableName_),
125  uPlusTable_(ptf.uPlusTable_)
126 {}
127 
128 
131 (
132  const fvPatch& p,
134  const dictionary& dict
135 )
136 :
138  uPlusTableName_(dict.get<word>("uPlusTable")),
139  uPlusTable_
140  (
141  IOobject
142  (
143  uPlusTableName_,
144  patch().boundaryMesh().mesh().time().constant(),
145  patch().boundaryMesh().mesh(),
148  false
149  ),
150  true
151  )
152 {}
153 
154 
157 (
159 )
160 :
162  uPlusTableName_(wfpsf.uPlusTableName_),
163  uPlusTable_(wfpsf.uPlusTable_)
164 {}
165 
166 
169 (
172 )
173 :
175  uPlusTableName_(wfpsf.uPlusTableName_),
176  uPlusTable_(wfpsf.uPlusTable_)
177 {}
178 
179 
180 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
181 
184 {
185  const label patchi = patch().index();
186 
187  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
188  (
190  (
192  internalField().group()
193  )
194  );
195  const scalarField& y = turbModel.y()[patchi];
196  const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
197  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
198  const tmp<scalarField> tnuw = turbModel.nu(patchi);
199  const scalarField& nuw = tnuw();
200  const scalarField Rey(magUp*y/nuw);
201 
202  return Rey/(calcUPlus(Rey) + ROOTVSMALL);
203 }
204 
205 
207 (
208  Ostream& os
209 ) const
210 {
212  os.writeEntry("uPlusTable", uPlusTableName_);
213  writeEntry("value", os);
214 }
215 
216 
217 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
218 
219 namespace Foam
220 {
222  (
225  );
226 }
227 
228 
229 // ************************************************************************* //
Foam::fvPatchField< vector >
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::nutUTabulatedWallFunctionFvPatchScalarField::calcUPlus
virtual tmp< scalarField > calcUPlus(const scalarField &Rey) const
Calculate wall u+ from table.
Definition: nutUTabulatedWallFunctionFvPatchScalarField.C:71
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::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::turbulenceModel::nu
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
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::nutUTabulatedWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the turbulent viscosity, i.e....
Definition: nutUTabulatedWallFunctionFvPatchScalarField.H:108
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::nutWallFunctionFvPatchScalarField::U
virtual const volVectorField & U(const turbulenceModel &turb) const
Definition: nutWallFunctionFvPatchScalarField.C:75
uPlus
scalar uPlus
Definition: evaluateNearWall.H:18
fvPatchFieldMapper.H
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
Foam::nutUTabulatedWallFunctionFvPatchScalarField::calcNut
virtual tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
Definition: nutUTabulatedWallFunctionFvPatchScalarField.C:39
Foam::turbulenceModel::y
const nearWallDist & y() const
Return the near wall distances.
Definition: turbulenceModel.H:159
magUp
scalar magUp
Definition: evaluateNearWall.H:10
Foam::nutUTabulatedWallFunctionFvPatchScalarField::yPlus
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
Definition: nutUTabulatedWallFunctionFvPatchScalarField.C:183
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::Field< scalar >
Foam::nutUTabulatedWallFunctionFvPatchScalarField::uPlusTable_
uniformInterpolationTable< scalar > uPlusTable_
u+ table
Definition: nutUTabulatedWallFunctionFvPatchScalarField.H:120
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::fvPatchField::patchInternalField
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:233
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
Foam::nutUTabulatedWallFunctionFvPatchScalarField::nutUTabulatedWallFunctionFvPatchScalarField
nutUTabulatedWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: nutUTabulatedWallFunctionFvPatchScalarField.C:91
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
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::nutUTabulatedWallFunctionFvPatchScalarField::write
virtual void write(Ostream &os) const
Write.
Definition: nutUTabulatedWallFunctionFvPatchScalarField.C:207
Foam::nutUTabulatedWallFunctionFvPatchScalarField::uPlusTableName_
word uPlusTableName_
Name of u+ table.
Definition: nutUTabulatedWallFunctionFvPatchScalarField.H:117
U
U
Definition: pEqn.H:72
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::IOobject::MUST_READ_IF_MODIFIED
Definition: IOobject.H:186
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
Rey
scalar Rey
Definition: evaluateNearWall.H:28
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:62
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
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)
nutUTabulatedWallFunctionFvPatchScalarField.H
Foam::nutWallFunctionFvPatchScalarField
The class nutWallFunction is a base class that parents the derived boundary conditions which provide ...
Definition: nutWallFunctionFvPatchScalarField.H:251
Foam::IOobject::NO_READ
Definition: IOobject.H:188
constant
constant condensation/saturation model.
turbulenceModel.H
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
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