FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

alphaSgsJayatillekeWallFunctionFvPatchScalarField.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 2008-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "alphaSgsJayatillekeWallFunctionFvPatchScalarField.H"
00027 #include <compressibleLESModels/LESModel.H>
00028 #include <finiteVolume/fvPatchFieldMapper.H>
00029 #include <finiteVolume/volFields.H>
00030 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00031 #include <finiteVolume/wallFvPatch.H>
00032 
00033 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00034 
00035 namespace Foam
00036 {
00037 namespace compressible
00038 {
00039 namespace LESModels
00040 {
00041 
00042 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00043 
00044 scalar alphaSgsJayatillekeWallFunctionFvPatchScalarField::maxExp_ = 50.0;
00045 scalar alphaSgsJayatillekeWallFunctionFvPatchScalarField::tolerance_ = 0.01;
00046 label alphaSgsJayatillekeWallFunctionFvPatchScalarField::maxIters_ = 10;
00047 
00048 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00049 
00050 void alphaSgsJayatillekeWallFunctionFvPatchScalarField::checkType()
00051 {
00052     if (!isA<wallFvPatch>(patch()))
00053     {
00054         FatalErrorIn
00055         (
00056             "alphaSgsJayatillekeWallFunctionFvPatchScalarField::checkType()"
00057         )
00058             << "Patch type for patch " << patch().name() << " must be wall\n"
00059             << "Current patch type is " << patch().type() << nl
00060             << exit(FatalError);
00061     }
00062 }
00063 
00064 
00065 scalar alphaSgsJayatillekeWallFunctionFvPatchScalarField::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 alphaSgsJayatillekeWallFunctionFvPatchScalarField::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 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00107 
00108 alphaSgsJayatillekeWallFunctionFvPatchScalarField::
00109 alphaSgsJayatillekeWallFunctionFvPatchScalarField
00110 (
00111     const fvPatch& p,
00112     const DimensionedField<scalar, volMesh>& iF
00113 )
00114 :
00115     fixedValueFvPatchScalarField(p, iF),
00116     Prt_(0.85),
00117     kappa_(0.41),
00118     E_(9.8),
00119     hsName_("hs")
00120 {
00121     checkType();
00122 }
00123 
00124 
00125 alphaSgsJayatillekeWallFunctionFvPatchScalarField::
00126 alphaSgsJayatillekeWallFunctionFvPatchScalarField
00127 (
00128     const alphaSgsJayatillekeWallFunctionFvPatchScalarField& 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     kappa_(ptf.kappa_),
00137     E_(ptf.E_),
00138     hsName_(ptf.hsName_)
00139 {}
00140 
00141 
00142 alphaSgsJayatillekeWallFunctionFvPatchScalarField::
00143 alphaSgsJayatillekeWallFunctionFvPatchScalarField
00144 (
00145     const fvPatch& p,
00146     const DimensionedField<scalar, volMesh>& iF,
00147     const dictionary& dict
00148 )
00149 :
00150     fixedValueFvPatchScalarField(p, iF, dict),
00151     Prt_(dict.lookupOrDefault<scalar>("Prt", 0.85)),
00152     kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
00153     E_(dict.lookupOrDefault<scalar>("E", 9.8)),
00154     hsName_(dict.lookupOrDefault<word>("hs", "hs"))
00155 {
00156     checkType();
00157 }
00158 
00159 
00160 alphaSgsJayatillekeWallFunctionFvPatchScalarField::
00161 alphaSgsJayatillekeWallFunctionFvPatchScalarField
00162 (
00163     const alphaSgsJayatillekeWallFunctionFvPatchScalarField& awfpsf
00164 )
00165 :
00166     fixedValueFvPatchScalarField(awfpsf),
00167     Prt_(awfpsf.Prt_),
00168     kappa_(awfpsf.kappa_),
00169     E_(awfpsf.E_),
00170     hsName_(awfpsf.hsName_)
00171 {
00172     checkType();
00173 }
00174 
00175 
00176 alphaSgsJayatillekeWallFunctionFvPatchScalarField::
00177 alphaSgsJayatillekeWallFunctionFvPatchScalarField
00178 (
00179     const alphaSgsJayatillekeWallFunctionFvPatchScalarField& awfpsf,
00180     const DimensionedField<scalar, volMesh>& iF
00181 )
00182 :
00183     fixedValueFvPatchScalarField(awfpsf, iF),
00184     Prt_(awfpsf.Prt_),
00185     kappa_(awfpsf.kappa_),
00186     E_(awfpsf.E_),
00187     hsName_(awfpsf.hsName_)
00188 {
00189     checkType();
00190 }
00191 
00192 
00193 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00194 
00195 void alphaSgsJayatillekeWallFunctionFvPatchScalarField::evaluate
00196 (
00197     const Pstream::commsTypes
00198 )
00199 {
00200     const LESModel& lesModel = db().lookupObject<LESModel>("LESProperties");
00201 
00202     // Field data
00203     const label patchI = patch().index();
00204 
00205     const scalarField& muw = lesModel.mu().boundaryField()[patchI];
00206     const scalarField muSgsw = lesModel.muSgs()().boundaryField()[patchI];
00207 
00208     const scalarField& alphaw = lesModel.alpha().boundaryField()[patchI];
00209     scalarField& alphaSgsw = *this;
00210 
00211     const fvPatchVectorField& Uw = lesModel.U().boundaryField()[patchI];
00212     const scalarField magUp = mag(Uw.patchInternalField() - Uw);
00213     const scalarField magGradUw = mag(Uw.snGrad());
00214 
00215     const scalarField& rhow = lesModel.rho().boundaryField()[patchI];
00216     const fvPatchScalarField& hw =
00217         patch().lookupPatchField<volScalarField, scalar>(hsName_);
00218 
00219     const scalarField& ry = patch().deltaCoeffs();
00220 
00221     // Heat flux [W/m2] - lagging alphaSgsw
00222     const scalarField qDot = (alphaw + alphaSgsw)*hw.snGrad();
00223 
00224     // Populate boundary values
00225     forAll(alphaSgsw, faceI)
00226     {
00227         // Calculate uTau using Newton-Raphson iteration
00228         scalar uTau =
00229             sqrt((muSgsw[faceI] + muw[faceI])/rhow[faceI]*magGradUw[faceI]);
00230 
00231         if (uTau > ROOTVSMALL)
00232         {
00233             label iter = 0;
00234             scalar err = GREAT;
00235 
00236             do
00237             {
00238                 scalar kUu = min(kappa_*magUp[faceI]/uTau, maxExp_);
00239                 scalar fkUu = exp(kUu) - 1.0 - kUu*(1.0 + 0.5*kUu);
00240 
00241                 scalar f =
00242                     - uTau/(ry[faceI]*muw[faceI]/rhow[faceI])
00243                     + magUp[faceI]/uTau
00244                     + 1.0/E_*(fkUu - 1.0/6.0*kUu*sqr(kUu));
00245 
00246                 scalar df =
00247                     - 1.0/(ry[faceI]*muw[faceI]/rhow[faceI])
00248                     - magUp[faceI]/sqr(uTau)
00249                     - 1.0/E_*kUu*fkUu/uTau;
00250 
00251                 scalar uTauNew = uTau - f/df;
00252                 err = mag((uTau - uTauNew)/uTau);
00253                 uTau = uTauNew;
00254 
00255             } while (uTau>VSMALL && err>tolerance_ && ++iter<maxIters_);
00256 
00257             scalar yPlus = uTau/ry[faceI]/(muw[faceI]/rhow[faceI]);
00258 
00259             // Molecular Prandtl number
00260             scalar Pr = muw[faceI]/alphaw[faceI];
00261 
00262             // Molecular-to-turbulenbt Prandtl number ratio
00263             scalar Prat = Pr/Prt_;
00264 
00265             // Thermal sublayer thickness
00266             scalar P = Psmooth(Prat);
00267             scalar yPlusTherm = this->yPlusTherm(P, Prat);
00268 
00269             // Evaluate new effective thermal diffusivity
00270             scalar alphaEff = 0.0;
00271             if (yPlus < yPlusTherm)
00272             {
00273                 scalar A = qDot[faceI]*rhow[faceI]*uTau/ry[faceI];
00274                 scalar B = qDot[faceI]*Pr*yPlus;
00275                 scalar C = Pr*0.5*rhow[faceI]*uTau*sqr(magUp[faceI]);
00276                 alphaEff = A/(B + C + VSMALL);
00277             }
00278             else
00279             {
00280                 scalar A = qDot[faceI]*rhow[faceI]*uTau/ry[faceI];
00281                 scalar B = qDot[faceI]*Prt_*(1.0/kappa_*log(E_*yPlus) + P);
00282                 scalar magUc = uTau/kappa_*log(E_*yPlusTherm) - mag(Uw[faceI]);
00283                 scalar C =
00284                     0.5*rhow[faceI]*uTau
00285                    *(Prt_*sqr(magUp[faceI]) + (Pr - Prt_)*sqr(magUc));
00286                 alphaEff = A/(B + C + VSMALL);
00287             }
00288 
00289             // Update turbulent thermal diffusivity
00290             alphaSgsw[faceI] = max(0.0, alphaEff - alphaw[faceI]);
00291 
00292             if (debug)
00293             {
00294                 Info<< "    uTau           = " << uTau << nl
00295                     << "    Pr             = " << Pr << nl
00296                     << "    Prt            = " << Prt_ << nl
00297                     << "    qDot           = " << qDot[faceI] << nl
00298                     << "    yPlus          = " << yPlus << nl
00299                     << "    yPlusTherm     = " << yPlusTherm << nl
00300                     << "    alphaEff       = " << alphaEff << nl
00301                     << "    alphaw         = " << alphaw[faceI] << nl
00302                     << "    alphaSgsw      = " << alphaSgsw[faceI] << nl
00303                     << endl;
00304             }
00305         }
00306         else
00307         {
00308             alphaSgsw[faceI] = 0.0;
00309         }
00310     }
00311 }
00312 
00313 
00314 void alphaSgsJayatillekeWallFunctionFvPatchScalarField::write(Ostream& os) const
00315 {
00316     fvPatchField<scalar>::write(os);
00317     os.writeKeyword("Prt") << Prt_ << token::END_STATEMENT << nl;
00318     os.writeKeyword("kappa") << kappa_ << token::END_STATEMENT << nl;
00319     os.writeKeyword("E") << E_ << token::END_STATEMENT << nl;
00320     os.writeKeyword("hs") << hsName_ << token::END_STATEMENT << nl;
00321     writeEntry("value", os);
00322 }
00323 
00324 
00325 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00326 
00327 makePatchTypeField
00328 (
00329     fvPatchScalarField,
00330     alphaSgsJayatillekeWallFunctionFvPatchScalarField
00331 );
00332 
00333 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00334 
00335 } // End namespace LESModels
00336 } // End namespace compressible
00337 } // End namespace Foam
00338 
00339 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines