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-2022 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 "distributionMap.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 momentumTransportModel& turbModel =
80  db().lookupObject<momentumTransportModel>
81  (
83  (
84  momentumTransportModel::typeName,
85  internalField().group()
86  )
87  );
88 
89  const scalarField& y = turbModel.y()[patchi];
90  const tmp<volScalarField> tk = turbModel.k();
91  const volScalarField& k = tk();
92  const tmp<scalarField> tnuw = turbModel.nu(patchi);
93  const scalarField& nuw = tnuw();
94 
95  const scalar Cmu25 = pow(Cmu_, 0.25);
96 
97  forAll(uTau, facei)
98  {
99  label faceCelli = patch().faceCells()[facei];
100 
101  scalar ut = Cmu25*sqrt(k[faceCelli]);
102 
103  scalar yPlus = y[facei]*ut/nuw[facei];
104 
105  scalar mStar = mDotFilmp[facei]/(y[facei]*ut);
106 
107  scalar factor = 0.0;
108  if (yPlus > yPlusCrit_)
109  {
110  scalar expTerm = exp(min(50.0, B_*mStar));
111  scalar powTerm = pow(yPlus, mStar/kappa_);
112  factor = mStar/(expTerm*powTerm - 1.0 + rootVSmall);
113  }
114  else
115  {
116  scalar expTerm = exp(min(50.0, mStar));
117  factor = mStar/(expTerm*yPlus - 1.0 + rootVSmall);
118  }
119 
120  uTau[facei] = sqrt(max(0, magGradU[facei]*ut*factor));
121  }
122 
123  return tuTau;
124 }
125 
126 
128 {
129  const label patchi = patch().index();
130 
131  const momentumTransportModel& turbModel =
132  db().lookupObject<momentumTransportModel>
133  (
135  (
136  momentumTransportModel::typeName,
137  internalField().group()
138  )
139  );
140 
141  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
142  const scalarField magGradU(mag(Uw.snGrad()));
143  const tmp<scalarField> tnuw = turbModel.nu(patchi);
144  const scalarField& nuw = tnuw();
145 
146  return max
147  (
148  scalar(0),
149  sqr(calcUTau(magGradU))/(magGradU + rootVSmall) - nuw
150  );
151 }
152 
153 
154 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
155 
157 (
158  const fvPatch& p,
160 )
161 :
163  B_(5.5),
164  yPlusCrit_(11.05)
165 {}
166 
167 
169 (
170  const fvPatch& p,
172  const dictionary& dict
173 )
174 :
176  B_(dict.lookupOrDefault("B", 5.5)),
177  yPlusCrit_(dict.lookupOrDefault("yPlusCrit", 11.05))
178 {}
179 
180 
182 (
184  const fvPatch& p,
186  const fvPatchFieldMapper& mapper
187 )
188 :
189  nutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
190  B_(5.5),
191  yPlusCrit_(11.05)
192 {}
193 
194 
196 (
199 )
200 :
202  B_(wfpsf.B_),
203  yPlusCrit_(wfpsf.yPlusCrit_)
204 {}
205 
206 
207 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
208 
210 {
211  const label patchi = patch().index();
212 
213  const momentumTransportModel& turbModel =
214  db().lookupObject<momentumTransportModel>
215  (
217  (
218  momentumTransportModel::typeName,
219  internalField().group()
220  )
221  );
222 
223  const scalarField& y = turbModel.y()[patchi];
224  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
225  const tmp<scalarField> tnuw = turbModel.nu(patchi);
226  const scalarField& nuw = tnuw();
227 
228  return y*calcUTau(mag(Uw.snGrad()))/nuw;
229 }
230 
231 
233 {
235  writeLocalEntries(os);
236  writeEntry(os, "B", B_);
237  writeEntry(os, "yPlusCrit", yPlusCrit_);
238  writeEntry(os, "value", *this);
239 }
240 
241 
242 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
243 
245 
246 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247 
248 } // End namespace RASModels
249 } // End namespace compressible
250 } // End namespace Foam
251 
252 // ************************************************************************* //
const char *const group
Group name for atomic constants.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:166
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:156
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
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.
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:37
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:243
const volScalarField::Boundary & y() const
Return the near wall distance.
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)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const volVectorField & U() const
Access function to velocity field.
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 > &)
virtual tmp< scalarField > nut() const
Calculate the turbulence viscosity.
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.
Base class for single-phase compressible turbulence models.
nutkWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Namespace for OpenFOAM.