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 | Copyright (C) 2011-2015 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 "surfaceFilmModel.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  filmRegionName_("surfaceFilmProperties"),
55  B_(5.5),
56  yPlusCrit_(11.05),
57  Cmu_(0.09),
58  kappa_(0.41),
59  Prt_(0.85)
60 {}
61 
62 
65 (
67  const fvPatch& p,
69  const fvPatchFieldMapper& mapper
70 )
71 :
72  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
73  filmRegionName_(ptf.filmRegionName_),
74  B_(ptf.B_),
75  yPlusCrit_(ptf.yPlusCrit_),
76  Cmu_(ptf.Cmu_),
77  kappa_(ptf.kappa_),
78  Prt_(ptf.Prt_)
79 {}
80 
81 
84 (
85  const fvPatch& p,
87  const dictionary& dict
88 )
89 :
90  fixedValueFvPatchScalarField(p, iF, dict),
91  filmRegionName_
92  (
93  dict.lookupOrDefault<word>("filmRegion", "surfaceFilmProperties")
94  ),
95  B_(dict.lookupOrDefault("B", 5.5)),
96  yPlusCrit_(dict.lookupOrDefault("yPlusCrit", 11.05)),
97  Cmu_(dict.lookupOrDefault("Cmu", 0.09)),
98  kappa_(dict.lookupOrDefault("kappa", 0.41)),
99  Prt_(dict.lookupOrDefault("Prt", 0.85))
100 {}
101 
102 
105 (
107 )
108 :
109  fixedValueFvPatchScalarField(fwfpsf),
110  filmRegionName_(fwfpsf.filmRegionName_),
111  B_(fwfpsf.B_),
112  yPlusCrit_(fwfpsf.yPlusCrit_),
113  Cmu_(fwfpsf.Cmu_),
114  kappa_(fwfpsf.kappa_),
115  Prt_(fwfpsf.Prt_)
116 {}
117 
118 
121 (
124 )
125 :
126  fixedValueFvPatchScalarField(fwfpsf, iF),
127  filmRegionName_(fwfpsf.filmRegionName_),
128  B_(fwfpsf.B_),
129  yPlusCrit_(fwfpsf.yPlusCrit_),
130  Cmu_(fwfpsf.Cmu_),
131  kappa_(fwfpsf.kappa_),
132  Prt_(fwfpsf.Prt_)
133 {}
134 
135 
136 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
137 
139 {
140  if (updated())
141  {
142  return;
143  }
144 
146 
147  // Since we're inside initEvaluate/evaluate there might be processor
148  // comms underway. Change the tag we use.
149  int oldTag = UPstream::msgType();
150  UPstream::msgType() = oldTag+1;
151 
152  bool foundFilm = db().time().foundObject<modelType>(filmRegionName_);
153 
154  if (!foundFilm)
155  {
156  // do nothing on construction - film model doesn't exist yet
157  return;
158  }
159 
160  const label patchi = patch().index();
161 
162  // Retrieve phase change mass from surface film model
163  const modelType& filmModel =
164  db().time().lookupObject<modelType>(filmRegionName_);
165 
166  const label filmPatchI = filmModel.regionPatchID(patchi);
167 
168  tmp<volScalarField> mDotFilm(filmModel.primaryMassTrans());
169  scalarField mDotFilmp = mDotFilm().boundaryField()[filmPatchI];
170  filmModel.toPrimary(filmPatchI, mDotFilmp);
171 
172  // Retrieve RAS turbulence model
173  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
174  (
176  (
179  )
180  );
181 
182  const scalarField& y = turbModel.y()[patchi];
183  const scalarField& rhow = turbModel.rho().boundaryField()[patchi];
184  const tmp<volScalarField> tk = turbModel.k();
185  const volScalarField& k = tk();
186  const tmp<scalarField> tmuw = turbModel.mu(patchi);
187  const scalarField& muw = tmuw();
188  const tmp<scalarField> talpha = turbModel.alpha(patchi);
189  const scalarField& alphaw = talpha();
190 
191  const scalar Cmu25 = pow(Cmu_, 0.25);
192 
193  // Populate alphat field values
194  scalarField& alphat = *this;
195  forAll(alphat, faceI)
196  {
197  label faceCellI = patch().faceCells()[faceI];
198 
199  scalar uTau = Cmu25*sqrt(k[faceCellI]);
200 
201  scalar yPlus = y[faceI]*uTau/(muw[faceI]/rhow[faceI]);
202 
203  scalar Pr = muw[faceI]/alphaw[faceI];
204 
205  scalar factor = 0.0;
206  scalar mStar = mDotFilmp[faceI]/(y[faceI]*uTau);
207  if (yPlus > yPlusCrit_)
208  {
209  scalar expTerm = exp(min(50.0, yPlusCrit_*mStar*Pr));
210  scalar yPlusRatio = yPlus/yPlusCrit_;
211  scalar powTerm = mStar*Prt_/kappa_;
212  factor =
213  mStar/(expTerm*(pow(yPlusRatio, powTerm)) - 1.0 + ROOTVSMALL);
214  }
215  else
216  {
217  scalar expTerm = exp(min(50.0, yPlus*mStar*Pr));
218  factor = mStar/(expTerm - 1.0 + ROOTVSMALL);
219  }
220 
221  scalar dx = patch().deltaCoeffs()[faceI];
222 
223  scalar alphaEff = dx*rhow[faceI]*uTau*factor;
224 
225  alphat[faceI] = max(alphaEff - alphaw[faceI], 0.0);
226  }
227 
228  // Restore tag
229  UPstream::msgType() = oldTag;
230 
231  fixedValueFvPatchScalarField::updateCoeffs();
232 }
233 
234 
235 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
236 
238 {
240  writeEntryIfDifferent<word>
241  (
242  os,
243  "filmRegion",
244  "surfaceFilmProperties",
246  );
247  os.writeKeyword("B") << B_ << token::END_STATEMENT << nl;
248  os.writeKeyword("yPlusCrit") << yPlusCrit_ << token::END_STATEMENT << nl;
249  os.writeKeyword("Cmu") << Cmu_ << token::END_STATEMENT << nl;
250  os.writeKeyword("kappa") << kappa_ << token::END_STATEMENT << nl;
251  os.writeKeyword("Prt") << Prt_ << token::END_STATEMENT << nl;
252  writeEntry("value", os);
253 }
254 
255 
256 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257 
259 (
262 );
263 
264 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
265 
266 } // End namespace RASModels
267 } // End namespace compressible
268 } // End namespace Foam
269 
270 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:349
Foam::token::END_STATEMENT
@ END_STATEMENT
Definition: token.H:99
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
mappedWallPolyPatch.H
Foam::compressible::RASModels::alphatFilmWallFunctionFvPatchScalarField::yPlusCrit_
scalar yPlusCrit_
y+ value for laminar -> turbulent transition (default = 11.05)
Definition: alphatFilmWallFunctionFvPatchScalarField.H:131
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
dimensionedInternalField
rDeltaT dimensionedInternalField()
alphat
alphat
Definition: TEqn.H:2
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:97
Foam::constant::atomic::group
const char *const group
Group name for atomic constants.
Definition: atomicConstants.C:40
fvPatchFieldMapper.H
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
Foam::compressible::RASModels::alphatFilmWallFunctionFvPatchScalarField::kappa_
scalar kappa_
Von-Karman constant (default = 0.41)
Definition: alphatFilmWallFunctionFvPatchScalarField.H:137
Foam::compressible::RASModels::alphatFilmWallFunctionFvPatchScalarField::alphatFilmWallFunctionFvPatchScalarField
alphatFilmWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: alphatFilmWallFunctionFvPatchScalarField.C:48
Foam::label
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
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::compressible::RASModels::alphatFilmWallFunctionFvPatchScalarField
This boundary condition provides a turbulent thermal diffusivity condition when using wall functions,...
Definition: alphatFilmWallFunctionFvPatchScalarField.H:116
surfaceFilmModel.H
Pr
dimensionedScalar Pr("Pr", dimless, laminarTransport)
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
uTau
scalar uTau
Definition: evaluateNearWall.H:14
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
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:140
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
alphaEff
const volScalarField & alphaEff
Definition: setAlphaEff.H:93
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::compressible::RASModels::alphatFilmWallFunctionFvPatchScalarField::B_
scalar B_
B Coefficient (default = 5.5)
Definition: alphatFilmWallFunctionFvPatchScalarField.H:128
Foam::UPstream::msgType
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:452
Foam::compressible::RASModels::alphatFilmWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: alphatFilmWallFunctionFvPatchScalarField.C:237
Foam::compressible::RASModels::alphatFilmWallFunctionFvPatchScalarField::filmRegionName_
word filmRegionName_
Name of film region.
Definition: alphatFilmWallFunctionFvPatchScalarField.H:125
mapDistribute.H
Foam::ThermalDiffusivity
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Definition: ThermalDiffusivity.H:48
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
compressible
bool compressible
Definition: pEqn.H:40
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
Foam::yPlus
This function object evaluates and outputs turbulence y+ for turbulence models. The field is stored o...
Definition: yPlus.H:110
patchi
label patchi
Definition: getPatchFieldScalar.H:1
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::regionModels::surfaceFilmModels::surfaceFilmModel
Base class for surface film models.
Definition: surfaceFilmModel.H:61
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
turbulentFluidThermoModel.H
Foam::ThermalDiffusivity::alpha
virtual tmp< volScalarField > alpha() const
Return the laminar thermal diffusivity for enthalpy [kg/m/s].
Definition: ThermalDiffusivity.H:123
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::compressible::RASModels::alphatFilmWallFunctionFvPatchScalarField::Cmu_
scalar Cmu_
Turbulent Cmu coefficient (default = 0.09)
Definition: alphatFilmWallFunctionFvPatchScalarField.H:134
Foam::compressible::RASModels::alphatFilmWallFunctionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: alphatFilmWallFunctionFvPatchScalarField.C:138