nutUSpaldingWallFunctionFvPatchScalarField.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 | Copyright (C) 2011-2015 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 
27 #include "turbulenceModel.H"
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38 
40 {
41  const label patchi = patch().index();
42 
43  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
44  (
46  (
49  )
50  );
51  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
52  const scalarField magGradU(mag(Uw.snGrad()));
53  const tmp<scalarField> tnuw = turbModel.nu(patchi);
54  const scalarField& nuw = tnuw();
55 
56  return max
57  (
58  scalar(0),
59  sqr(calcUTau(magGradU))/(magGradU + ROOTVSMALL) - nuw
60  );
61 }
62 
63 
65 (
66  const scalarField& magGradU
67 ) const
68 {
69  const label patchi = patch().index();
70 
71  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
72  (
74  (
77  )
78  );
79  const scalarField& y = turbModel.y()[patchi];
80 
81  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
82  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
83 
84  const tmp<scalarField> tnuw = turbModel.nu(patchi);
85  const scalarField& nuw = tnuw();
86 
87  const scalarField& nutw = *this;
88 
89  tmp<scalarField> tuTau(new scalarField(patch().size(), 0.0));
90  scalarField& uTau = tuTau();
91 
92  forAll(uTau, faceI)
93  {
94  scalar ut = sqrt((nutw[faceI] + nuw[faceI])*magGradU[faceI]);
95 
96  if (ut > ROOTVSMALL)
97  {
98  int iter = 0;
99  scalar err = GREAT;
100 
101  do
102  {
103  scalar kUu = min(kappa_*magUp[faceI]/ut, 50);
104  scalar fkUu = exp(kUu) - 1 - kUu*(1 + 0.5*kUu);
105 
106  scalar f =
107  - ut*y[faceI]/nuw[faceI]
108  + magUp[faceI]/ut
109  + 1/E_*(fkUu - 1.0/6.0*kUu*sqr(kUu));
110 
111  scalar df =
112  y[faceI]/nuw[faceI]
113  + magUp[faceI]/sqr(ut)
114  + 1/E_*kUu*fkUu/ut;
115 
116  scalar uTauNew = ut + f/df;
117  err = mag((ut - uTauNew)/ut);
118  ut = uTauNew;
119 
120  } while (ut > ROOTVSMALL && err > 0.01 && ++iter < 10);
121 
122  uTau[faceI] = max(0.0, ut);
123  }
124  }
125 
126  return tuTau;
127 }
128 
129 
130 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
131 
134 (
135  const fvPatch& p,
137 )
138 :
140 {}
141 
142 
145 (
147  const fvPatch& p,
149  const fvPatchFieldMapper& mapper
150 )
151 :
152  nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper)
153 {}
154 
155 
158 (
159  const fvPatch& p,
161  const dictionary& dict
162 )
163 :
165 {}
166 
167 
170 (
172 )
173 :
175 {}
176 
177 
180 (
183 )
184 :
186 {}
187 
188 
189 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
190 
192 {
193  const label patchi = patch().index();
194 
195  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
196  (
198  (
201  )
202  );
203  const scalarField& y = turbModel.y()[patchi];
204  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
205  const tmp<scalarField> tnuw = turbModel.nu(patchi);
206  const scalarField& nuw = tnuw();
207 
208  return y*calcUTau(mag(Uw.snGrad()))/nuw;
209 }
210 
211 
213 {
215  writeLocalEntries(os);
216  writeEntry("value", os);
217 }
218 
219 
220 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
221 
223 (
226 );
227 
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229 
230 } // End namespace Foam
231 
232 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:349
Foam::fvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:200
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
Foam::turbulenceModel::nu
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
p
p
Definition: pEqn.H:62
Foam::nutUSpaldingWallFunctionFvPatchScalarField::yPlus
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
Definition: nutUSpaldingWallFunctionFvPatchScalarField.C:191
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
dimensionedInternalField
rDeltaT dimensionedInternalField()
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:97
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::constant::atomic::group
const char *const group
Group name for atomic constants.
Definition: atomicConstants.C:40
fvPatchFieldMapper.H
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
Foam::turbulenceModel::y
const nearWallDist & y() const
Return the near wall distances.
Definition: turbulenceModel.H:157
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
magUp
scalar magUp
Definition: evaluateNearWall.H:10
Foam::nutUSpaldingWallFunctionFvPatchScalarField::calcUTau
virtual tmp< scalarField > calcUTau(const scalarField &magGradU) const
Calculate the friction velocity.
Definition: nutUSpaldingWallFunctionFvPatchScalarField.C:65
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::fvPatchField::patchInternalField
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:208
uTau
scalar uTau
Definition: evaluateNearWall.H:14
Foam::nutUSpaldingWallFunctionFvPatchScalarField::write
virtual void write(Ostream &os) const
Write.
Definition: nutUSpaldingWallFunctionFvPatchScalarField.C:212
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:60
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::nutUSpaldingWallFunctionFvPatchScalarField
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions ...
Definition: nutUSpaldingWallFunctionFvPatchScalarField.H:91
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::makePatchTypeField
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::nutUSpaldingWallFunctionFvPatchScalarField::calcNut
virtual tmp< scalarField > calcNut() const
Calculate the turbulence viscosity.
Definition: nutUSpaldingWallFunctionFvPatchScalarField.C:39
f
labelList f(nPoints)
Foam::nutUSpaldingWallFunctionFvPatchScalarField::nutUSpaldingWallFunctionFvPatchScalarField
nutUSpaldingWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: nutUSpaldingWallFunctionFvPatchScalarField.C:134
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
patchi
label patchi
Definition: getPatchFieldScalar.H:1
nutUSpaldingWallFunctionFvPatchScalarField.H
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::nutWallFunctionFvPatchScalarField
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions,...
Definition: nutWallFunctionFvPatchScalarField.H:101
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
turbulenceModel.H
Foam::nutWallFunctionFvPatchScalarField::writeLocalEntries
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
Definition: nutWallFunctionFvPatchScalarField.C:58
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
Foam::turbulenceModel::U
const volVectorField & U() const
Access function to velocity field.
Definition: turbulenceModel.H:142
y
scalar y
Definition: LISASMDCalcMethod1.H:14