nutkFilmWallFunctionFvPatchScalarField.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 
30 #include "fvPatchFieldMapper.H"
31 #include "volFields.H"
34 #include "surfaceFilmRegionModel.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 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
48 
50 (
51  const scalarField& magGradU
52 ) const
53 {
54  tmp<scalarField> tuTau(new scalarField(patch().size(), Zero));
55  scalarField& uTau = tuTau.ref();
56 
57  const auto* filmModelPtr = db().time().findObject
59  (filmRegionName_);
60 
61  if (!filmModelPtr)
62  {
63  // Do nothing on construction - film model doesn't exist yet
64  return tuTau;
65  }
66 
67  const label patchi = patch().index();
68 
69  // Retrieve phase change mass from surface film model
70  const auto& filmModel = *filmModelPtr;
71 
72  const label filmPatchi = filmModel.regionPatchID(patchi);
73 
74  tmp<volScalarField> mDotFilm(filmModel.primaryMassTrans());
75  scalarField mDotFilmp = mDotFilm().boundaryField()[filmPatchi];
76  filmModel.toPrimary(filmPatchi, mDotFilmp);
77 
78 
79  // Retrieve RAS turbulence model
80  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
81  (
83  (
85  internalField().group()
86  )
87  );
88 
89  const scalarField& y = turbModel.y()[patchi];
90  const tmp<volScalarField> tk = turbModel.k();
91  const volScalarField& k = tk();
92  const tmp<scalarField> tnuw = turbModel.nu(patchi);
93  const scalarField& nuw = tnuw();
94 
95  const scalar Cmu25 = pow(Cmu_, 0.25);
96 
97  forAll(uTau, facei)
98  {
99  label faceCelli = patch().faceCells()[facei];
100 
101  scalar ut = Cmu25*sqrt(k[faceCelli]);
102 
103  scalar yPlus = y[facei]*ut/nuw[facei];
104 
105  scalar mStar = mDotFilmp[facei]/(y[facei]*ut);
106 
107  scalar factor = 0.0;
108  if (yPlus > yPlusCrit_)
109  {
110  scalar expTerm = exp(min(50.0, B_*mStar));
111  scalar powTerm = pow(yPlus, mStar/kappa_);
112  factor = mStar/(expTerm*powTerm - 1.0 + ROOTVSMALL);
113  }
114  else
115  {
116  scalar expTerm = exp(min(50.0, mStar));
117  factor = mStar/(expTerm*yPlus - 1.0 + ROOTVSMALL);
118  }
119 
120  uTau[facei] = sqrt(max(0, magGradU[facei]*ut*factor));
121  }
122 
123  return tuTau;
124 }
125 
126 
128 {
129  const label patchi = patch().index();
130 
131  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
132  (
134  (
136  internalField().group()
137  )
138  );
139 
140  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
141  const scalarField magGradU(mag(Uw.snGrad()));
142  const tmp<scalarField> tnuw = turbModel.nu(patchi);
143  const scalarField& nuw = tnuw();
144 
145  return max
146  (
147  scalar(0),
148  sqr(calcUTau(magGradU))/(magGradU + ROOTVSMALL) - nuw
149  );
150 }
151 
152 
153 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
154 
156 (
157  const fvPatch& p,
159 )
160 :
162  filmRegionName_("surfaceFilmProperties"),
163  B_(5.5),
164  yPlusCrit_(11.05)
165 {}
166 
167 
169 (
171  const fvPatch& p,
173  const fvPatchFieldMapper& mapper
174 )
175 :
176  nutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
177  filmRegionName_(ptf.filmRegionName_),
178  B_(5.5),
179  yPlusCrit_(11.05)
180 {}
181 
182 
184 (
185  const fvPatch& p,
187  const dictionary& dict
188 )
189 :
191  filmRegionName_
192  (
193  dict.getOrDefault<word>("filmRegion", "surfaceFilmProperties")
194  ),
195  B_(dict.getOrDefault("B", 5.5)),
196  yPlusCrit_(dict.getOrDefault("yPlusCrit", 11.05))
197 {}
198 
199 
201 (
203 )
204 :
206  filmRegionName_(wfpsf.filmRegionName_),
207  B_(wfpsf.B_),
208  yPlusCrit_(wfpsf.yPlusCrit_)
209 {}
210 
211 
213 (
216 )
217 :
219  filmRegionName_(wfpsf.filmRegionName_),
220  B_(wfpsf.B_),
221  yPlusCrit_(wfpsf.yPlusCrit_)
222 {}
223 
224 
225 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
226 
228 {
229  const label patchi = patch().index();
230 
231  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
232  (
234  (
236  internalField().group()
237  )
238  );
239 
240  const scalarField& y = turbModel.y()[patchi];
241  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
242  const tmp<scalarField> tnuw = turbModel.nu(patchi);
243  const scalarField& nuw = tnuw();
244 
245  return y*calcUTau(mag(Uw.snGrad()))/nuw;
246 }
247 
248 
250 {
252  os.writeEntryIfDifferent<word>
253  (
254  "filmRegion",
255  "surfaceFilmProperties",
257  );
259  os.writeEntry("B", B_);
260  os.writeEntry("yPlusCrit", yPlusCrit_);
261  writeEntry("value", os);
262 }
263 
264 
265 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
266 
268 
269 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
270 
271 } // End namespace RASModels
272 } // End namespace compressible
273 } // End namespace Foam
274 
275 // ************************************************************************* //
Foam::fvPatchField< vector >
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::fvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:225
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
p
volScalarField & p
Definition: createFieldRefs.H:8
mappedWallPolyPatch.H
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::write
virtual void write(Ostream &os) const
Write.
Definition: nutkFilmWallFunctionFvPatchScalarField.C:249
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::yPlusCrit_
scalar yPlusCrit_
y+ value for laminar -> turbulent transition (default = 11.05)
Definition: nutkFilmWallFunctionFvPatchScalarField.H:87
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
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::yPlus
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
Definition: nutkFilmWallFunctionFvPatchScalarField.C:227
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::calcNut
virtual tmp< scalarField > calcNut() const
Calculate the turbulence viscosity.
Definition: nutkFilmWallFunctionFvPatchScalarField.C:127
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::B_
scalar B_
B Coefficient (default = 5.5)
Definition: nutkFilmWallFunctionFvPatchScalarField.H:84
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
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::nutkFilmWallFunctionFvPatchScalarField
nutkFilmWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: nutkFilmWallFunctionFvPatchScalarField.C:156
fvPatchFieldMapper.H
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
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::tmp::ref
T & ref() const
Definition: tmpI.H:227
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
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::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::calcUTau
virtual tmp< scalarField > calcUTau(const scalarField &magGradU) const
Calculate the friction velocity.
Definition: nutkFilmWallFunctionFvPatchScalarField.C:50
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField::filmRegionName_
word filmRegionName_
Name of film region.
Definition: nutkFilmWallFunctionFvPatchScalarField.H:81
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
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
nutkFilmWallFunctionFvPatchScalarField.H
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
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
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::compressible::RASModels::nutkFilmWallFunctionFvPatchScalarField
This boundary condition provides a turbulent viscosity condition when using wall functions,...
Definition: nutkFilmWallFunctionFvPatchScalarField.H:72
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::nutWallFunctionFvPatchScalarField::writeLocalEntries
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
Definition: nutWallFunctionFvPatchScalarField.C:91
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::nutkWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the turbulent viscosity, i.e....
Definition: nutkWallFunctionFvPatchScalarField.H:94