alphatFilmWallFunctionFvPatchScalarField.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 
28 #include "surfaceFilmRegionModel.H"
29 #include "fvPatchFieldMapper.H"
30 #include "volFields.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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
48 (
49  const fvPatch& p,
51 )
52 :
53  fixedValueFvPatchScalarField(p, iF),
54  B_(5.5),
55  yPlusCrit_(11.05),
56  Cmu_(0.09),
57  kappa_(0.41),
58  Prt_(0.85)
59 {}
60 
61 
64 (
65  const fvPatch& p,
67  const dictionary& dict
68 )
69 :
70  fixedValueFvPatchScalarField(p, iF, dict),
71  B_(dict.lookupOrDefault("B", 5.5)),
72  yPlusCrit_(dict.lookupOrDefault("yPlusCrit", 11.05)),
73  Cmu_(dict.lookupOrDefault("Cmu", 0.09)),
74  kappa_(dict.lookupOrDefault("kappa", 0.41)),
75  Prt_(dict.lookupOrDefault("Prt", 0.85))
76 {}
77 
78 
81 (
83  const fvPatch& p,
85  const fvPatchFieldMapper& mapper
86 )
87 :
88  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
89  B_(ptf.B_),
90  yPlusCrit_(ptf.yPlusCrit_),
91  Cmu_(ptf.Cmu_),
92  kappa_(ptf.kappa_),
93  Prt_(ptf.Prt_)
94 {}
95 
96 
99 (
102 )
103 :
104  fixedValueFvPatchScalarField(fwfpsf, iF),
105  B_(fwfpsf.B_),
106  yPlusCrit_(fwfpsf.yPlusCrit_),
107  Cmu_(fwfpsf.Cmu_),
108  kappa_(fwfpsf.kappa_),
109  Prt_(fwfpsf.Prt_)
110 {}
111 
112 
113 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
114 
116 {
117  if (updated())
118  {
119  return;
120  }
121 
123 
124  // Since we're inside initEvaluate/evaluate there might be processor
125  // comms underway. Change the tag we use.
126  int oldTag = UPstream::msgType();
127  UPstream::msgType() = oldTag+1;
128 
129  bool foundFilm =
130  db().time().foundObject<modelType>("surfaceFilmProperties");
131 
132  if (!foundFilm)
133  {
134  // Do nothing on construction - film model doesn't exist yet
135  return;
136  }
137 
138  const label patchi = patch().index();
139 
140  // Retrieve phase change mass from surface film model
141  const modelType& filmModel =
142  db().time().lookupObject<modelType>("surfaceFilmProperties");
143 
144  const label filmPatchi = filmModel.regionPatchID(patchi);
145 
146  tmp<volScalarField> mDotFilm(filmModel.primaryMassTrans());
147  scalarField mDotFilmp = mDotFilm().boundaryField()[filmPatchi];
148  filmModel.toPrimary(filmPatchi, mDotFilmp);
149 
150  const thermophysicalTransportModel& ttm =
151  db().lookupObject<thermophysicalTransportModel>
152  (
154  (
155  thermophysicalTransportModel::typeName,
156  internalField().group()
157  )
158  );
159 
160  const compressibleMomentumTransportModel& turbModel =
161  ttm.momentumTransport();
162 
163  const scalarField& y = turbModel.y()[patchi];
164  const scalarField& rhow = turbModel.rho().boundaryField()[patchi];
165  const tmp<volScalarField> tk = turbModel.k();
166  const volScalarField& k = tk();
167  const tmp<scalarField> tnuw = turbModel.nu(patchi);
168  const scalarField& nuw = tnuw();
169 
170  const tmp<scalarField> talphaw
171  (
172  ttm.thermo().kappa().boundaryField()[patchi]
173  /ttm.thermo().Cp().boundaryField()[patchi]
174  );
175  const scalarField& alphaw = talphaw();
176 
177  const scalar Cmu25 = pow(Cmu_, 0.25);
178 
179  // Populate alphat field values
180  scalarField& alphat = *this;
181  forAll(alphat, facei)
182  {
183  label faceCelli = patch().faceCells()[facei];
184 
185  scalar uTau = Cmu25*sqrt(k[faceCelli]);
186 
187  scalar yPlus = y[facei]*uTau/nuw[facei];
188 
189  scalar Pr = rhow[facei]*nuw[facei]/alphaw[facei];
190 
191  scalar factor = 0.0;
192  scalar mStar = mDotFilmp[facei]/(y[facei]*uTau);
193  if (yPlus > yPlusCrit_)
194  {
195  scalar expTerm = exp(min(50.0, yPlusCrit_*mStar*Pr));
196  scalar yPlusRatio = yPlus/yPlusCrit_;
197  scalar powTerm = mStar*Prt_/kappa_;
198  factor =
199  mStar/(expTerm*(pow(yPlusRatio, powTerm)) - 1.0 + rootVSmall);
200  }
201  else
202  {
203  scalar expTerm = exp(min(50.0, yPlus*mStar*Pr));
204  factor = mStar/(expTerm - 1.0 + rootVSmall);
205  }
206 
207  scalar dx = patch().deltaCoeffs()[facei];
208 
209  scalar alphaEff = dx*rhow[facei]*uTau*factor;
210 
211  alphat[facei] = max(alphaEff - alphaw[facei], 0.0);
212  }
213 
214  // Restore tag
215  UPstream::msgType() = oldTag;
216 
217  fixedValueFvPatchScalarField::updateCoeffs();
218 }
219 
220 
221 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
222 
224 {
226  writeEntry(os, "B", B_);
227  writeEntry(os, "yPlusCrit", yPlusCrit_);
228  writeEntry(os, "Cmu", Cmu_);
229  writeEntry(os, "kappa", kappa_);
230  writeEntry(os, "Prt", Prt_);
231  writeEntry(os, "value", *this);
232 }
233 
234 
235 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236 
238 (
241 );
242 
243 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244 
245 } // End namespace RASModels
246 } // End namespace compressible
247 } // End namespace Foam
248 
249 // ************************************************************************* //
const char *const group
Group name for atomic constants.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual const fluidThermo & thermo() const =0
Access function to incompressible transport model.
const compressibleMomentumTransportModel & momentumTransport() const
const Boundary & boundaryField() const
Return const-reference to the boundary field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
scalar yPlusCrit_
y+ value for laminar -> turbulent transition (default = 11.05)
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
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:476
label k
Boltzmann constant.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:243
alphatFilmWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Macros for easy insertion into run-time selection tables.
virtual const volScalarField & Cp() const =0
Heat capacity at constant pressure [J/kg/K].
scalar uTau
scalar y
dimensionedScalar exp(const dimensionedScalar &ds)
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
virtual const volScalarField & kappa() const =0
Thermal conductivity of mixture [W/m/K].
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
Abstract base class for thermophysical transport models (RAS, LES and laminar).
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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...
A class for managing temporary objects.
Definition: PtrList.H:53
Base class for single-phase compressible turbulence models.
Namespace for OpenFOAM.
This boundary condition provides a turbulent thermal diffusivity condition when using wall functions...