kLowReWallFunctionFvPatchScalarField.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) 2012-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 #include "wallFvPatch.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39 
41 {
42  if (!isA<wallFvPatch>(patch()))
43  {
45  << "Invalid wall function specification" << nl
46  << " Patch type for patch " << patch().name()
47  << " must be wall" << nl
48  << " Current patch type is " << patch().type() << nl << endl
49  << abort(FatalError);
50  }
51 }
52 
53 
54 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
55 
57 (
58  const scalar kappa,
59  const scalar E
60 )
61 {
62  scalar ypl = 11.0;
63 
64  for (int i=0; i<10; i++)
65  {
66  ypl = log(max(E*ypl, 1))/kappa;
67  }
68 
69  return ypl;
70 }
71 
72 
73 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
74 
76 (
77  const fvPatch& p,
79 )
80 :
82  Cmu_(0.09),
83  kappa_(0.41),
84  E_(9.8),
85  Ceps2_(1.9),
86  yPlusLam_(yPlusLam(kappa_, E_))
87 {
88  checkType();
89 }
90 
91 
93 (
95  const fvPatch& p,
97  const fvPatchFieldMapper& mapper
98 )
99 :
100  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
101  Cmu_(ptf.Cmu_),
102  kappa_(ptf.kappa_),
103  E_(ptf.E_),
104  Ceps2_(ptf.Ceps2_),
105  yPlusLam_(ptf.yPlusLam_)
106 {
107  checkType();
108 }
109 
110 
112 (
113  const fvPatch& p,
115  const dictionary& dict
116 )
117 :
119  Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
120  kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
121  E_(dict.lookupOrDefault<scalar>("E", 9.8)),
122  Ceps2_(dict.lookupOrDefault<scalar>("Ceps2", 1.9)),
123  yPlusLam_(yPlusLam(kappa_, E_))
124 {
125  checkType();
126 }
127 
128 
130 (
132 )
133 :
135  Cmu_(kwfpsf.Cmu_),
136  kappa_(kwfpsf.kappa_),
137  E_(kwfpsf.E_),
138  Ceps2_(kwfpsf.Ceps2_),
139  yPlusLam_(kwfpsf.yPlusLam_)
140 {
141  checkType();
142 }
143 
144 
146 (
149 )
150 :
151  fixedValueFvPatchField<scalar>(kwfpsf, iF),
152  Cmu_(kwfpsf.Cmu_),
153  kappa_(kwfpsf.kappa_),
154  E_(kwfpsf.E_),
155  Ceps2_(kwfpsf.Ceps2_),
156  yPlusLam_(kwfpsf.yPlusLam_)
157 {
158  checkType();
159 }
160 
161 
162 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
163 
165 {
166  if (updated())
167  {
168  return;
169  }
170 
171  const label patchi = patch().index();
172 
173  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
174  (
176  (
179  )
180  );
181  const scalarField& y = turbModel.y()[patchi];
182 
183  const tmp<volScalarField> tk = turbModel.k();
184  const volScalarField& k = tk();
185 
186  const tmp<scalarField> tnuw = turbModel.nu(patchi);
187  const scalarField& nuw = tnuw();
188 
189  const scalar Cmu25 = pow025(Cmu_);
190 
191  scalarField& kw = *this;
192 
193  // Set k wall values
194  forAll(kw, faceI)
195  {
196  label faceCellI = patch().faceCells()[faceI];
197 
198  scalar uTau = Cmu25*sqrt(k[faceCellI]);
199 
200  scalar yPlus = uTau*y[faceI]/nuw[faceI];
201 
202  if (yPlus > yPlusLam_)
203  {
204  scalar Ck = -0.416;
205  scalar Bk = 8.366;
206  kw[faceI] = Ck/kappa_*log(yPlus) + Bk;
207  }
208  else
209  {
210  scalar C = 11.0;
211  scalar Cf = (1.0/sqr(yPlus + C) + 2.0*yPlus/pow3(C) - 1.0/sqr(C));
212  kw[faceI] = 2400.0/sqr(Ceps2_)*Cf;
213  }
214 
215  kw[faceI] *= sqr(uTau);
216  }
217 
218  // Limit kw to avoid failure of the turbulence model due to division by kw
219  kw = max(kw, SMALL);
220 
222 
223  // TODO: perform averaging for cells sharing more than one boundary face
224 }
225 
226 
228 (
229  const Pstream::commsTypes commsType
230 )
231 {
233 }
234 
235 
237 {
238  os.writeKeyword("Cmu") << Cmu_ << token::END_STATEMENT << nl;
239  os.writeKeyword("kappa") << kappa_ << token::END_STATEMENT << nl;
240  os.writeKeyword("E") << E_ << token::END_STATEMENT << nl;
241  os.writeKeyword("Ceps2") << Ceps2_ << token::END_STATEMENT << nl;
243 }
244 
245 
246 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247 
249 (
252 );
253 
254 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
255 
256 } // End namespace Foam
257 
258 // ************************************************************************* //
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::token::END_STATEMENT
@ END_STATEMENT
Definition: token.H:99
Foam::turbulenceModel::nu
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
Foam::fvPatchField< Type >::evaluate
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::blocking)
Evaluate the patch field, sets Updated to false.
Definition: fvPatchField.C:318
p
p
Definition: pEqn.H:62
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
Foam::kLowReWallFunctionFvPatchScalarField::Cmu_
scalar Cmu_
Cmu coefficient.
Definition: kLowReWallFunctionFvPatchScalarField.H:110
Foam::fvPatchField< scalar >::dimensionedInternalField
const DimensionedField< Type, volMesh > & dimensionedInternalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:307
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:97
wallFvPatch.H
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Definition: dictionaryTemplates.C:33
Foam::kLowReWallFunctionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: kLowReWallFunctionFvPatchScalarField.C:164
Foam::kLowReWallFunctionFvPatchScalarField::kLowReWallFunctionFvPatchScalarField
kLowReWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: kLowReWallFunctionFvPatchScalarField.C:76
Foam::constant::atomic::group
const char *const group
Group name for atomic constants.
Definition: atomicConstants.C:40
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
fvPatchFieldMapper.H
Foam::kLowReWallFunctionFvPatchScalarField::kappa_
scalar kappa_
Von Karman constant.
Definition: kLowReWallFunctionFvPatchScalarField.H:113
Foam::turbulenceModel::y
const nearWallDist & y() const
Return the near wall distances.
Definition: turbulenceModel.H:157
Foam::constant::electromagnetic::kappa
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Foam::fvPatchField< scalar >::updated
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:348
Foam::fixedValueFvPatchField< scalar >
Foam::fixedValueFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fixedValueFvPatchField.C:142
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:131
Foam::fvPatch::name
const word & name() const
Return name.
Definition: fvPatch.H:149
Foam::kLowReWallFunctionFvPatchScalarField::evaluate
virtual void evaluate(const Pstream::commsTypes)
Evaluate the patchField.
Definition: kLowReWallFunctionFvPatchScalarField.C:228
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::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:87
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::kLowReWallFunctionFvPatchScalarField::Ceps2_
scalar Ceps2_
Ceps2 coefficient.
Definition: kLowReWallFunctionFvPatchScalarField.H:119
uTau
scalar uTau
Definition: evaluateNearWall.H:14
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:60
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::fvPatch::index
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:179
Foam::kLowReWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: kLowReWallFunctionFvPatchScalarField.C:236
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:253
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::turbulenceModel::k
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::makePatchTypeField
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
Foam::fvPatchField< Type >::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:296
Foam::kLowReWallFunctionFvPatchScalarField
This boundary condition provides a turbulence kinetic energy wall function condition for low- and hig...
Definition: kLowReWallFunctionFvPatchScalarField.H:101
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::kLowReWallFunctionFvPatchScalarField::checkType
virtual void checkType()
Check the type of the patch.
Definition: kLowReWallFunctionFvPatchScalarField.C:40
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:64
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::fvPatchField< scalar >::db
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:181
kLowReWallFunctionFvPatchScalarField.H
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
Foam::yPlus
This function object evaluates and outputs turbulence y+ for turbulence models. The field is stored o...
Definition: yPlus.H:110
patchi
label patchi
Definition: getPatchFieldScalar.H:1
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::kLowReWallFunctionFvPatchScalarField::yPlusLam_
scalar yPlusLam_
Y+ at the edge of the laminar sublayer.
Definition: kLowReWallFunctionFvPatchScalarField.H:122
Foam::fvPatch::faceCells
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:93
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::fvPatchField< scalar >::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:300
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Definition: objectRegistryTemplates.C:165
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::C
Graphite solid properties.
Definition: C.H:57
Foam::kLowReWallFunctionFvPatchScalarField::E_
scalar E_
E coefficient.
Definition: kLowReWallFunctionFvPatchScalarField.H:116
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::kLowReWallFunctionFvPatchScalarField::yPlusLam
scalar yPlusLam(const scalar kappa, const scalar E)
Calculate the Y+ at the edge of the laminar sublayer.
Definition: kLowReWallFunctionFvPatchScalarField.C:57
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:51
y
scalar y
Definition: LISASMDCalcMethod1.H:14