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-2019 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 "mapDistribute.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 (
66  const fvPatch& p,
68  const fvPatchFieldMapper& mapper
69 )
70 :
71  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
72  B_(ptf.B_),
73  yPlusCrit_(ptf.yPlusCrit_),
74  Cmu_(ptf.Cmu_),
75  kappa_(ptf.kappa_),
76  Prt_(ptf.Prt_)
77 {}
78 
79 
82 (
83  const fvPatch& p,
85  const dictionary& dict
86 )
87 :
88  fixedValueFvPatchScalarField(p, iF, dict),
89  B_(dict.lookupOrDefault("B", 5.5)),
90  yPlusCrit_(dict.lookupOrDefault("yPlusCrit", 11.05)),
91  Cmu_(dict.lookupOrDefault("Cmu", 0.09)),
92  kappa_(dict.lookupOrDefault("kappa", 0.41)),
93  Prt_(dict.lookupOrDefault("Prt", 0.85))
94 {}
95 
96 
99 (
101 )
102 :
103  fixedValueFvPatchScalarField(fwfpsf),
104  B_(fwfpsf.B_),
105  yPlusCrit_(fwfpsf.yPlusCrit_),
106  Cmu_(fwfpsf.Cmu_),
107  kappa_(fwfpsf.kappa_),
108  Prt_(fwfpsf.Prt_)
109 {}
110 
111 
114 (
117 )
118 :
119  fixedValueFvPatchScalarField(fwfpsf, iF),
120  B_(fwfpsf.B_),
121  yPlusCrit_(fwfpsf.yPlusCrit_),
122  Cmu_(fwfpsf.Cmu_),
123  kappa_(fwfpsf.kappa_),
124  Prt_(fwfpsf.Prt_)
125 {}
126 
127 
128 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
129 
131 {
132  if (updated())
133  {
134  return;
135  }
136 
138 
139  // Since we're inside initEvaluate/evaluate there might be processor
140  // comms underway. Change the tag we use.
141  int oldTag = UPstream::msgType();
142  UPstream::msgType() = oldTag+1;
143 
144  bool foundFilm =
145  db().time().foundObject<modelType>("surfaceFilmProperties");
146 
147  if (!foundFilm)
148  {
149  // Do nothing on construction - film model doesn't exist yet
150  return;
151  }
152 
153  const label patchi = patch().index();
154 
155  // Retrieve phase change mass from surface film model
156  const modelType& filmModel =
157  db().time().lookupObject<modelType>("surfaceFilmProperties");
158 
159  const label filmPatchi = filmModel.regionPatchID(patchi);
160 
161  tmp<volScalarField> mDotFilm(filmModel.primaryMassTrans());
162  scalarField mDotFilmp = mDotFilm().boundaryField()[filmPatchi];
163  filmModel.toPrimary(filmPatchi, mDotFilmp);
164 
165  // Retrieve RAS turbulence model
166  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
167  (
169  (
171  internalField().group()
172  )
173  );
174 
175  const scalarField& y = turbModel.y()[patchi];
176  const scalarField& rhow = turbModel.rho().boundaryField()[patchi];
177  const tmp<volScalarField> tk = turbModel.k();
178  const volScalarField& k = tk();
179  const tmp<scalarField> tmuw = turbModel.mu(patchi);
180  const scalarField& muw = tmuw();
181  const tmp<scalarField> talpha = turbModel.alpha(patchi);
182  const scalarField& alphaw = talpha();
183 
184  const scalar Cmu25 = pow(Cmu_, 0.25);
185 
186  // Populate alphat field values
187  scalarField& alphat = *this;
188  forAll(alphat, facei)
189  {
190  label faceCelli = patch().faceCells()[facei];
191 
192  scalar uTau = Cmu25*sqrt(k[faceCelli]);
193 
194  scalar yPlus = y[facei]*uTau/(muw[facei]/rhow[facei]);
195 
196  scalar Pr = muw[facei]/alphaw[facei];
197 
198  scalar factor = 0.0;
199  scalar mStar = mDotFilmp[facei]/(y[facei]*uTau);
200  if (yPlus > yPlusCrit_)
201  {
202  scalar expTerm = exp(min(50.0, yPlusCrit_*mStar*Pr));
203  scalar yPlusRatio = yPlus/yPlusCrit_;
204  scalar powTerm = mStar*Prt_/kappa_;
205  factor =
206  mStar/(expTerm*(pow(yPlusRatio, powTerm)) - 1.0 + rootVSmall);
207  }
208  else
209  {
210  scalar expTerm = exp(min(50.0, yPlus*mStar*Pr));
211  factor = mStar/(expTerm - 1.0 + rootVSmall);
212  }
213 
214  scalar dx = patch().deltaCoeffs()[facei];
215 
216  scalar alphaEff = dx*rhow[facei]*uTau*factor;
217 
218  alphat[facei] = max(alphaEff - alphaw[facei], 0.0);
219  }
220 
221  // Restore tag
222  UPstream::msgType() = oldTag;
223 
224  fixedValueFvPatchScalarField::updateCoeffs();
225 }
226 
227 
228 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
229 
231 {
233  writeEntry(os, "B", B_);
234  writeEntry(os, "yPlusCrit", yPlusCrit_);
235  writeEntry(os, "Cmu", Cmu_);
236  writeEntry(os, "kappa", kappa_);
237  writeEntry(os, "Prt", Prt_);
238  writeEntry(os, "value", *this);
239 }
240 
241 
242 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
243 
245 (
248 );
249 
250 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251 
252 } // End namespace RASModels
253 } // End namespace compressible
254 } // End namespace Foam
255 
256 // ************************************************************************* //
const char *const group
Group name for atomic constants.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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:61
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.
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:280
alphatFilmWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Macros for easy insertion into run-time selection tables.
scalar uTau
scalar y
dimensionedScalar exp(const dimensionedScalar &ds)
static const word propertiesName
Default name of the turbulence properties dictionary.
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
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
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
virtual tmp< volScalarField > alpha() const
Thermal diffusivity for enthalpy of mixture [kg/m/s].
Namespace for OpenFOAM.
This boundary condition provides a turbulent thermal diffusivity condition when using wall functions...