00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include "kappatJayatillekeWallFunctionFvPatchScalarField.H"
00027 #include <incompressibleRASModels/RASModel.H>
00028 #include <finiteVolume/fvPatchFieldMapper.H>
00029 #include <finiteVolume/volFields.H>
00030 #include <finiteVolume/wallFvPatch.H>
00031 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00032
00033
00034
00035 namespace Foam
00036 {
00037 namespace incompressible
00038 {
00039 namespace RASModels
00040 {
00041
00042
00043
00044 scalar kappatJayatillekeWallFunctionFvPatchScalarField::tolerance_ = 0.01;
00045 label kappatJayatillekeWallFunctionFvPatchScalarField::maxIters_ = 10;
00046
00047
00048
00049 void kappatJayatillekeWallFunctionFvPatchScalarField::checkType()
00050 {
00051 if (!isA<wallFvPatch>(patch()))
00052 {
00053 FatalErrorIn
00054 (
00055 "kappatJayatillekeWallFunctionFvPatchScalarField::checkType()"
00056 ) << "Invalid wall function specification" << nl
00057 << " Patch type for patch " << patch().name()
00058 << " must be wall" << nl
00059 << " Current patch type is " << patch().type() << nl << endl
00060 << abort(FatalError);
00061 }
00062 }
00063
00064
00065 scalar kappatJayatillekeWallFunctionFvPatchScalarField::Psmooth
00066 (
00067 const scalar Prat
00068 ) const
00069 {
00070 return 9.24*(pow(Prat, 0.75) - 1.0)*(1.0 + 0.28*exp(-0.007*Prat));
00071 }
00072
00073
00074 scalar kappatJayatillekeWallFunctionFvPatchScalarField::yPlusTherm
00075 (
00076 const scalar P,
00077 const scalar Prat
00078 ) const
00079 {
00080 scalar ypt = 11.0;
00081
00082 for (int i=0; i<maxIters_; i++)
00083 {
00084 scalar f = ypt - (log(E_*ypt)/kappa_ + P)/Prat;
00085 scalar df = 1.0 - 1.0/(ypt*kappa_*Prat);
00086 scalar yptNew = ypt - f/df;
00087
00088 if (yptNew < VSMALL)
00089 {
00090 return 0;
00091 }
00092 else if (mag(yptNew - ypt) < tolerance_)
00093 {
00094 return yptNew;
00095 }
00096 else
00097 {
00098 ypt = yptNew;
00099 }
00100 }
00101
00102 return ypt;
00103 }
00104
00105
00106
00107
00108 kappatJayatillekeWallFunctionFvPatchScalarField::
00109 kappatJayatillekeWallFunctionFvPatchScalarField
00110 (
00111 const fvPatch& p,
00112 const DimensionedField<scalar, volMesh>& iF
00113 )
00114 :
00115 fixedValueFvPatchScalarField(p, iF),
00116 Prt_(0.85),
00117 Cmu_(0.09),
00118 kappa_(0.41),
00119 E_(9.8)
00120 {
00121 checkType();
00122 }
00123
00124
00125 kappatJayatillekeWallFunctionFvPatchScalarField::
00126 kappatJayatillekeWallFunctionFvPatchScalarField
00127 (
00128 const kappatJayatillekeWallFunctionFvPatchScalarField& ptf,
00129 const fvPatch& p,
00130 const DimensionedField<scalar, volMesh>& iF,
00131 const fvPatchFieldMapper& mapper
00132 )
00133 :
00134 fixedValueFvPatchScalarField(ptf, p, iF, mapper),
00135 Prt_(ptf.Prt_),
00136 Cmu_(ptf.Cmu_),
00137 kappa_(ptf.kappa_),
00138 E_(ptf.E_)
00139 {
00140 checkType();
00141 }
00142
00143
00144 kappatJayatillekeWallFunctionFvPatchScalarField::
00145 kappatJayatillekeWallFunctionFvPatchScalarField
00146 (
00147 const fvPatch& p,
00148 const DimensionedField<scalar, volMesh>& iF,
00149 const dictionary& dict
00150 )
00151 :
00152 fixedValueFvPatchScalarField(p, iF, dict),
00153 Prt_(readScalar(dict.lookup("Prt"))),
00154 Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
00155 kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
00156 E_(dict.lookupOrDefault<scalar>("E", 9.8))
00157 {
00158 checkType();
00159 }
00160
00161
00162 kappatJayatillekeWallFunctionFvPatchScalarField::
00163 kappatJayatillekeWallFunctionFvPatchScalarField
00164 (
00165 const kappatJayatillekeWallFunctionFvPatchScalarField& wfpsf
00166 )
00167 :
00168 fixedValueFvPatchScalarField(wfpsf),
00169 Prt_(wfpsf.Prt_),
00170 Cmu_(wfpsf.Cmu_),
00171 kappa_(wfpsf.kappa_),
00172 E_(wfpsf.E_)
00173 {
00174 checkType();
00175 }
00176
00177
00178 kappatJayatillekeWallFunctionFvPatchScalarField::
00179 kappatJayatillekeWallFunctionFvPatchScalarField
00180 (
00181 const kappatJayatillekeWallFunctionFvPatchScalarField& wfpsf,
00182 const DimensionedField<scalar, volMesh>& iF
00183 )
00184 :
00185 fixedValueFvPatchScalarField(wfpsf, iF),
00186 Prt_(wfpsf.Prt_),
00187 Cmu_(wfpsf.Cmu_),
00188 kappa_(wfpsf.kappa_),
00189 E_(wfpsf.E_)
00190 {
00191 checkType();
00192 }
00193
00194
00195
00196
00197 void kappatJayatillekeWallFunctionFvPatchScalarField::updateCoeffs()
00198 {
00199 if (updated())
00200 {
00201 return;
00202 }
00203
00204 const label patchI = patch().index();
00205
00206
00207 const RASModel& rasModel = db().lookupObject<RASModel>("RASProperties");
00208 const scalar Cmu25 = pow(Cmu_, 0.25);
00209 const scalarField& y = rasModel.y()[patchI];
00210 const scalarField& nuw = rasModel.nu().boundaryField()[patchI];
00211 const tmp<volScalarField> tk = rasModel.k();
00212 const volScalarField& k = tk();
00213
00214
00215 const scalar
00216 Pr(dimensionedScalar(rasModel.transport().lookup("Pr")).value());
00217
00218
00219 scalarField& kappatw = *this;
00220 forAll(kappatw, faceI)
00221 {
00222 label faceCellI = patch().faceCells()[faceI];
00223
00224
00225 scalar yPlus = Cmu25*sqrt(k[faceCellI])*y[faceI]/nuw[faceI];
00226
00227
00228 scalar Prat = Pr/Prt_;
00229
00230
00231 scalar P = Psmooth(Prat);
00232 scalar yPlusTherm = this->yPlusTherm(P, Prat);
00233
00234
00235 if (yPlus > yPlusTherm)
00236 {
00237 scalar nu = nuw[faceI];
00238 scalar kt = nu*(yPlus/(Prt_*(log(E_*yPlus)/kappa_ + P)) - 1/Pr);
00239 kappatw[faceI] = max(0.0, kt);
00240 }
00241 else
00242 {
00243 kappatw[faceI] = 0.0;
00244 }
00245 }
00246
00247 fixedValueFvPatchField<scalar>::updateCoeffs();
00248 }
00249
00250
00251 void kappatJayatillekeWallFunctionFvPatchScalarField::write(Ostream& os) const
00252 {
00253 fvPatchField<scalar>::write(os);
00254 os.writeKeyword("Prt") << Prt_ << token::END_STATEMENT << nl;
00255 os.writeKeyword("Cmu") << Cmu_ << token::END_STATEMENT << nl;
00256 os.writeKeyword("kappa") << kappa_ << token::END_STATEMENT << nl;
00257 os.writeKeyword("E") << E_ << token::END_STATEMENT << nl;
00258 writeEntry("value", os);
00259 }
00260
00261
00262
00263
00264 makePatchTypeField(fvPatchScalarField, kappatJayatillekeWallFunctionFvPatchScalarField);
00265
00266
00267
00268 }
00269 }
00270 }
00271
00272