v2WallFunctionFvPatchScalarField.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 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),
93  yPlusLam_(yPlusLam(kappa_, E_))
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)),
128  yPlusLam_(yPlusLam(kappa_, E_))
129 {
130  checkType();
131 }
132 
133 
135 (
136  const v2WallFunctionFvPatchScalarField& 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 v2WallFunctionFvPatchScalarField& 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  (
182  )
183  );
184  const scalarField& y = turbModel.y()[patchi];
185 
186  const tmp<volScalarField> tk = turbModel.k();
187  const volScalarField& k = tk();
188 
189  const tmp<scalarField> tnuw = turbModel.nu(patchi);
190  const scalarField& nuw = tnuw();
191 
192  const scalar Cmu25 = pow025(Cmu_);
193 
194  scalarField& v2 = *this;
195 
196  // Set v2 wall values
197  forAll(v2, faceI)
198  {
199  label faceCellI = patch().faceCells()[faceI];
200 
201  scalar uTau = Cmu25*sqrt(k[faceCellI]);
202 
203  scalar yPlus = uTau*y[faceI]/nuw[faceI];
204 
205  if (yPlus > yPlusLam_)
206  {
207  scalar Cv2 = 0.193;
208  scalar Bv2 = -0.94;
209  v2[faceI] = Cv2/kappa_*log(yPlus) + Bv2;
210  }
211  else
212  {
213  scalar Cv2 = 0.193;
214  v2[faceI] = Cv2*pow4(yPlus);
215  }
216 
217  v2[faceI] *= sqr(uTau);
218  }
219 
221 
222  // TODO: perform averaging for cells sharing more than one boundary face
223 }
224 
225 
227 (
228  const Pstream::commsTypes commsType
229 )
230 {
232 }
233 
234 
236 {
237  writeLocalEntries(os);
239 }
240 
241 
242 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
243 
245 (
248 );
249 
250 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251 
252 } // End namespace RASModels
253 } // End namespace Foam
254 
255 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Foam::RASModels::v2WallFunctionFvPatchScalarField::kappa_
scalar kappa_
Von Karman constant.
Definition: v2WallFunctionFvPatchScalarField.H:111
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::RASModels::v2WallFunctionFvPatchScalarField::writeLocalEntries
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
Definition: v2WallFunctionFvPatchScalarField.C:56
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
Foam::RASModels::v2WallFunctionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: v2WallFunctionFvPatchScalarField.C:167
p
p
Definition: pEqn.H:62
Foam::RASModels::v2WallFunctionFvPatchScalarField::E_
scalar E_
E coefficient.
Definition: v2WallFunctionFvPatchScalarField.H:114
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::RASModels::v2WallFunctionFvPatchScalarField::Cmu_
scalar Cmu_
Cmu coefficient.
Definition: v2WallFunctionFvPatchScalarField.H:108
Foam::fvPatchField< scalar >::dimensionedInternalField
const DimensionedField< Type, volMesh > & dimensionedInternalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:307
Foam::RASModels::v2WallFunctionFvPatchScalarField::v2WallFunctionFvPatchScalarField
v2WallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: v2WallFunctionFvPatchScalarField.C:84
Foam::RASModels::v2WallFunctionFvPatchScalarField::checkType
virtual void checkType()
Check the type of the patch.
Definition: v2WallFunctionFvPatchScalarField.C:42
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:97
wallFvPatch.H
Foam::RASModels::v2WallFunctionFvPatchScalarField::yPlusLam_
scalar yPlusLam_
Y+ at the edge of the laminar sublayer.
Definition: v2WallFunctionFvPatchScalarField.H:117
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::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::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:98
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::nl
static const char nl
Definition: Ostream.H:260
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
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::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::fvPatchField< Type >::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:296
Foam::RASModels::v2WallFunctionFvPatchScalarField::yPlusLam
scalar yPlusLam(const scalar kappa, const scalar E)
Calculate the Y+ at the edge of the laminar sublayer.
Definition: v2WallFunctionFvPatchScalarField.C:65
Foam::RASModels::v2WallFunctionFvPatchScalarField::evaluate
virtual void evaluate(const Pstream::commsTypes)
Evaluate the patchField.
Definition: v2WallFunctionFvPatchScalarField.C:227
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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::RASModels::makePatchTypeField
makePatchTypeField(fvPatchScalarField, fWallFunctionFvPatchScalarField)
Foam::fvPatchField< scalar >::db
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:181
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::RASModels::v2WallFunctionFvPatchScalarField
This boundary condition provides a turbulence stress normal to streamlines wall function condition fo...
Definition: v2WallFunctionFvPatchScalarField.H:99
Foam::fvPatch::faceCells
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:93
Foam::RASModels::v2WallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: v2WallFunctionFvPatchScalarField.C:235
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
v2WallFunctionFvPatchScalarField.H
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::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
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