nutkFilmWallFunctionFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 "surfaceFilmRegionModel.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.ref();
53 
55 
56  bool foundFilm =
57  db().time().foundObject<modelType>("surfaceFilmProperties");
58 
59  if (!foundFilm)
60  {
61  // Do nothing on construction - film model doesn't exist yet
62  return tuTau;
63  }
64 
65  const label patchi = patch().index();
66 
67  // Retrieve phase change mass from surface film model
68  const modelType& filmModel =
69  db().time().lookupObject<modelType>("surfaceFilmProperties");
70 
71  const label filmPatchi = filmModel.regionPatchID(patchi);
72 
73  tmp<volScalarField> mDotFilm(filmModel.primaryMassTrans());
74  scalarField mDotFilmp = mDotFilm().boundaryField()[filmPatchi];
75  filmModel.toPrimary(filmPatchi, mDotFilmp);
76 
77 
78  // Retrieve RAS turbulence model
79  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
80  (
82  (
84  internalField().group()
85  )
86  );
87 
88  const scalarField& y = turbModel.y()[patchi];
89  const tmp<volScalarField> tk = turbModel.k();
90  const volScalarField& k = tk();
91  const tmp<scalarField> tnuw = turbModel.nu(patchi);
92  const scalarField& nuw = tnuw();
93 
94  const scalar Cmu25 = pow(Cmu_, 0.25);
95 
96  forAll(uTau, facei)
97  {
98  label faceCelli = patch().faceCells()[facei];
99 
100  scalar ut = Cmu25*sqrt(k[faceCelli]);
101 
102  scalar yPlus = y[facei]*ut/nuw[facei];
103 
104  scalar mStar = mDotFilmp[facei]/(y[facei]*ut);
105 
106  scalar factor = 0.0;
107  if (yPlus > yPlusCrit_)
108  {
109  scalar expTerm = exp(min(50.0, B_*mStar));
110  scalar powTerm = pow(yPlus, mStar/kappa_);
111  factor = mStar/(expTerm*powTerm - 1.0 + rootVSmall);
112  }
113  else
114  {
115  scalar expTerm = exp(min(50.0, mStar));
116  factor = mStar/(expTerm*yPlus - 1.0 + rootVSmall);
117  }
118 
119  uTau[facei] = sqrt(max(0, magGradU[facei]*ut*factor));
120  }
121 
122  return tuTau;
123 }
124 
125 
127 {
128  const label patchi = patch().index();
129 
130  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
131  (
133  (
135  internalField().group()
136  )
137  );
138 
139  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
140  const scalarField magGradU(mag(Uw.snGrad()));
141  const tmp<scalarField> tnuw = turbModel.nu(patchi);
142  const scalarField& nuw = tnuw();
143 
144  return max
145  (
146  scalar(0),
147  sqr(calcUTau(magGradU))/(magGradU + rootVSmall) - nuw
148  );
149 }
150 
151 
152 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
153 
155 (
156  const fvPatch& p,
158 )
159 :
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  B_(5.5),
176  yPlusCrit_(11.05)
177 {}
178 
179 
181 (
182  const fvPatch& p,
184  const dictionary& dict
185 )
186 :
188  B_(dict.lookupOrDefault("B", 5.5)),
189  yPlusCrit_(dict.lookupOrDefault("yPlusCrit", 11.05))
190 {}
191 
192 
194 (
196 )
197 :
199  B_(wfpsf.B_),
200  yPlusCrit_(wfpsf.yPlusCrit_)
201 {}
202 
203 
205 (
208 )
209 :
211  B_(wfpsf.B_),
212  yPlusCrit_(wfpsf.yPlusCrit_)
213 {}
214 
215 
216 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
217 
219 {
220  const label patchi = patch().index();
221 
222  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
223  (
225  (
227  internalField().group()
228  )
229  );
230 
231  const scalarField& y = turbModel.y()[patchi];
232  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
233  const tmp<scalarField> tnuw = turbModel.nu(patchi);
234  const scalarField& nuw = tnuw();
235 
236  return y*calcUTau(mag(Uw.snGrad()))/nuw;
237 }
238 
239 
241 {
243  writeLocalEntries(os);
244  os.writeKeyword("B") << B_ << token::END_STATEMENT << nl;
245  os.writeKeyword("yPlusCrit") << yPlusCrit_ << token::END_STATEMENT << nl;
246  writeEntry("value", os);
247 }
248 
249 
250 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251 
253 
254 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
255 
256 } // End namespace RASModels
257 } // End namespace compressible
258 } // End namespace Foam
259 
260 // ************************************************************************* //
const char *const group
Group name for atomic constants.
#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
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:216
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
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 > &)
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
bool foundObject(const word &name) const
Is the named Type found?
makePatchTypeField(fvPatchScalarField, alphatFilmWallFunctionFvPatchScalarField)
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
label k
Boltzmann constant.
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:37
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:362
Macros for easy insertion into run-time selection tables.
scalar yPlusCrit_
y+ value for laminar -> turbulent transition (default = 11.05)
scalar uTau
scalar y
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
dimensionedScalar exp(const dimensionedScalar &ds)
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual tmp< scalarField > calcNut() const
Calculate the turbulence viscosity.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:265
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
nutkFilmWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
scalar yPlus
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< scalarField > calcUTau(const scalarField &magGradU) const
Calculate the friction velocity.
This boundary condition provides a turbulent viscosity condition when using wall functions, based on turbulence kinetic energy, for use with surface film models.
nutkWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Namespace for OpenFOAM.