alphatFilmWallFunctionFvPatchScalarField.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 
31 #include "surfaceFilmRegionModel.H"
32 #include "fvPatchFieldMapper.H"
33 #include "volFields.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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
51 (
52  const fvPatch& p,
54 )
55 :
56  fixedValueFvPatchScalarField(p, iF),
57  filmRegionName_("surfaceFilmProperties"),
58  B_(5.5),
59  yPlusCrit_(11.05),
60  Cmu_(0.09),
61  kappa_(0.41),
62  Prt_(0.85)
63 {}
64 
65 
68 (
70  const fvPatch& p,
72  const fvPatchFieldMapper& mapper
73 )
74 :
75  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
76  filmRegionName_(ptf.filmRegionName_),
77  B_(ptf.B_),
78  yPlusCrit_(ptf.yPlusCrit_),
79  Cmu_(ptf.Cmu_),
80  kappa_(ptf.kappa_),
81  Prt_(ptf.Prt_)
82 {}
83 
84 
87 (
88  const fvPatch& p,
90  const dictionary& dict
91 )
92 :
93  fixedValueFvPatchScalarField(p, iF, dict),
94  filmRegionName_
95  (
96  dict.getOrDefault<word>("filmRegion", "surfaceFilmProperties")
97  ),
98  B_(dict.getOrDefault("B", 5.5)),
99  yPlusCrit_(dict.getOrDefault("yPlusCrit", 11.05)),
100  Cmu_(dict.getOrDefault("Cmu", 0.09)),
101  kappa_(dict.getOrDefault("kappa", 0.41)),
102  Prt_(dict.getOrDefault("Prt", 0.85))
103 {}
104 
105 
108 (
110 )
111 :
112  fixedValueFvPatchScalarField(fwfpsf),
113  filmRegionName_(fwfpsf.filmRegionName_),
114  B_(fwfpsf.B_),
115  yPlusCrit_(fwfpsf.yPlusCrit_),
116  Cmu_(fwfpsf.Cmu_),
117  kappa_(fwfpsf.kappa_),
118  Prt_(fwfpsf.Prt_)
119 {}
120 
121 
124 (
127 )
128 :
129  fixedValueFvPatchScalarField(fwfpsf, iF),
130  filmRegionName_(fwfpsf.filmRegionName_),
131  B_(fwfpsf.B_),
132  yPlusCrit_(fwfpsf.yPlusCrit_),
133  Cmu_(fwfpsf.Cmu_),
134  kappa_(fwfpsf.kappa_),
135  Prt_(fwfpsf.Prt_)
136 {}
137 
138 
139 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
140 
142 {
143  if (updated())
144  {
145  return;
146  }
147 
148  const auto* filmModelPtr = db().time().findObject
150  (filmRegionName_);
151 
152  if (!filmModelPtr)
153  {
154  // Do nothing on construction - film model doesn't exist yet
155  return;
156  }
157 
158  const auto& filmModel = *filmModelPtr;
159 
160 
161  // Since we're inside initEvaluate/evaluate there might be processor
162  // comms underway. Change the tag we use.
163  const int oldTag = UPstream::incrMsgType();
164 
165  const label patchi = patch().index();
166 
167  // Retrieve phase change mass from surface film model
168  const label filmPatchi = filmModel.regionPatchID(patchi);
169 
170  tmp<volScalarField> mDotFilm = filmModel.primaryMassTrans();
171  scalarField mDotFilmp = mDotFilm().boundaryField()[filmPatchi];
172  filmModel.toPrimary(filmPatchi, mDotFilmp);
173 
174  // Retrieve RAS turbulence model
175  const auto& turbModel = db().lookupObject<turbulenceModel>
176  (
178  (
180  internalField().group()
181  )
182  );
183 
184  const scalarField& y = turbModel.y()[patchi];
185 
186  const scalarField& rhow = turbModel.rho().boundaryField()[patchi];
187 
188  const tmp<volScalarField> tk = turbModel.k();
189  const volScalarField& k = tk();
190 
191  const tmp<scalarField> tmuw = turbModel.mu(patchi);
192  const scalarField& muw = tmuw();
193 
194  const tmp<scalarField> talpha = turbModel.alpha(patchi);
195  const scalarField& alphaw = talpha();
196 
197  const scalar Cmu25 = pow025(Cmu_);
198 
199  // Populate alphat field values
200  scalarField& alphat = *this;
201  forAll(alphat, facei)
202  {
203  const label faceCelli = patch().faceCells()[facei];
204 
205  const scalar uTau = Cmu25*sqrt(k[faceCelli]);
206 
207  const scalar yPlus = y[facei]*uTau/(muw[facei]/rhow[facei]);
208 
209  const scalar Pr = muw[facei]/alphaw[facei];
210 
211  scalar factor = 0;
212  const scalar mStar = mDotFilmp[facei]/(y[facei]*uTau);
213  if (yPlus > yPlusCrit_)
214  {
215  const scalar expTerm = exp(min(scalar(50), yPlusCrit_*mStar*Pr));
216  const scalar yPlusRatio = yPlus/yPlusCrit_;
217  const scalar powTerm = mStar*Prt_/kappa_;
218 
219  factor =
220  mStar/(expTerm*(pow(yPlusRatio, powTerm)) - 1.0 + ROOTVSMALL);
221  }
222  else
223  {
224  const scalar expTerm = exp(min(scalar(50), yPlus*mStar*Pr));
225 
226  factor = mStar/(expTerm - 1.0 + ROOTVSMALL);
227  }
228 
229  const scalar dx = patch().deltaCoeffs()[facei];
230 
231  const scalar alphaEff = dx*rhow[facei]*uTau*factor;
232 
233  alphat[facei] = max(alphaEff - alphaw[facei], scalar(0));
234  }
235 
236  UPstream::msgType(oldTag); // Restore tag
238  fixedValueFvPatchScalarField::updateCoeffs();
239 }
240 
241 
242 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
243 
245 {
247  os.writeEntryIfDifferent<word>
248  (
249  "filmRegion",
250  "surfaceFilmProperties",
252  );
253  os.writeEntryIfDifferent<scalar>("B", 5.5, B_);
254  os.writeEntryIfDifferent<scalar>("yPlusCrit", 11.05, yPlusCrit_);
255  os.writeEntryIfDifferent<scalar>("Cmu", 0.09, Cmu_);
256  os.writeEntryIfDifferent<scalar>("kappa", 0.41, kappa_);
257  os.writeEntryIfDifferent<scalar>("Prt", 0.85, Prt_);
258 
260 }
261 
262 
263 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
264 
266 (
268  alphatFilmWallFunctionFvPatchScalarField
269 );
270 
271 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
272 
273 } // End namespace RASModels
274 } // End namespace compressible
275 } // End namespace Foam
276 
277 // ************************************************************************* //
dimensionedScalar Pr("Pr", dimless, laminarTransport)
dictionary dict
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition: UPstream.H:1251
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
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)
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:1229
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
label k
Boltzmann constant.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
Abstract base class for turbulence models (RAS, LES and laminar).
alphatFilmWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
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.
scalar uTau
scalar y
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
dimensionedScalar exp(const dimensionedScalar &ds)
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
A FieldMapper for finite-volume patch fields.
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)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
scalar yPlus
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.
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Namespace for OpenFOAM.
This boundary condition provides a turbulent thermal diffusivity condition when using wall functions...