nutUBlendedWallFunctionFvPatchScalarField.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) 2016-2020 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
29 #include "turbulenceModel.H"
30 #include "fvPatchFieldMapper.H"
31 #include "volFields.H"
33 
34 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
35 
38 {
39  const label patchi = patch().index();
40 
41  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
42  (
44  (
46  internalField().group()
47  )
48  );
49  const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
50  const scalarField magGradU(mag(Uw.snGrad()));
51  const tmp<scalarField> tnuw = turbModel.nu(patchi);
52  const scalarField& nuw = tnuw();
53 
54  return max
55  (
56  scalar(0),
57  sqr(calcUTau(magGradU))/(magGradU + ROOTVSMALL) - nuw
58  );
59 }
60 
61 
64 (
65  const scalarField& magGradU
66 ) const
67 {
68  const label patchi = patch().index();
69 
70  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
71  (
73  (
75  internalField().group()
76  )
77  );
78 
79  const scalarField& y = turbModel.y()[patchi];
80 
81  const tmp<scalarField> tnuw = turbModel.nu(patchi);
82  const scalarField& nuw = tnuw();
83 
84  const vectorField n(patch().nf());
85  const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
86  vectorField Up(Uw.patchInternalField() - Uw);
87  Up -= n*(n & Up);
88  const scalarField magUp(mag(Up));
89 
90  tmp<scalarField> tuTaup(new scalarField(patch().size(), Zero));
91  scalarField& uTaup = tuTaup.ref();
92 
93  const scalarField& nutw = *this;
94 
95  forAll(uTaup, facei)
96  {
97  scalar ut = sqrt((nutw[facei] + nuw[facei])*magGradU[facei]);
98  if (mag(ut) > ROOTVSMALL)
99  {
100  scalar error = GREAT;
101  label iter = 0;
102  while (iter++ < 10 && error > 0.001)
103  {
104  const scalar yPlus = y[facei]*ut/nuw[facei];
105  const scalar uTauVis = magUp[facei]/yPlus;
106  const scalar uTauLog = kappa_*magUp[facei]/log(E_*yPlus);
107 
108  const scalar utNew =
109  pow(pow(uTauVis, n_) + pow(uTauLog, n_), 1.0/n_);
110  error = mag(ut - utNew)/(ut + ROOTVSMALL);
111  ut = 0.5*(ut + utNew);
112  }
113  }
114  uTaup[facei] = ut;
115  }
116 
117  return tuTaup;
118 }
119 
120 
121 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
122 
125 (
126  const fvPatch& p,
128 )
129 :
131  n_(4)
132 {}
133 
134 
137 (
139  const fvPatch& p,
141  const fvPatchFieldMapper& mapper
142 )
143 :
144  nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
145  n_(ptf.n_)
146 {}
147 
148 
151 (
152  const fvPatch& p,
154  const dictionary& dict
155 )
156 :
158  n_(dict.getOrDefault<scalar>("n", 4.0))
159 {}
160 
161 
164 (
166 )
167 :
169  n_(wfpsf.n_)
170 {}
171 
172 
175 (
178 )
179 :
181  n_(wfpsf.n_)
182 {}
183 
184 
185 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
186 
189 {
190  const label patchi = patch().index();
191 
192  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
193  (
195  (
197  internalField().group()
198  )
199  );
200  const scalarField& y = turbModel.y()[patchi];
201  const tmp<scalarField> tnuw = turbModel.nu(patchi);
202  const scalarField& nuw = tnuw();
203  const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
204  const scalarField magGradU(mag(Uw.snGrad()));
205 
206  return y*calcUTau(magGradU)/nuw;
207 }
208 
209 
211 (
212  Ostream& os
213 ) const
214 {
216  writeLocalEntries(os);
217  os.writeEntry("n", n_);
218  writeEntry("value", os);
219 }
220 
221 
222 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223 
224 namespace Foam
225 {
227  (
230  );
231 }
232 
233 
234 // ************************************************************************* //
Foam::fvPatchField< vector >
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
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::nutUBlendedWallFunctionFvPatchScalarField::write
virtual void write(Ostream &os) const
Write.
Definition: nutUBlendedWallFunctionFvPatchScalarField.C:211
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::nutUBlendedWallFunctionFvPatchScalarField::nutUBlendedWallFunctionFvPatchScalarField
nutUBlendedWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: nutUBlendedWallFunctionFvPatchScalarField.C:125
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
fvPatchFieldMapper.H
Foam::turbulenceModel::y
const nearWallDist & y() const
Return the near wall distances.
Definition: turbulenceModel.H:159
magUp
scalar magUp
Definition: evaluateNearWall.H:10
Foam::nutUBlendedWallFunctionFvPatchScalarField::yPlus
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
Definition: nutUBlendedWallFunctionFvPatchScalarField.C:188
Foam::nutUBlendedWallFunctionFvPatchScalarField::calcUTau
virtual tmp< scalarField > calcUTau(const scalarField &magGradU) const
Calculate the friction velocity.
Definition: nutUBlendedWallFunctionFvPatchScalarField.C:64
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::nutUBlendedWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the turbulent viscosity, i.e....
Definition: nutUBlendedWallFunctionFvPatchScalarField.H:147
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::Field< scalar >
nutUBlendedWallFunctionFvPatchScalarField.H
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::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
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::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
U
U
Definition: pEqn.H:72
Foam::nutUBlendedWallFunctionFvPatchScalarField::n_
scalar n_
Model coefficient; default = 4.
Definition: nutUBlendedWallFunctionFvPatchScalarField.H:156
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::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
yPlus
scalar yPlus
Definition: evaluateNearWall.H:16
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)
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::nutWallFunctionFvPatchScalarField
The class nutWallFunction is a base class that parents the derived boundary conditions which provide ...
Definition: nutWallFunctionFvPatchScalarField.H:251
Foam::nutUBlendedWallFunctionFvPatchScalarField::calcNut
virtual tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
Definition: nutUBlendedWallFunctionFvPatchScalarField.C:37
Foam::error
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Definition: error.H:73
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