fWallFunctionFvPatchScalarField.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-2016 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 "fvPatchFieldMapper.H"
28 #include "volFields.H"
29 #include "wallFvPatch.H"
30 #include "v2f.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace RASModels
38 {
39 
40 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41 
43 {
44  if (!isA<wallFvPatch>(patch()))
45  {
47  << "Invalid wall function specification" << nl
48  << " Patch type for patch " << patch().name()
49  << " must be wall" << nl
50  << " Current patch type is " << patch().type() << nl << endl
51  << abort(FatalError);
52  }
53 }
54 
55 
57 {
58  os.writeKeyword("Cmu") << Cmu_ << token::END_STATEMENT << nl;
59  os.writeKeyword("kappa") << kappa_ << token::END_STATEMENT << nl;
60  os.writeKeyword("E") << E_ << token::END_STATEMENT << nl;
61 }
62 
63 
65 (
66  const scalar kappa,
67  const scalar E
68 )
69 {
70  scalar ypl = 11.0;
71 
72  for (int i=0; i<10; i++)
73  {
74  ypl = log(max(E*ypl, 1))/kappa;
75  }
76 
77  return ypl;
78 }
79 
80 
81 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
82 
84 (
85  const fvPatch& p,
87 )
88 :
90  Cmu_(0.09),
91  kappa_(0.41),
92  E_(9.8),
94 {
95  checkType();
96 }
97 
98 
100 (
102  const fvPatch& p,
104  const fvPatchFieldMapper& mapper
105 )
106 :
107  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
108  Cmu_(ptf.Cmu_),
109  kappa_(ptf.kappa_),
110  E_(ptf.E_),
111  yPlusLam_(ptf.yPlusLam_)
112 {
113  checkType();
114 }
115 
116 
118 (
119  const fvPatch& p,
121  const dictionary& dict
122 )
123 :
125  Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
126  kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
127  E_(dict.lookupOrDefault<scalar>("E", 9.8)),
129 {
130  checkType();
131 }
132 
133 
135 (
136  const fWallFunctionFvPatchScalarField& v2wfpsf
137 )
138 :
140  Cmu_(v2wfpsf.Cmu_),
141  kappa_(v2wfpsf.kappa_),
142  E_(v2wfpsf.E_),
143  yPlusLam_(v2wfpsf.yPlusLam_)
144 {
145  checkType();
146 }
147 
148 
150 (
151  const fWallFunctionFvPatchScalarField& v2wfpsf,
153 )
154 :
155  fixedValueFvPatchField<scalar>(v2wfpsf, iF),
156  Cmu_(v2wfpsf.Cmu_),
157  kappa_(v2wfpsf.kappa_),
158  E_(v2wfpsf.E_),
159  yPlusLam_(v2wfpsf.yPlusLam_)
160 {
161  checkType();
162 }
163 
164 
165 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
166 
168 {
169  if (updated())
170  {
171  return;
172  }
173 
174  const label patchi = patch().index();
175 
176  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
177  (
179  (
181  internalField().group()
182  )
183  );
184  const v2fBase& v2fModel = refCast<const v2fBase>(turbModel);
185 
186  const scalarField& y = turbModel.y()[patchi];
187 
188  const tmp<volScalarField> tk = turbModel.k();
189  const volScalarField& k = tk();
190 
191  const tmp<volScalarField> tepsilon = turbModel.epsilon();
192  const volScalarField& epsilon = tepsilon();
193 
194  const tmp<volScalarField> tv2 = v2fModel.v2();
195  const volScalarField& v2 = tv2();
196 
197  const tmp<scalarField> tnuw = turbModel.nu(patchi);
198  const scalarField& nuw = tnuw();
199 
200  const scalar Cmu25 = pow025(Cmu_);
201 
202  scalarField& f = *this;
203 
204  // Set f wall values
205  forAll(f, facei)
206  {
207  label faceCelli = patch().faceCells()[facei];
208 
209  scalar uTau = Cmu25*sqrt(k[faceCelli]);
210 
211  scalar yPlus = uTau*y[facei]/nuw[facei];
212 
213  if (yPlus > yPlusLam_)
214  {
215  scalar N = 6.0;
216  scalar v2c = v2[faceCelli];
217  scalar epsc = epsilon[faceCelli];
218  scalar kc = k[faceCelli];
219 
220  f[facei] = N*v2c*epsc/(sqr(kc) + ROOTVSMALL);
221  f[facei] /= sqr(uTau) + ROOTVSMALL;
222  }
223  else
224  {
225  f[facei] = 0.0;
226  }
227  }
228 
230 
231  // TODO: perform averaging for cells sharing more than one boundary face
232 }
233 
234 
236 (
237  const Pstream::commsTypes commsType
238 )
239 {
241 }
242 
243 
245 {
246  writeLocalEntries(os);
248 }
249 
250 
251 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
252 
254 (
257 );
258 
259 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
260 
261 } // End namespace RASModels
262 } // End namespace Foam
263 
264 // ************************************************************************* //
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:339
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::blocking)
Evaluate the patch field, sets Updated to false.
Definition: fvPatchField.C:332
const char *const group
Group name for atomic constants.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
dimensionedScalar log(const dimensionedScalar &ds)
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual void write(Ostream &) const
Write.
scalar yPlusLam(const scalar kappa, const scalar E)
Calculate the Y+ at the edge of the laminar sublayer.
commsTypes
Types of communications.
Definition: UPstream.H:64
virtual void evaluate(const Pstream::commsTypes)
Evaluate the patchField.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
const word & name() const
Return name.
Definition: fvPatch.H:149
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
dimensionedScalar pow025(const dimensionedScalar &ds)
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
label k
Boltzmann constant.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
scalar uTau
scalar yPlusLam_
Y+ at the edge of the laminar sublayer.
Abstract base-class for v2-f models to provide BCs access to the v2 and f fields. ...
Definition: v2fBase.H:57
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
scalar y
static const word propertiesName
Default name of the turbulence properties dictionary.
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:370
errorManip< error > abort(error &err)
Definition: errorManip.H:131
makePatchTypeField(fvPatchScalarField, fWallFunctionFvPatchScalarField)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:93
static const char nl
Definition: Ostream.H:262
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
labelList f(nPoints)
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:198
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:179
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
fWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
label patchi
This boundary condition provides a turbulence damping function, f, wall function condition for low- a...
scalar yPlus
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:312
scalar epsilon
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
label N
Definition: createFields.H:22
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:54
virtual void checkType()
Check the type of the patch.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
const nearWallDist & y() const
Return the near wall distances.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:346
Namespace for OpenFOAM.