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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2020-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "fvPatchFieldMapper.H"
31 #include "volFields.H"
34 #include "surfaceFilmRegionModel.H"
35 #include "mappedWallPolyPatch.H"
36 #include "mapDistribute.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 namespace compressible
43 {
44 namespace RASModels
45 {
46 
47 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
48 
50 (
51  const scalarField& magGradU
52 ) const
53 {
54  auto tuTau = tmp<scalarField>::New(patch().size(), Zero);
55  auto& uTau = tuTau.ref();
56 
57  const auto* filmModelPtr = db().time().findObject
60 
61  if (!filmModelPtr)
62  {
63  // Do nothing on construction - film model doesn't exist yet
64  return tuTau;
65  }
66 
67  const label patchi = patch().index();
68 
69  // Retrieve phase change mass from surface film model
70  const auto& filmModel = *filmModelPtr;
71 
72  const label filmPatchi = filmModel.regionPatchID(patchi);
73 
74  tmp<volScalarField> mDotFilm = filmModel.primaryMassTrans();
75  scalarField mDotFilmp = mDotFilm().boundaryField()[filmPatchi];
76  filmModel.toPrimary(filmPatchi, mDotFilmp);
77 
78 
79  // Retrieve RAS turbulence model
80  const auto& turbModel = db().lookupObject<turbulenceModel>
81  (
83  (
85  internalField().group()
86  )
87  );
88 
89  const scalarField& y = turbModel.y()[patchi];
90 
91  const tmp<volScalarField> tk = turbModel.k();
92  const volScalarField& k = tk();
93 
94  const tmp<scalarField> tnuw = turbModel.nu(patchi);
95  const scalarField& nuw = tnuw();
96 
97  const scalar Cmu25 = pow025(wallCoeffs_.Cmu());
98  const scalar kappa = wallCoeffs_.kappa();
99 
100  forAll(uTau, facei)
101  {
102  const label faceCelli = patch().faceCells()[facei];
103 
104  const scalar ut = Cmu25*sqrt(k[faceCelli]);
105 
106  const scalar yPlus = y[facei]*ut/nuw[facei];
107 
108  const scalar mStar = mDotFilmp[facei]/(y[facei]*ut);
109 
110  scalar factor = 0;
111  if (yPlus > yPlusCrit_)
112  {
113  const scalar expTerm = exp(min(scalar(50), B_*mStar));
114  const scalar powTerm = pow(yPlus, mStar/kappa);
115  factor = mStar/(expTerm*powTerm - 1.0 + ROOTVSMALL);
116  }
117  else
118  {
119  const scalar expTerm = exp(min(scalar(50), mStar));
120 
121  factor = mStar/(expTerm*yPlus - 1.0 + ROOTVSMALL);
122  }
123 
124  uTau[facei] = sqrt(max(scalar(0), magGradU[facei]*ut*factor));
125  }
126 
127  return tuTau;
128 }
129 
130 
132 {
133  const label patchi = patch().index();
134 
135  const auto& turbModel = db().lookupObject<turbulenceModel>
136  (
138  (
140  internalField().group()
141  )
142  );
143 
144  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
145  const scalarField magGradU(mag(Uw.snGrad()));
146 
147  const tmp<scalarField> tnuw = turbModel.nu(patchi);
148  const scalarField& nuw = tnuw();
149 
150  return max
151  (
152  scalar(0),
153  sqr(calcUTau(magGradU))/(magGradU + ROOTVSMALL) - nuw
154  );
155 }
156 
157 
159 (
160  Ostream& os
161 ) const
162 {
163  os.writeEntryIfDifferent<word>
164  (
165  "filmRegion",
166  "surfaceFilmProperties",
168  );
169  os.writeEntryIfDifferent<scalar>("B", 5.5, B_);
170  os.writeEntryIfDifferent<scalar>("yPlusCrit", 11.05, yPlusCrit_);
171 }
172 
173 
174 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
175 
177 (
178  const fvPatch& p,
180 )
181 :
183  filmRegionName_("surfaceFilmProperties"),
184  B_(5.5),
185  yPlusCrit_(11.05)
186 {}
187 
188 
190 (
192  const fvPatch& p,
194  const fvPatchFieldMapper& mapper
195 )
196 :
197  nutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
198  filmRegionName_(ptf.filmRegionName_),
199  B_(5.5),
200  yPlusCrit_(11.05)
201 {}
202 
203 
205 (
206  const fvPatch& p,
208  const dictionary& dict
209 )
210 :
212  filmRegionName_
213  (
214  dict.getOrDefault<word>("filmRegion", "surfaceFilmProperties")
215  ),
216  B_(dict.getOrDefault("B", 5.5)),
217  yPlusCrit_(dict.getOrDefault("yPlusCrit", 11.05))
218 {}
219 
220 
222 (
224 )
225 :
227  filmRegionName_(wfpsf.filmRegionName_),
228  B_(wfpsf.B_),
229  yPlusCrit_(wfpsf.yPlusCrit_)
230 {}
231 
232 
234 (
237 )
238 :
240  filmRegionName_(wfpsf.filmRegionName_),
241  B_(wfpsf.B_),
242  yPlusCrit_(wfpsf.yPlusCrit_)
243 {}
244 
245 
246 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
247 
249 {
250  const label patchi = patch().index();
251 
252  const auto& turbModel = db().lookupObject<turbulenceModel>
253  (
255  (
257  internalField().group()
258  )
259  );
260 
261  const scalarField& y = turbModel.y()[patchi];
262 
263  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
264 
265  const tmp<scalarField> tnuw = turbModel.nu(patchi);
266  const scalarField& nuw = tnuw();
267 
268  return y*calcUTau(mag(Uw.snGrad()))/nuw;
269 }
270 
271 
273 {
277 }
278 
279 
280 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281 
282 makePatchTypeField(fvPatchScalarField, nutkFilmWallFunctionFvPatchScalarField);
283 
284 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
285 
286 } // End namespace RASModels
287 } // End namespace compressible
288 } // End namespace Foam
289 
290 // ************************************************************************* //
dictionary dict
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:213
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
dimensionedScalar pow025(const dimensionedScalar &ds)
makePatchTypeField(fvPatchScalarField, alphatFilmWallFunctionFvPatchScalarField)
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
label k
Boltzmann constant.
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
constexpr const char *const group
Group name for atomic constants.
void writeLocalEntries(Ostream &) const
Write local wall function variables.
scalar uTau
scalar y
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
dimensionedScalar exp(const dimensionedScalar &ds)
This boundary condition provides a wall function for the turbulent viscosity (i.e. nut) based on the turbulent kinetic energy, (i.e. k) for for low- and high-Reynolds number applications.
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual tmp< scalarField > calcNut() const
Calculate the turbulence viscosity.
A class for handling words, derived from Foam::string.
Definition: word.H:63
A FieldMapper for finite-volume patch fields.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
bool compressible
Definition: pEqn.H:2
OBJstream os(runTime.globalPath()/outputName)
scalar Cmu() const noexcept
Return the object: Cmu.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
nutkFilmWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
wallFunctionCoefficients wallCoeffs_
Wall-function coefficients.
scalar kappa() const noexcept
Return the object: kappa.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
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.
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127