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 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  int oldTag = UPstream::msgType();
164  UPstream::msgType() = oldTag+1;
165 
166  const label patchi = patch().index();
167 
168  // Retrieve phase change mass from surface film model
169  const label filmPatchi = filmModel.regionPatchID(patchi);
170 
171  tmp<volScalarField> mDotFilm(filmModel.primaryMassTrans());
172  scalarField mDotFilmp = mDotFilm().boundaryField()[filmPatchi];
173  filmModel.toPrimary(filmPatchi, mDotFilmp);
174 
175  // Retrieve RAS turbulence model
176  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
177  (
179  (
181  internalField().group()
182  )
183  );
184 
185  const scalarField& y = turbModel.y()[patchi];
186  const scalarField& rhow = turbModel.rho().boundaryField()[patchi];
187  const tmp<volScalarField> tk = turbModel.k();
188  const volScalarField& k = tk();
189  const tmp<scalarField> tmuw = turbModel.mu(patchi);
190  const scalarField& muw = tmuw();
191  const tmp<scalarField> talpha = turbModel.alpha(patchi);
192  const scalarField& alphaw = talpha();
193 
194  const scalar Cmu25 = pow(Cmu_, 0.25);
195 
196  // Populate alphat field values
197  scalarField& alphat = *this;
198  forAll(alphat, facei)
199  {
200  label faceCelli = patch().faceCells()[facei];
201 
202  scalar uTau = Cmu25*sqrt(k[faceCelli]);
203 
204  scalar yPlus = y[facei]*uTau/(muw[facei]/rhow[facei]);
205 
206  scalar Pr = muw[facei]/alphaw[facei];
207 
208  scalar factor = 0.0;
209  scalar mStar = mDotFilmp[facei]/(y[facei]*uTau);
210  if (yPlus > yPlusCrit_)
211  {
212  scalar expTerm = exp(min(50.0, yPlusCrit_*mStar*Pr));
213  scalar yPlusRatio = yPlus/yPlusCrit_;
214  scalar powTerm = mStar*Prt_/kappa_;
215  factor =
216  mStar/(expTerm*(pow(yPlusRatio, powTerm)) - 1.0 + ROOTVSMALL);
217  }
218  else
219  {
220  scalar expTerm = exp(min(50.0, yPlus*mStar*Pr));
221  factor = mStar/(expTerm - 1.0 + ROOTVSMALL);
222  }
223 
224  scalar dx = patch().deltaCoeffs()[facei];
225 
226  scalar alphaEff = dx*rhow[facei]*uTau*factor;
227 
228  alphat[facei] = max(alphaEff - alphaw[facei], 0.0);
229  }
230 
231  // Restore tag
232  UPstream::msgType() = oldTag;
233 
234  fixedValueFvPatchScalarField::updateCoeffs();
235 }
236 
237 
238 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
239 
241 {
243  os.writeEntryIfDifferent<word>
244  (
245  "filmRegion",
246  "surfaceFilmProperties",
248  );
249  os.writeEntry("B", B_);
250  os.writeEntry("yPlusCrit", yPlusCrit_);
251  os.writeEntry("Cmu", Cmu_);
252  os.writeEntry("kappa", kappa_);
253  os.writeEntry("Prt", Prt_);
254  writeEntry("value", os);
255 }
256 
257 
258 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
259 
261 (
264 );
265 
266 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
267 
268 } // End namespace RASModels
269 } // End namespace compressible
270 } // End namespace Foam
271 
272 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
p
volScalarField & p
Definition: createFieldRefs.H:8
mappedWallPolyPatch.H
Foam::compressible::RASModels::alphatFilmWallFunctionFvPatchScalarField::yPlusCrit_
scalar yPlusCrit_
y+ value for laminar -> turbulent transition (default = 11.05)
Definition: alphatFilmWallFunctionFvPatchScalarField.H:132
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::constant::atomic::group
constexpr const char *const group
Group name for atomic constants.
Definition: atomicConstants.H:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
alphatFilmWallFunctionFvPatchScalarField.H
Foam::compressible::RASModels::makePatchTypeField
makePatchTypeField(fvPatchScalarField, alphatFilmWallFunctionFvPatchScalarField)
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
fvPatchFieldMapper.H
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::compressible::RASModels::alphatFilmWallFunctionFvPatchScalarField::kappa_
scalar kappa_
Von-Karman constant (default = 0.41)
Definition: alphatFilmWallFunctionFvPatchScalarField.H:138
surfaceFilmRegionModel.H
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::compressible::RASModels::alphatFilmWallFunctionFvPatchScalarField::alphatFilmWallFunctionFvPatchScalarField
alphatFilmWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: alphatFilmWallFunctionFvPatchScalarField.C:51
Foam::Field< scalar >
Foam::regionModels::surfaceFilmModels::surfaceFilmRegionModel
Base class for surface film models.
Definition: surfaceFilmRegionModel.H:55
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Pr
dimensionedScalar Pr("Pr", dimless, laminarTransport)
Foam::compressible::RASModels::alphatFilmWallFunctionFvPatchScalarField
This boundary condition provides a turbulent thermal diffusivity condition when using wall functions,...
Definition: alphatFilmWallFunctionFvPatchScalarField.H:117
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
uTau
scalar uTau
Definition: evaluateNearWall.H:14
compressible
bool compressible
Definition: pEqn.H:2
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::compressible::RASModels::alphatFilmWallFunctionFvPatchScalarField::Prt_
scalar Prt_
Turbulent Prandtl number (default = 0.85)
Definition: alphatFilmWallFunctionFvPatchScalarField.H:141
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::compressible::RASModels::alphatFilmWallFunctionFvPatchScalarField::B_
scalar B_
B Coefficient (default = 5.5)
Definition: alphatFilmWallFunctionFvPatchScalarField.H:129
Foam::compressible::RASModels::alphatFilmWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: alphatFilmWallFunctionFvPatchScalarField.C:240
Foam::UPstream::msgType
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:540
Foam::compressible::RASModels::alphatFilmWallFunctionFvPatchScalarField::filmRegionName_
word filmRegionName_
Name of film region.
Definition: alphatFilmWallFunctionFvPatchScalarField.H:126
mapDistribute.H
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::ThermalDiffusivity
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Definition: phaseCompressibleTurbulenceModelFwd.H:47
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
alphaEff
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
yPlus
scalar yPlus
Definition: evaluateNearWall.H:16
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::GeometricField< scalar, fvPatchField, volMesh >
turbulentFluidThermoModel.H
Foam::ThermalDiffusivity::alpha
virtual tmp< volScalarField > alpha() const
Return the laminar thermal diffusivity for enthalpy [kg/m/s].
Definition: ThermalDiffusivity.H:125
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::compressible::RASModels::alphatFilmWallFunctionFvPatchScalarField::Cmu_
scalar Cmu_
Turbulent Cmu coefficient (default = 0.09)
Definition: alphatFilmWallFunctionFvPatchScalarField.H:135
Foam::compressible::RASModels::alphatFilmWallFunctionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: alphatFilmWallFunctionFvPatchScalarField.C:141