1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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.
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
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "alphatFilmWallFunctionFvPatchScalarField.H"
28 #include "surfaceFilmModel.H"
29 #include "fvPatchFieldMapper.H"
30 #include "volFields.H"
31 #include "addToRunTimeSelectionTable.H"
32 #include "directMappedWallPolyPatch.H"
33 #include "mapDistribute.H"
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 namespace compressible
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 alphatFilmWallFunctionFvPatchScalarField::
47 alphatFilmWallFunctionFvPatchScalarField
50 const DimensionedField<scalar, volMesh>& iF
53 fixedValueFvPatchScalarField(p, iF),
62 alphatFilmWallFunctionFvPatchScalarField::
63 alphatFilmWallFunctionFvPatchScalarField
65 const alphatFilmWallFunctionFvPatchScalarField& ptf,
67 const DimensionedField<scalar, volMesh>& iF,
68 const fvPatchFieldMapper& mapper
71 fixedValueFvPatchScalarField(ptf, p, iF, mapper),
73 yPlusCrit_(ptf.yPlusCrit_),
80 alphatFilmWallFunctionFvPatchScalarField::
81 alphatFilmWallFunctionFvPatchScalarField
84 const DimensionedField<scalar, volMesh>& iF,
85 const dictionary& dict
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))
97 alphatFilmWallFunctionFvPatchScalarField::
98 alphatFilmWallFunctionFvPatchScalarField
100 const alphatFilmWallFunctionFvPatchScalarField& fwfpsf
103 fixedValueFvPatchScalarField(fwfpsf),
105 yPlusCrit_(fwfpsf.yPlusCrit_),
107 kappa_(fwfpsf.kappa_),
112 alphatFilmWallFunctionFvPatchScalarField::
113 alphatFilmWallFunctionFvPatchScalarField
115 const alphatFilmWallFunctionFvPatchScalarField& fwfpsf,
116 const DimensionedField<scalar, volMesh>& iF
119 fixedValueFvPatchScalarField(fwfpsf, iF),
121 yPlusCrit_(fwfpsf.yPlusCrit_),
123 kappa_(fwfpsf.kappa_),
128 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
130 void alphatFilmWallFunctionFvPatchScalarField::updateCoeffs()
137 typedef regionModels::surfaceFilmModels::surfaceFilmModel modelType;
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;
145 db().objectRegistry::foundObject<modelType>("surfaceFilmProperties");
149 // do nothing on construction - film model doesn't exist yet
153 const label patchI = patch().index();
155 // Retrieve phase change mass from surface film model
156 const modelType& filmModel =
157 db().objectRegistry::lookupObject<modelType>("surfaceFilmProperties");
159 const label filmPatchI = filmModel.regionPatchID(patchI);
161 const mapDistribute& distMap = filmModel.mappedPatches()[filmPatchI].map();
163 tmp<volScalarField> mDotFilm(filmModel.primaryMassTrans());
164 scalarField mDotFilmp = mDotFilm().boundaryField()[filmPatchI];
165 distMap.distribute(mDotFilmp);
167 // Retrieve RAS turbulence model
168 const RASModel& rasModel = db().lookupObject<RASModel>("RASProperties");
170 const scalarField& y = rasModel.y()[patchI];
171 const scalarField& rhow = rasModel.rho().boundaryField()[patchI];
172 const tmp<volScalarField> tk = rasModel.k();
173 const volScalarField& k = tk();
174 const scalarField& muw = rasModel.mu().boundaryField()[patchI];
175 const scalarField& alphaw = rasModel.alpha().boundaryField()[patchI];
177 const scalar Cmu25 = pow(Cmu_, 0.25);
179 // Populate alphat field values
180 scalarField& alphat = *this;
181 forAll(alphat, faceI)
183 label faceCellI = patch().faceCells()[faceI];
185 scalar uTau = Cmu25*sqrt(k[faceCellI]);
187 scalar yPlus = y[faceI]*uTau/(muw[faceI]/rhow[faceI]);
189 scalar Pr = muw[faceI]/alphaw[faceI];
192 scalar mStar = mDotFilmp[faceI]/(y[faceI]*uTau);
193 if (yPlus > yPlusCrit_)
195 scalar expTerm = exp(min(50.0, yPlusCrit_*mStar*Pr));
196 scalar yPlusRatio = yPlus/yPlusCrit_;
197 scalar powTerm = mStar*Prt_/kappa_;
199 mStar/(expTerm*(pow(yPlusRatio, powTerm)) - 1.0 + ROOTVSMALL);
203 scalar expTerm = exp(min(50.0, yPlus*mStar*Pr));
204 factor = mStar/(expTerm - 1.0 + ROOTVSMALL);
207 scalar dx = patch().deltaCoeffs()[faceI];
209 scalar alphaEff = dx*rhow[faceI]*uTau*factor;
211 alphat[faceI] = max(alphaEff - alphaw[faceI], 0.0);
215 UPstream::msgType() = oldTag;
217 fixedValueFvPatchScalarField::updateCoeffs();
222 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
224 void alphatFilmWallFunctionFvPatchScalarField::write(Ostream& os) const
226 fvPatchField<scalar>::write(os);
227 os.writeKeyword("B") << B_ << token::END_STATEMENT << nl;
228 os.writeKeyword("yPlusCrit") << yPlusCrit_ << token::END_STATEMENT << nl;
229 os.writeKeyword("Cmu") << Cmu_ << token::END_STATEMENT << nl;
230 os.writeKeyword("kappa") << kappa_ << token::END_STATEMENT << nl;
231 os.writeKeyword("Prt") << Prt_ << token::END_STATEMENT << nl;
232 writeEntry("value", os);
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241 alphatFilmWallFunctionFvPatchScalarField
244 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
246 } // End namespace RASModels
247 } // End namespace compressible
248 } // End namespace Foam
250 // ************************************************************************* //