nutkFilmWallFunctionFvPatchScalarField.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 "fvPatchFieldMapper.H"
28 #include "volFields.H"
31 #include "surfaceFilmModel.H"
32 #include "mappedWallPolyPatch.H"
33 #include "mapDistribute.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace compressible
40 {
41 namespace RASModels
42 {
43 
44 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
45 
47 (
48  const scalarField& magGradU
49 ) const
50 {
51  tmp<scalarField> tuTau(new scalarField(patch().size(), 0.0));
52  scalarField& uTau = tuTau();
53 
55 
56  bool foundFilm = db().time().foundObject<modelType>(filmRegionName_);
57 
58  if (!foundFilm)
59  {
60  // do nothing on construction - film model doesn't exist yet
61  return tuTau;
62  }
63 
64  const label patchi = patch().index();
65 
66  // Retrieve phase change mass from surface film model
67  const modelType& filmModel =
68  db().time().lookupObject<modelType>(filmRegionName_);
69 
70  const label filmPatchI = filmModel.regionPatchID(patchi);
71 
72  tmp<volScalarField> mDotFilm(filmModel.primaryMassTrans());
73  scalarField mDotFilmp = mDotFilm().boundaryField()[filmPatchI];
74  filmModel.toPrimary(filmPatchI, mDotFilmp);
75 
76 
77  // Retrieve RAS turbulence model
78  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
79  (
81  (
84  )
85  );
86 
87  const scalarField& y = turbModel.y()[patchi];
88  const tmp<volScalarField> tk = turbModel.k();
89  const volScalarField& k = tk();
90  const tmp<scalarField> tnuw = turbModel.nu(patchi);
91  const scalarField& nuw = tnuw();
92 
93  const scalar Cmu25 = pow(Cmu_, 0.25);
94 
95  forAll(uTau, faceI)
96  {
97  label faceCellI = patch().faceCells()[faceI];
98 
99  scalar ut = Cmu25*sqrt(k[faceCellI]);
100 
101  scalar yPlus = y[faceI]*ut/nuw[faceI];
102 
103  scalar mStar = mDotFilmp[faceI]/(y[faceI]*ut);
104 
105  scalar factor = 0.0;
106  if (yPlus > yPlusCrit_)
107  {
108  scalar expTerm = exp(min(50.0, B_*mStar));
109  scalar powTerm = pow(yPlus, mStar/kappa_);
110  factor = mStar/(expTerm*powTerm - 1.0 + ROOTVSMALL);
111  }
112  else
113  {
114  scalar expTerm = exp(min(50.0, mStar));
115  factor = mStar/(expTerm*yPlus - 1.0 + ROOTVSMALL);
116  }
117 
118  uTau[faceI] = sqrt(max(0, magGradU[faceI]*ut*factor));
119  }
120 
121  return tuTau;
122 }
123 
124 
126 {
127  const label patchi = patch().index();
128 
129  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
130  (
132  (
135  )
136  );
137 
138  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
139  const scalarField magGradU(mag(Uw.snGrad()));
140  const tmp<scalarField> tnuw = turbModel.nu(patchi);
141  const scalarField& nuw = tnuw();
142 
143  return max
144  (
145  scalar(0),
146  sqr(calcUTau(magGradU))/(magGradU + ROOTVSMALL) - nuw
147  );
148 }
149 
150 
151 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
152 
154 (
155  const fvPatch& p,
157 )
158 :
160  filmRegionName_("surfaceFilmProperties"),
161  B_(5.5),
162  yPlusCrit_(11.05)
163 {}
164 
165 
167 (
169  const fvPatch& p,
171  const fvPatchFieldMapper& mapper
172 )
173 :
174  nutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
175  filmRegionName_(ptf.filmRegionName_),
176  B_(5.5),
177  yPlusCrit_(11.05)
178 {}
179 
180 
182 (
183  const fvPatch& p,
185  const dictionary& dict
186 )
187 :
189  filmRegionName_
190  (
191  dict.lookupOrDefault<word>("filmRegion", "surfaceFilmProperties")
192  ),
193  B_(dict.lookupOrDefault("B", 5.5)),
194  yPlusCrit_(dict.lookupOrDefault("yPlusCrit", 11.05))
195 {}
196 
197 
199 (
201 )
202 :
204  filmRegionName_(wfpsf.filmRegionName_),
205  B_(wfpsf.B_),
206  yPlusCrit_(wfpsf.yPlusCrit_)
207 {}
208 
209 
211 (
214 )
215 :
217  filmRegionName_(wfpsf.filmRegionName_),
218  B_(wfpsf.B_),
219  yPlusCrit_(wfpsf.yPlusCrit_)
220 {}
221 
222 
223 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
224 
226 {
227  const label patchi = patch().index();
228 
229  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
230  (
232  (
235  )
236  );
237 
238  const scalarField& y = turbModel.y()[patchi];
239  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
240  const tmp<scalarField> tnuw = turbModel.nu(patchi);
241  const scalarField& nuw = tnuw();
242 
243  return y*calcUTau(mag(Uw.snGrad()))/nuw;
244 }
245 
246 
248 {
250  writeEntryIfDifferent<word>
251  (
252  os,
253  "filmRegion",
254  "surfaceFilmProperties",
256  );
257  writeLocalEntries(os);
258  os.writeKeyword("B") << B_ << token::END_STATEMENT << nl;
259  os.writeKeyword("yPlusCrit") << yPlusCrit_ << token::END_STATEMENT << nl;
260  writeEntry("value", os);
261 }
262 
263 
264 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
265 
267 
268 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
269 
270 } // End namespace RASModels
271 } // End namespace compressible
272 } // End namespace Foam
273 
274 // ************************************************************************* //
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::token::END_STATEMENT
@ END_STATEMENT
Definition: token.H:99
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
mappedWallPolyPatch.H
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::write
virtual void write(Ostream &os) const
Write.
Definition: nutkFilmWallFunctionFvPatchScalarField.C:247
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::yPlusCrit_
scalar yPlusCrit_
y+ value for laminar -> turbulent transition (default = 11.05)
Definition: nutkFilmWallFunctionFvPatchScalarField.H:86
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::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::yPlus
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
Definition: nutkFilmWallFunctionFvPatchScalarField.C:225
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::calcNut
virtual tmp< scalarField > calcNut() const
Calculate the turbulence viscosity.
Definition: nutkFilmWallFunctionFvPatchScalarField.C:125
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::B_
scalar B_
B Coefficient (default = 5.5)
Definition: nutkFilmWallFunctionFvPatchScalarField.H:83
dimensionedInternalField
rDeltaT dimensionedInternalField()
Foam::compressible::RASModels::makePatchTypeField
makePatchTypeField(fvPatchScalarField, alphatFilmWallFunctionFvPatchScalarField)
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:97
Foam::regionModels::regionModel::time
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:37
Foam::constant::atomic::group
const char *const group
Group name for atomic constants.
Definition: atomicConstants.C:40
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::nutkFilmWallFunctionFvPatchScalarField
nutkFilmWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: nutkFilmWallFunctionFvPatchScalarField.C:154
fvPatchFieldMapper.H
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
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
surfaceFilmModel.H
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
uTau
scalar uTau
Definition: evaluateNearWall.H:14
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::calcUTau
virtual tmp< scalarField > calcUTau(const scalarField &magGradU) const
Calculate the friction velocity.
Definition: nutkFilmWallFunctionFvPatchScalarField.C:47
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::filmRegionName_
word filmRegionName_
Name of film region.
Definition: nutkFilmWallFunctionFvPatchScalarField.H:80
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::objectRegistry::foundObject
bool foundObject(const word &name) const
Is the named Type found?
Definition: objectRegistryTemplates.C:142
nutkFilmWallFunctionFvPatchScalarField.H
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
mapDistribute.H
Foam::ThermalDiffusivity
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Definition: ThermalDiffusivity.H:48
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
compressible
bool compressible
Definition: pEqn.H:40
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::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField
This boundary condition provides a turbulent viscosity condition when using wall functions,...
Definition: nutkFilmWallFunctionFvPatchScalarField.H:71
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::regionModels::surfaceFilmModels::surfaceFilmModel
Base class for surface film models.
Definition: surfaceFilmModel.H:61
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
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
turbulentFluidThermoModel.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
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::nutkWallFunctionFvPatchScalarField
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions,...
Definition: nutkWallFunctionFvPatchScalarField.H:66